libpappsomspp
Library for mass spectrometry
Loading...
Searching...
No Matches
pappso::PsmFeatures Class Reference

#include <psmfeatures.h>

Public Member Functions

 PsmFeatures (PrecisionPtr ms2precision, double minimumMz)
 compute psm features
 
 ~PsmFeatures ()
 
void setPeptideSpectrumCharge (const pappso::PeptideSp peptideSp, const MassSpectrum *p_spectrum, unsigned int parent_charge, unsigned int max_isotope_number)
 
double getIntensityOfMatchedIon (Enums::PeptideIon ion_type)
 get the sum of intensity of a specific ion
 
double getTotalIntensity () const
 sum of all peak intensities (matched or not)
 
double getTotalIntensityOfMatchedIons () const
 sum of matched peak intensities
 
std::size_t getNumberOfMatchedIons () const
 number of matched ions (peaks)
 
std::size_t countMatchedIonComplementPairs () const
 count the number of matched ion complement
 
double getTotalIntensityOfMatchedIonComplementPairs () const
 intensity of matched ion complement
 
const std::vector< std::pair< pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch > > & getPeakIonPairs () const
 
double getMatchedMzDiffMean () const
 get mean deviation of matched peak mass delta
 
double getMatchedMzDiffSd () const
 get standard deviation of matched peak mass delta
 
double getMaxIntensityMatchedIonComplementPairPrecursorMassDelta () const
 get the precursor mass delta of the maximum intensity pair of complement ions
 
std::size_t getMaxConsecutiveIon (Enums::PeptideIon ion_type)
 get the maximum consecutive fragments of one ion type
 
std::size_t getAaSequenceCoverage (Enums::PeptideIon ion_type)
 number of amino acid covered by matched ions
 
std::size_t getComplementPairsAaSequenceCoverage ()
 number of amino acid covered by matched complement pairs of ions
 
double getMaxIntensityPeakIonMatch (Enums::PeptideIon ion_type) const
 
double getIonPairPrecursorMassDelta (const std::pair< pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch > &ion_pair) const
 
LinearRegression getIonIsotopeLinearRegression () const
 

Private Member Functions

void findComplementIonPairs (const pappso::PeptideSp &peptideSp)
 

Private Attributes

std::shared_ptr< FilterResampleKeepGreatermsp_filterKeepGreater
 
std::shared_ptr< PeptideIsotopeSpectrumMatchmsp_peptideSpectrumMatch
 
pappso::PeptideSp msp_peptide
 
PrecisionPtr m_ms2precision
 
std::list< Enums::PeptideIonm_ionList
 
double m_spectrumSumIntensity
 
double m_precursorTheoreticalMz
 
double m_precursorTheoreticalMass
 
unsigned int m_parentCharge = 1
 
std::vector< std::pair< pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch > > m_peakIonPairs
 
double m_matchedMzDiffMean = 0
 
double m_matchedMzDiffMedian = 0
 
double m_matchedMzDiffSd = 0
 

Detailed Description

Todo
write docs

Definition at line 45 of file psmfeatures.h.

Constructor & Destructor Documentation

◆ PsmFeatures()

PsmFeatures::PsmFeatures ( PrecisionPtr ms2precision,
double minimumMz )

compute psm features

Parameters
ms2precisionprecision of mass measurements for MS2 fragments
minimumMzignore mz values under this threshold

Definition at line 35 of file psmfeatures.cpp.

36{
37
38 m_ms2precision = ms2precision;
39
42
43
44 msp_filterKeepGreater = std::make_shared<FilterResampleKeepGreater>(minimumMz);
45
46 // msp_filterChargeDeconvolution =
47 // std::make_shared<FilterChargeDeconvolution>(m_ms2precision);
48 // msp_filterMzExclusion = std::make_shared<FilterMzExclusion>(
49 // PrecisionFactory::getPrecisionPtrFractionInstance(m_ms2precision, 0.5));
50}
std::shared_ptr< FilterResampleKeepGreater > msp_filterKeepGreater
std::list< Enums::PeptideIon > m_ionList
PrecisionPtr m_ms2precision
@ y
Cter amino ions.
Definition types.h:295
@ b
Nter acylium ions.
Definition types.h:287

References pappso::Enums::b, m_ionList, m_ms2precision, msp_filterKeepGreater, and pappso::Enums::y.

◆ ~PsmFeatures()

PsmFeatures::~PsmFeatures ( )

Destructor

Definition at line 52 of file psmfeatures.cpp.

53{
54}

Member Function Documentation

◆ countMatchedIonComplementPairs()

std::size_t pappso::PsmFeatures::countMatchedIonComplementPairs ( ) const

count the number of matched ion complement

matched ion complement are ions with a sum compatible to the precursor mass

Definition at line 158 of file psmfeatures.cpp.

159{
160 return m_peakIonPairs.size();
161}
std::vector< std::pair< pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch > > m_peakIonPairs

References m_peakIonPairs.

◆ findComplementIonPairs()

void pappso::PsmFeatures::findComplementIonPairs ( const pappso::PeptideSp & peptideSp)
private

Definition at line 343 of file psmfeatures.cpp.

344{
345 std::size_t peptide_size = peptideSp.get()->size();
346 std::vector<PeakIonIsotopeMatch> ion_isotope_list(
347 msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList().begin(),
348 msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList().end());
349 for(const pappso::PeakIonIsotopeMatch &peak_ion_ext : ion_isotope_list)
350 {
351 if(peptideIonIsNter(peak_ion_ext.getPeptideIonType()))
352 {
353 auto it = findComplementIonType(
354 ion_isotope_list.begin(), ion_isotope_list.end(), peak_ion_ext, peptide_size);
355 if(it != ion_isotope_list.end())
356 { // contains the complementary ion
357 m_peakIonPairs.push_back({peak_ion_ext, *it});
358 }
359 }
360 }
361}
unsigned int size() const override
Definition peptide.cpp:215
std::shared_ptr< PeptideIsotopeSpectrumMatch > msp_peptideSpectrumMatch
bool peptideIonIsNter(Enums::PeptideIon ion_type)
tells if an ion is Nter
Definition peptide.cpp:85
std::vector< PeakIonIsotopeMatch >::iterator findComplementIonType(std::vector< PeakIonIsotopeMatch >::iterator begin, std::vector< PeakIonIsotopeMatch >::iterator end, const PeakIonIsotopeMatch &peak_ion, std::size_t peptide_size)
find the first element containing the complementary ion complementary ion of y1 is b(n-1) for instanc...

References pappso::findComplementIonType(), m_peakIonPairs, msp_peptideSpectrumMatch, pappso::peptideIonIsNter(), and pappso::Peptide::size().

Referenced by setPeptideSpectrumCharge().

◆ getAaSequenceCoverage()

std::size_t pappso::PsmFeatures::getAaSequenceCoverage ( Enums::PeptideIon ion_type)

number of amino acid covered by matched ions

Definition at line 254 of file psmfeatures.cpp.

255{
256 std::vector<bool> covered;
257 covered.resize(msp_peptide.get()->size(), false);
258
259 for(auto &peak : msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList())
260 {
261 if(peak.getPeptideNaturalIsotopeAverageSp().get()->getIsotopeNumber() == 0)
262 { // only count with monotisotopic
263 if(peak.getPeptideIonType() == ion_type)
264 {
265 covered[peak.getPeptideFragmentIonSp().get()->size() - 1] = true;
266 }
267 }
268 }
269 return std::count(covered.begin(), covered.end(), true);
270}
pappso::PeptideSp msp_peptide

References msp_peptide, and msp_peptideSpectrumMatch.

◆ getComplementPairsAaSequenceCoverage()

std::size_t pappso::PsmFeatures::getComplementPairsAaSequenceCoverage ( )

number of amino acid covered by matched complement pairs of ions

Definition at line 273 of file psmfeatures.cpp.

274{
275
276 std::vector<bool> covered;
277 covered.resize(msp_peptide.get()->size(), false);
278
279 for(auto &peak_pair : m_peakIonPairs)
280 {
281 std::size_t pos = peak_pair.first.getPeptideFragmentIonSp().get()->size() - 1;
282 covered[pos] = true;
283 covered[pos + 1] = true;
284 }
285 return std::count(covered.begin(), covered.end(), true);
286}

References m_peakIonPairs, and msp_peptide.

◆ getIntensityOfMatchedIon()

double PsmFeatures::getIntensityOfMatchedIon ( Enums::PeptideIon ion_type)

get the sum of intensity of a specific ion

Parameters
ion_typeion species (y, b, ...)

Definition at line 128 of file psmfeatures.cpp.

129{
130 double sum = 0;
131 for(const PeakIonMatch &peak_ion_match : *msp_peptideSpectrumMatch.get())
132 {
133 if(peak_ion_match.getPeptideIonType() == ion_type)
134 {
135 sum += peak_ion_match.getPeak().y;
136 }
137 }
138 return sum;
139}
@ sum
sum of intensities
Definition types.h:279

References msp_peptideSpectrumMatch.

◆ getIonIsotopeLinearRegression()

LinearRegression pappso::PsmFeatures::getIonIsotopeLinearRegression ( ) const

experimental use with caution

Definition at line 365 of file psmfeatures.cpp.

366{
367 std::size_t peptide_size = msp_peptide.get()->size();
368 std::vector<PeakIonIsotopeMatch> ion_isotope_list(
369 msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList().begin(),
370 msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList().end());
371
372 Trace scaterplot;
373
374 for(Enums::PeptideIon ion_type : m_ionList)
375 {
376 std::vector<double> mono_th_intensities(peptide_size, 0);
377 std::vector<double> isotope_th_intensities(peptide_size, 0);
378
379 std::vector<double> mono_exp_intensities(peptide_size, 0);
380 std::vector<double> isotope_exp_intensities(peptide_size, 0);
381 for(const PeakIonIsotopeMatch &peak_ion_match : ion_isotope_list)
382 {
383 if(peak_ion_match.getPeptideIonType() == ion_type)
384 {
385 std::size_t vector_position =
386 peak_ion_match.getPeptideFragmentIonSp().get()->size() - 1;
387 PeptideNaturalIsotopeAverageSp iso_average_sp =
388 peak_ion_match.getPeptideNaturalIsotopeAverageSp();
389 if(iso_average_sp.get()->getIsotopeNumber() == 0)
390 {
391 mono_th_intensities[vector_position] = iso_average_sp.get()->getIntensityRatio();
392 mono_exp_intensities[vector_position] = peak_ion_match.getPeak().y;
393 }
394 else if(iso_average_sp.get()->getIsotopeNumber() == 1)
395 {
396 isotope_th_intensities[vector_position] =
397 iso_average_sp.get()->getIntensityRatio();
398 isotope_exp_intensities[vector_position] = peak_ion_match.getPeak().y;
399 }
400 }
401 }
402
403 for(std::size_t i = 0; i < mono_th_intensities.size(); i++)
404 {
405 if((mono_th_intensities[i] != 0) && (isotope_th_intensities[i] != 0))
406 {
407 DataPoint xy(mono_th_intensities[i] / isotope_th_intensities[i],
408 mono_exp_intensities[i] / isotope_exp_intensities[i]);
409 scaterplot.push_back(xy);
410 }
411 }
412 }
413
414 scaterplot.sortX();
415
416 LinearRegression linear_regression(scaterplot);
417
418 return linear_regression;
419}
void sortX(Enums::SortOrder sort_order=Enums::SortOrder::ascending)
Definition trace.cpp:1039
PeptideIon
Enums::PeptideIon enum defines all types of ions (Nter or Cter)
Definition types.h:286
std::shared_ptr< const PeptideNaturalIsotopeAverage > PeptideNaturalIsotopeAverageSp

References pappso::PeptideNaturalIsotopeAverage::getIntensityRatio(), pappso::PeptideNaturalIsotopeAverage::getIsotopeNumber(), m_ionList, msp_peptide, msp_peptideSpectrumMatch, and pappso::Trace::sortX().

◆ getIonPairPrecursorMassDelta()

double pappso::PsmFeatures::getIonPairPrecursorMassDelta ( const std::pair< pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch > & ion_pair) const

Definition at line 327 of file psmfeatures.cpp.

329{
330 qDebug() << m_precursorTheoreticalMz << " " << ion_pair.first.getPeak().x << " "
331 << ion_pair.second.getPeak().x << " " << ion_pair.second.getCharge() << " "
332 << ion_pair.first.getCharge() << " " << m_parentCharge;
333 double diff = (m_precursorTheoreticalMass + (MHPLUS * ion_pair.first.getCharge())) /
334 ion_pair.first.getCharge();
335
336
337 return (diff - (ion_pair.first.getPeak().x + ion_pair.second.getPeak().x -
338 ((MHPLUS * ion_pair.first.getCharge())) / ion_pair.first.getCharge()));
339}
double m_precursorTheoreticalMz
unsigned int m_parentCharge
double m_precursorTheoreticalMass
const pappso_double MHPLUS(1.007276466879)

References m_parentCharge, m_precursorTheoreticalMass, m_precursorTheoreticalMz, and pappso::MHPLUS().

Referenced by getMaxIntensityMatchedIonComplementPairPrecursorMassDelta().

◆ getMatchedMzDiffMean()

double pappso::PsmFeatures::getMatchedMzDiffMean ( ) const

get mean deviation of matched peak mass delta

Definition at line 189 of file psmfeatures.cpp.

190{
191 return m_matchedMzDiffMean;
192}

References m_matchedMzDiffMean.

◆ getMatchedMzDiffSd()

double pappso::PsmFeatures::getMatchedMzDiffSd ( ) const

get standard deviation of matched peak mass delta

Definition at line 183 of file psmfeatures.cpp.

184{
185 return m_matchedMzDiffSd;
186}

References m_matchedMzDiffSd.

◆ getMaxConsecutiveIon()

std::size_t pappso::PsmFeatures::getMaxConsecutiveIon ( Enums::PeptideIon ion_type)

get the maximum consecutive fragments of one ion type

Parameters
ion_typeion species (y, b, ...)

Definition at line 202 of file psmfeatures.cpp.

203{
204 std::size_t max = 0;
205 auto peak_ion_match_list = msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList();
206
207 peak_ion_match_list.erase(
208 std::remove_if(peak_ion_match_list.begin(),
209 peak_ion_match_list.end(),
210 [ion_type](const PeakIonIsotopeMatch &a) {
211 if(a.getPeptideIonType() != ion_type)
212 return true;
213 if(a.getPeptideNaturalIsotopeAverageSp().get()->getIsotopeNumber() > 0)
214 return true;
215 return false;
216 }),
217 peak_ion_match_list.end());
218
219 peak_ion_match_list.sort([](const PeakIonIsotopeMatch &a, const PeakIonIsotopeMatch &b) {
220 if(a.getCharge() < b.getCharge())
221 return true;
222 if(a.getPeptideIonType() < b.getPeptideIonType())
223 return true;
224 if(a.getPeptideFragmentIonSp().get()->size() < b.getPeptideFragmentIonSp().get()->size())
225 return true;
226 return false;
227 });
228
229 unsigned int charge = 0;
230 std::size_t size = 0;
231 std::size_t count = 0;
232 for(std::list<PeakIonIsotopeMatch>::iterator it = peak_ion_match_list.begin();
233 it != peak_ion_match_list.end();
234 it++)
235 {
236 qDebug() << it->toString() << max << " " << it->getPeak().x << " "
237 << it->getPeptideNaturalIsotopeAverageSp().get()->getIsotopeNumber();
238 count++;
239 if((charge != it->getCharge()) || (size != (it->getPeptideFragmentIonSp().get()->size() - 1)))
240 {
241 count = 1;
242 charge = it->getCharge();
243 }
244 if(max < count)
245 max = count;
246
247 size = it->getPeptideFragmentIonSp().get()->size();
248 }
249
250 return max;
251}
@ max
maximum of intensities
Definition types.h:280

References pappso::a, pappso::b, and msp_peptideSpectrumMatch.

◆ getMaxIntensityMatchedIonComplementPairPrecursorMassDelta()

double pappso::PsmFeatures::getMaxIntensityMatchedIonComplementPairPrecursorMassDelta ( ) const

get the precursor mass delta of the maximum intensity pair of complement ions

Definition at line 309 of file psmfeatures.cpp.

310{
311 auto peak_it = std::max_element(
312 m_peakIonPairs.begin(),
313 m_peakIonPairs.end(),
314 [](const std::pair<pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch> &a,
315 const std::pair<pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch> &b) {
316 return ((a.first.getPeak().y + a.second.getPeak().y) <
317 (b.first.getPeak().y + b.second.getPeak().y));
318 });
319
320 if(peak_it == m_peakIonPairs.end())
321 return 0;
322
323 return getIonPairPrecursorMassDelta(*peak_it);
324}
double getIonPairPrecursorMassDelta(const std::pair< pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch > &ion_pair) const

References pappso::a, pappso::b, getIonPairPrecursorMassDelta(), and m_peakIonPairs.

◆ getMaxIntensityPeakIonMatch()

double pappso::PsmFeatures::getMaxIntensityPeakIonMatch ( Enums::PeptideIon ion_type) const

Definition at line 290 of file psmfeatures.cpp.

291{
292 std::list<pappso::PeakIonIsotopeMatch> peak_ion_type =
293 msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList();
294
295 peak_ion_type.remove_if(
296 [ion_type](const PeakIonIsotopeMatch &a) { return (a.getPeptideIonType() != ion_type); });
297 auto peak_it = std::max_element(peak_ion_type.begin(),
298 peak_ion_type.end(),
299 [](const PeakIonIsotopeMatch &a, const PeakIonIsotopeMatch &b) {
300 return (a.getPeak().y < b.getPeak().y);
301 });
302
303 if(peak_it == peak_ion_type.end())
304 return 0;
305 return peak_it->getPeak().y;
306}

References pappso::a, pappso::b, and msp_peptideSpectrumMatch.

◆ getNumberOfMatchedIons()

std::size_t pappso::PsmFeatures::getNumberOfMatchedIons ( ) const

number of matched ions (peaks)

Definition at line 196 of file psmfeatures.cpp.

197{
198 return msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList().size();
199}

References msp_peptideSpectrumMatch.

◆ getPeakIonPairs()

const std::vector< std::pair< pappso::PeakIonIsotopeMatch, pappso::PeakIonIsotopeMatch > > & pappso::PsmFeatures::getPeakIonPairs ( ) const

Definition at line 164 of file psmfeatures.cpp.

165{
166 return m_peakIonPairs;
167}

References m_peakIonPairs.

◆ getTotalIntensity()

double PsmFeatures::getTotalIntensity ( ) const

sum of all peak intensities (matched or not)

Definition at line 152 of file psmfeatures.cpp.

153{
155}
double m_spectrumSumIntensity

References m_spectrumSumIntensity.

◆ getTotalIntensityOfMatchedIonComplementPairs()

double pappso::PsmFeatures::getTotalIntensityOfMatchedIonComplementPairs ( ) const

intensity of matched ion complement

Definition at line 170 of file psmfeatures.cpp.

171{
172
173 double sum = 0;
174 for(auto &peak_pairs : m_peakIonPairs)
175 {
176 sum += peak_pairs.first.getPeak().y;
177 sum += peak_pairs.second.getPeak().y;
178 }
179 return sum;
180}

References m_peakIonPairs.

◆ getTotalIntensityOfMatchedIons()

double PsmFeatures::getTotalIntensityOfMatchedIons ( ) const

sum of matched peak intensities

Definition at line 141 of file psmfeatures.cpp.

142{
143 double sum = 0;
144 for(const PeakIonMatch &peak_ion_match : *msp_peptideSpectrumMatch.get())
145 {
146 sum += peak_ion_match.getPeak().y;
147 }
148 return sum;
149}

References msp_peptideSpectrumMatch.

◆ setPeptideSpectrumCharge()

void PsmFeatures::setPeptideSpectrumCharge ( const pappso::PeptideSp peptideSp,
const MassSpectrum * p_spectrum,
unsigned int parent_charge,
unsigned int max_isotope_number )

Definition at line 57 of file psmfeatures.cpp.

61{
62 m_peakIonPairs.clear();
63 msp_peptide = peptideSp;
64 MassSpectrum spectrum(*p_spectrum);
65 msp_filterKeepGreater.get()->filter(spectrum);
66 // msp_filterChargeDeconvolution.get()->filter(spectrum);
67 // msp_filterMzExclusion.get()->filter(spectrum);
68
69 msp_peptideSpectrumMatch = std::make_shared<PeptideIsotopeSpectrumMatch>(
70 spectrum, peptideSp, parent_charge, m_ms2precision, m_ionList, max_isotope_number, 1);
71
72 msp_peptideSpectrumMatch.get()->dropPeaksLackingMonoisotope();
73 m_spectrumSumIntensity = spectrum.sumY();
74
75
76 qDebug() << " accumulate";
77 std::vector<double> delta_list;
78
79
80 // TODO compute number of matched complementary peaks having m/z compatible
81 // with the precursor
82
83 m_precursorTheoreticalMz = peptideSp.get()->getMz(parent_charge);
84 m_precursorTheoreticalMass = peptideSp.get()->getMass();
85 m_parentCharge = parent_charge;
86
87
88 findComplementIonPairs(peptideSp);
89
90
91 for(const pappso::PeakIonIsotopeMatch &peak_ion :
92 msp_peptideSpectrumMatch.get()->getPeakIonIsotopeMatchList())
93 {
94 delta_list.push_back(peak_ion.getPeptideFragmentIonSp().get()->getMz(peak_ion.getCharge()) -
95 peak_ion.getPeak().x);
96 }
97 pappso::pappso_double sum = std::accumulate(delta_list.begin(), delta_list.end(), 0);
98
99 qDebug() << " delta_list.size()=" << delta_list.size();
103 if(delta_list.size() > 0)
104 {
105 m_matchedMzDiffMean = sum / ((pappso::pappso_double)delta_list.size());
106
107 std::sort(delta_list.begin(), delta_list.end());
108 m_matchedMzDiffMedian = delta_list[(delta_list.size() / 2)];
109
110
111 qDebug() << " sd";
113 for(pappso::pappso_double val : delta_list)
114 {
115 // sd = sd + ((val - mean) * (val - mean));
116 m_matchedMzDiffSd += std::pow((val - m_matchedMzDiffMean), 2);
117 }
118 m_matchedMzDiffSd = m_matchedMzDiffSd / delta_list.size();
120 }
121 else
122 {
123 }
124}
virtual pappso_double getMz(unsigned int charge) const final
Definition ion.cpp:45
pappso_double getMass()
Definition peptide.cpp:322
void findComplementIonPairs(const pappso::PeptideSp &peptideSp)
double pappso_double
A type definition for doubles.
Definition types.h:61

References findComplementIonPairs(), pappso::Peptide::getMass(), pappso::Ion::getMz(), m_ionList, m_matchedMzDiffMean, m_matchedMzDiffMedian, m_matchedMzDiffSd, m_ms2precision, m_parentCharge, m_peakIonPairs, m_precursorTheoreticalMass, m_precursorTheoreticalMz, m_spectrumSumIntensity, msp_filterKeepGreater, msp_peptide, msp_peptideSpectrumMatch, and pappso::Trace::sumY().

Member Data Documentation

◆ m_ionList

std::list<Enums::PeptideIon> pappso::PsmFeatures::m_ionList
private

◆ m_matchedMzDiffMean

double pappso::PsmFeatures::m_matchedMzDiffMean = 0
private

Definition at line 159 of file psmfeatures.h.

Referenced by getMatchedMzDiffMean(), and setPeptideSpectrumCharge().

◆ m_matchedMzDiffMedian

double pappso::PsmFeatures::m_matchedMzDiffMedian = 0
private

Definition at line 160 of file psmfeatures.h.

Referenced by setPeptideSpectrumCharge().

◆ m_matchedMzDiffSd

double pappso::PsmFeatures::m_matchedMzDiffSd = 0
private

Definition at line 161 of file psmfeatures.h.

Referenced by getMatchedMzDiffSd(), and setPeptideSpectrumCharge().

◆ m_ms2precision

PrecisionPtr pappso::PsmFeatures::m_ms2precision
private

Definition at line 148 of file psmfeatures.h.

Referenced by PsmFeatures(), and setPeptideSpectrumCharge().

◆ m_parentCharge

unsigned int pappso::PsmFeatures::m_parentCharge = 1
private

Definition at line 155 of file psmfeatures.h.

Referenced by getIonPairPrecursorMassDelta(), and setPeptideSpectrumCharge().

◆ m_peakIonPairs

◆ m_precursorTheoreticalMass

double pappso::PsmFeatures::m_precursorTheoreticalMass
private

Definition at line 154 of file psmfeatures.h.

Referenced by getIonPairPrecursorMassDelta(), and setPeptideSpectrumCharge().

◆ m_precursorTheoreticalMz

double pappso::PsmFeatures::m_precursorTheoreticalMz
private

Definition at line 153 of file psmfeatures.h.

Referenced by getIonPairPrecursorMassDelta(), and setPeptideSpectrumCharge().

◆ m_spectrumSumIntensity

double pappso::PsmFeatures::m_spectrumSumIntensity
private

Definition at line 151 of file psmfeatures.h.

Referenced by getTotalIntensity(), and setPeptideSpectrumCharge().

◆ msp_filterKeepGreater

std::shared_ptr<FilterResampleKeepGreater> pappso::PsmFeatures::msp_filterKeepGreater
private

Definition at line 143 of file psmfeatures.h.

Referenced by PsmFeatures(), and setPeptideSpectrumCharge().

◆ msp_peptide

◆ msp_peptideSpectrumMatch


The documentation for this class was generated from the following files: