Home  · Classes  · Annotated Classes  · Modules  · Members  · Namespaces  · Related Pages
Static Public Member Functions | List of all members
OPXLHelper Class Reference

The OPXLHelper class contains functions needed by OpenPepXL and OpenPepXLLF to reduce duplicated code. More...

#include <OpenMS/ANALYSIS/XLMS/OPXLHelper.h>

Static Public Member Functions

static std::vector< OPXLDataStructs::XLPrecursorenumerateCrossLinksAndMasses (const std::vector< OPXLDataStructs::AASeqWithMass > &peptides, double cross_link_mass_light, const DoubleList &cross_link_mass_mono_link, const StringList &cross_link_residue1, const StringList &cross_link_residue2, std::vector< double > &spectrum_precursors, double precursor_mass_tolerance, bool precursor_mass_tolerance_unit_ppm)
 Enumerates precursor masses for all candidates in an XL-MS search. More...
 
static std::vector< ResidueModificationgetModificationsFromStringList (StringList modNames)
 A helper function, that turns a StringList with modification names into a vector of ResidueModifications. More...
 
static std::vector< OPXLDataStructs::AASeqWithMassdigestDatabase (std::vector< FASTAFile::FASTAEntry > fasta_db, EnzymaticDigestion digestor, Size min_peptide_length, StringList cross_link_residue1, StringList cross_link_residue2, std::vector< ResidueModification > fixed_modifications, std::vector< ResidueModification > variable_modifications, Size max_variable_mods_per_peptide, Size count_proteins=0, Size count_peptides=0, bool n_term_linker=false, bool c_term_linker=false)
 Digests a database with the given EnzymaticDigestion settings and precomputes masses for all peptides. More...
 
static std::vector< OPXLDataStructs::ProteinProteinCrossLinkbuildCandidates (const std::vector< OPXLDataStructs::XLPrecursor > &candidates, const std::vector< OPXLDataStructs::AASeqWithMass > &peptide_masses, const StringList &cross_link_residue1, const StringList &cross_link_residue2, double cross_link_mass, const DoubleList &cross_link_mass_mono_link, double precursor_mass, double allowed_error, String cross_link_name, bool n_term_linker, bool c_term_linker)
 Builds specific cross-link candidates with all possible combinations of linked positions from peptide pairs. Used to build candidates for the precursor mass window of a single MS2 spectrum. More...
 
static void buildFragmentAnnotations (std::vector< PeptideHit::PeakAnnotation > &frag_annotations, const std::vector< std::pair< Size, Size > > &matching, const PeakSpectrum &theoretical_spectrum, const PeakSpectrum &experiment_spectrum)
 Fills up the given FragmentAnnotation vector with annotations from a theoretical spectrum4. More...
 
static void buildPeptideIDs (std::vector< PeptideIdentification > &peptide_ids, const std::vector< OPXLDataStructs::CrossLinkSpectrumMatch > &top_csms_spectrum, std::vector< std::vector< OPXLDataStructs::CrossLinkSpectrumMatch > > &all_top_csms, Size all_top_csms_current_index, const PeakMap &spectra, Size scan_index, Size scan_index_heavy)
 Builds PeptideIdentifications and PeptideHits. More...
 

Detailed Description

The OPXLHelper class contains functions needed by OpenPepXL and OpenPepXLLF to reduce duplicated code.

Member Function Documentation

◆ buildCandidates()

static std::vector<OPXLDataStructs::ProteinProteinCrossLink> buildCandidates ( const std::vector< OPXLDataStructs::XLPrecursor > &  candidates,
const std::vector< OPXLDataStructs::AASeqWithMass > &  peptide_masses,
const StringList cross_link_residue1,
const StringList cross_link_residue2,
double  cross_link_mass,
const DoubleList cross_link_mass_mono_link,
double  precursor_mass,
double  allowed_error,
String  cross_link_name,
bool  n_term_linker,
bool  c_term_linker 
)
static

Builds specific cross-link candidates with all possible combinations of linked positions from peptide pairs. Used to build candidates for the precursor mass window of a single MS2 spectrum.

Parameters
candidatesThe XLPrecursors containing indices of two peptides
peptide_massesThe digested peptide database, that the indices in the XLPrecursors refer to
cross_link_residue1A list of residues, to which the first side of the linker can react
cross_link_residue2A list of residues, to which the second side of the linker can react
cross_link_massmass of the cross-linker, only the light one if a labeled linker is used
cross_link_mass_mono_linkA list of possible masses for the cross-link, if it is attached to a peptide on one side
precursor_massThe precursor mass of the experimental spectrum (used to filter out certain candidates, e.g. mono- and loop-links have a different mass)
allowed_errorThe maximal precursor mass error in Da
cross_link_nameThe name of the cross-linker
n_term_linkerTrue, if the cross-linker can react with the N-terminal of a protein
c_term_linkerTrue, if the cross-linker can react with the C-terminal of a protein
Returns
A vector of ProteinProteinCrossLink candidates containing all necessary information to generate theoretical spectra

◆ buildFragmentAnnotations()

static void buildFragmentAnnotations ( std::vector< PeptideHit::PeakAnnotation > &  frag_annotations,
const std::vector< std::pair< Size, Size > > &  matching,
const PeakSpectrum theoretical_spectrum,
const PeakSpectrum experiment_spectrum 
)
static

Fills up the given FragmentAnnotation vector with annotations from a theoretical spectrum4.

This function takes an alignment of a theoretical spectrum with meta information and an experimental spectrum and builds annotations taking the MZ and intensity values from the experimental spectrum and the ion names and charges from the theoretical spectrum to annotate matched experimental peaks.

Parameters
frag_annotationsThe vector to fill. Does not have to be empty, as annotations from several alignments can just be added on to the same vector.
matchingThe alignment between the two spectra
theoretical_spectrumThe theoretical spectrum with meta information
experiment_spectrumThe experimental specrum

◆ buildPeptideIDs()

static void buildPeptideIDs ( std::vector< PeptideIdentification > &  peptide_ids,
const std::vector< OPXLDataStructs::CrossLinkSpectrumMatch > &  top_csms_spectrum,
std::vector< std::vector< OPXLDataStructs::CrossLinkSpectrumMatch > > &  all_top_csms,
Size  all_top_csms_current_index,
const PeakMap spectra,
Size  scan_index,
Size  scan_index_heavy 
)
static

Builds PeptideIdentifications and PeptideHits.

Parameters
peptide_idsThe vector of PeptideIdentifications for the whole experiment. The created PepIds will be pushed on this one.
top_csms_spectrumAll CrossLinkSpectrumMatches from the current spectrum to be written out
all_top_csmsA vector of all CrossLinkSpectrumMatches of the experiment, that is also extended in this function
all_top_csms_current_indexThe index of the current spectrum in all_top_csms (some spectra have no matches, so this is not equal to the spectrum index)
spectraThe searched spectra as a PeakMap
scan_indexThe index of the current spectrum
scan_index_heavyThe index of the heavy spectrum in a spectrum pair with labeled linkers

◆ digestDatabase()

static std::vector<OPXLDataStructs::AASeqWithMass> digestDatabase ( std::vector< FASTAFile::FASTAEntry fasta_db,
EnzymaticDigestion  digestor,
Size  min_peptide_length,
StringList  cross_link_residue1,
StringList  cross_link_residue2,
std::vector< ResidueModification fixed_modifications,
std::vector< ResidueModification variable_modifications,
Size  max_variable_mods_per_peptide,
Size  count_proteins = 0,
Size  count_peptides = 0,
bool  n_term_linker = false,
bool  c_term_linker = false 
)
static

Digests a database with the given EnzymaticDigestion settings and precomputes masses for all peptides.

Also keeps track of the peptides at protein terminals and builds peptide candidates with all possible modification patterns according to the parameters.

Parameters
fasta_dbThe protein database
digestorThe object containing the digestion settings, e.g. the enzyme
min_peptide_lengthThe minimal peptide length for the digestion
cross_link_residue1A list of residues, to which the first side of the linker can react
cross_link_residue2A list of residues, to which the second side of the linker can react
fixed_modificationsThe list of fixed modifications
variable_modificationsThe list of variable modifications
max_variable_mods_per_peptideThe maximal number of variable modifications per peptide
count_proteinsA variable to keep track of the number of proteins in the database. Should be an externally declared variable and = 0 when calling this function.
count_peptidesA variable to keep track of the number of peptides after digestion. Should be an externally declared variable and = 0 when calling this function.
n_term_linkerTrue, if the cross-linker can react with the N-terminal of a protein
c_term_linkerTrue, if the cross-linker can react with the C-terminal of a protein
Returns
A vector of AASeqWithMass containing the peptides, their masses and information about terminal peptides

◆ enumerateCrossLinksAndMasses()

static std::vector<OPXLDataStructs::XLPrecursor> enumerateCrossLinksAndMasses ( const std::vector< OPXLDataStructs::AASeqWithMass > &  peptides,
double  cross_link_mass_light,
const DoubleList cross_link_mass_mono_link,
const StringList cross_link_residue1,
const StringList cross_link_residue2,
std::vector< double > &  spectrum_precursors,
double  precursor_mass_tolerance,
bool  precursor_mass_tolerance_unit_ppm 
)
static

Enumerates precursor masses for all candidates in an XL-MS search.

Parameters
peptidesThe peptides with precomputed masses from the digestDatabase function
cross_link_mass_lightmass of the cross-linker, only the light one if a labeled linker is used
cross_link_mass_mono_linkA list of possible masses for the cross-link, if it is attached to a peptide on one side
cross_link_residue1A list of residues, to which the first side of the linker can react
cross_link_residue2A list of residues, to which the second side of the linker can react
spectrum_precursorsA vector of all MS2 precursor masses of the searched spectra. Used to filter out candidates.
precursor_mass_toleranceThe precursor mass tolerance
precursor_mass_tolerance_unit_ppmThe unit of the precursor mass tolerance ("Da" or "ppm")
Returns
A vector of XLPrecursors containing all possible candidate cross-links

◆ getModificationsFromStringList()

static std::vector<ResidueModification> getModificationsFromStringList ( StringList  modNames)
static

A helper function, that turns a StringList with modification names into a vector of ResidueModifications.

Parameters
modNamesThe list of modification names
Returns
A vector of modifications

OpenMS / TOPP release 2.3.0 Documentation generated on Tue Jan 9 2018 18:22:11 using doxygen 1.8.13