OpenMS  2.8.0
Classes | Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Friends | List of all members
IdentificationData Class Reference

Representation of spectrum identification results and associated data. More...

#include <OpenMS/METADATA/ID/IdentificationData.h>

Inheritance diagram for IdentificationData:
[legend]
Collaboration diagram for IdentificationData:
[legend]

Classes

struct  ModifyMultiIndexAddProcessingStep
 Helper functor for adding processing steps to elements in a @t boost::multi_index_container structure. More...
 
struct  ModifyMultiIndexAddScore
 Helper functor for adding scores to elements in a boost::multi_index_container structure. More...
 
struct  ModifyMultiIndexRemoveParentMatches
 Helper functor for removing invalid parent matches from elements in a boost::multi_index_container structure. More...
 
struct  RefTranslator
 structure that maps references of corresponding objects after copying More...
 

Public Types

using MoleculeType = IdentificationDataInternal::MoleculeType
 
using MassType = IdentificationDataInternal::MassType
 
using InputFile = IdentificationDataInternal::InputFile
 
using InputFiles = IdentificationDataInternal::InputFiles
 
using InputFileRef = IdentificationDataInternal::InputFileRef
 
using ProcessingSoftware = IdentificationDataInternal::ProcessingSoftware
 
using ProcessingSoftwares = IdentificationDataInternal::ProcessingSoftwares
 
using ProcessingSoftwareRef = IdentificationDataInternal::ProcessingSoftwareRef
 
using ProcessingStep = IdentificationDataInternal::ProcessingStep
 
using ProcessingSteps = IdentificationDataInternal::ProcessingSteps
 
using ProcessingStepRef = IdentificationDataInternal::ProcessingStepRef
 
using DBSearchParam = IdentificationDataInternal::DBSearchParam
 
using DBSearchParams = IdentificationDataInternal::DBSearchParams
 
using SearchParamRef = IdentificationDataInternal::SearchParamRef
 
using DBSearchSteps = IdentificationDataInternal::DBSearchSteps
 
using ScoreType = IdentificationDataInternal::ScoreType
 
using ScoreTypes = IdentificationDataInternal::ScoreTypes
 
using ScoreTypeRef = IdentificationDataInternal::ScoreTypeRef
 
using ScoredProcessingResult = IdentificationDataInternal::ScoredProcessingResult
 
using AppliedProcessingStep = IdentificationDataInternal::AppliedProcessingStep
 
using AppliedProcessingSteps = IdentificationDataInternal::AppliedProcessingSteps
 
using Observation = IdentificationDataInternal::Observation
 
using Observations = IdentificationDataInternal::Observations
 
using ObservationRef = IdentificationDataInternal::ObservationRef
 
using ParentSequence = IdentificationDataInternal::ParentSequence
 
using ParentSequences = IdentificationDataInternal::ParentSequences
 
using ParentSequenceRef = IdentificationDataInternal::ParentSequenceRef
 
using ParentMatch = IdentificationDataInternal::ParentMatch
 
using ParentMatches = IdentificationDataInternal::ParentMatches
 
using IdentifiedPeptide = IdentificationDataInternal::IdentifiedPeptide
 
using IdentifiedPeptides = IdentificationDataInternal::IdentifiedPeptides
 
using IdentifiedPeptideRef = IdentificationDataInternal::IdentifiedPeptideRef
 
using IdentifiedCompound = IdentificationDataInternal::IdentifiedCompound
 
using IdentifiedCompounds = IdentificationDataInternal::IdentifiedCompounds
 
using IdentifiedCompoundRef = IdentificationDataInternal::IdentifiedCompoundRef
 
using IdentifiedOligo = IdentificationDataInternal::IdentifiedOligo
 
using IdentifiedOligos = IdentificationDataInternal::IdentifiedOligos
 
using IdentifiedOligoRef = IdentificationDataInternal::IdentifiedOligoRef
 
using IdentifiedMolecule = IdentificationDataInternal::IdentifiedMolecule
 
using PeakAnnotations = IdentificationDataInternal::PeakAnnotations
 
using Adducts = IdentificationDataInternal::Adducts
 
using AdductRef = IdentificationDataInternal::AdductRef
 
using AdductOpt = IdentificationDataInternal::AdductOpt
 
using ObservationMatch = IdentificationDataInternal::ObservationMatch
 
using ObservationMatches = IdentificationDataInternal::ObservationMatches
 
using ObservationMatchRef = IdentificationDataInternal::ObservationMatchRef
 
using ObservationMatchGroup = IdentificationDataInternal::ObservationMatchGroup
 
using ObservationMatchGroups = IdentificationDataInternal::ObservationMatchGroups
 
using MatchGroupRef = IdentificationDataInternal::MatchGroupRef
 
using ParentGroup = IdentificationDataInternal::ParentGroup
 
using ParentGroups = IdentificationDataInternal::ParentGroups
 
using ParentGroupRef = IdentificationDataInternal::ParentGroupRef
 
using ParentGroupSet = IdentificationDataInternal::ParentGroupSet
 
using ParentGroupSets = IdentificationDataInternal::ParentGroupSets
 
using AddressLookup = std::unordered_set< uintptr_t >
 

Public Member Functions

 IdentificationData ()
 Default constructor. More...
 
 IdentificationData (const IdentificationData &other)
 Copy constructor. More...
 
 IdentificationData (IdentificationData &&other) noexcept
 Move constructor. More...
 
InputFileRef registerInputFile (const InputFile &file)
 Register an input file. More...
 
ProcessingSoftwareRef registerProcessingSoftware (const ProcessingSoftware &software)
 Register data processing software. More...
 
SearchParamRef registerDBSearchParam (const DBSearchParam &param)
 Register database search parameters. More...
 
ProcessingStepRef registerProcessingStep (const ProcessingStep &step)
 Register a data processing step. More...
 
ProcessingStepRef registerProcessingStep (const ProcessingStep &step, SearchParamRef search_ref)
 Register a database search step with associated parameters. More...
 
ScoreTypeRef registerScoreType (const ScoreType &score)
 Register a score type. More...
 
ObservationRef registerObservation (const Observation &obs)
 Register an observation (e.g. MS2 spectrum or feature) More...
 
ParentSequenceRef registerParentSequence (const ParentSequence &parent)
 Register a parent sequence (e.g. protein or intact RNA) More...
 
void registerParentGroupSet (const ParentGroupSet &groups)
 Register a grouping of parent sequences (e.g. protein inference result) More...
 
IdentifiedPeptideRef registerIdentifiedPeptide (const IdentifiedPeptide &peptide)
 Register an identified peptide. More...
 
IdentifiedCompoundRef registerIdentifiedCompound (const IdentifiedCompound &compound)
 Register an identified compound (small molecule) More...
 
IdentifiedOligoRef registerIdentifiedOligo (const IdentifiedOligo &oligo)
 Register an identified RNA oligonucleotide. More...
 
AdductRef registerAdduct (const AdductInfo &adduct)
 Register an adduct. More...
 
ObservationMatchRef registerObservationMatch (const ObservationMatch &match)
 Register an observation match (e.g. peptide-spectrum match) More...
 
MatchGroupRef registerObservationMatchGroup (const ObservationMatchGroup &group)
 Register a group of observation matches that belong together. More...
 
const InputFilesgetInputFiles () const
 Return the registered input files (immutable) More...
 
const ProcessingSoftwaresgetProcessingSoftwares () const
 Return the registered data processing software (immutable) More...
 
const ProcessingStepsgetProcessingSteps () const
 Return the registered data processing steps (immutable) More...
 
const DBSearchParamsgetDBSearchParams () const
 Return the registered database search parameters (immutable) More...
 
const DBSearchStepsgetDBSearchSteps () const
 Return the registered database search steps (immutable) More...
 
const ScoreTypesgetScoreTypes () const
 Return the registered score types (immutable) More...
 
const ObservationsgetObservations () const
 Return the registered observations (immutable) More...
 
const ParentSequencesgetParentSequences () const
 Return the registered parent sequences (immutable) More...
 
const ParentGroupSetsgetParentGroupSets () const
 Return the registered parent sequence groupings (immutable) More...
 
const IdentifiedPeptidesgetIdentifiedPeptides () const
 Return the registered identified peptides (immutable) More...
 
const IdentifiedCompoundsgetIdentifiedCompounds () const
 Return the registered compounds (immutable) More...
 
const IdentifiedOligosgetIdentifiedOligos () const
 Return the registered identified oligonucleotides (immutable) More...
 
const AdductsgetAdducts () const
 Return the registered adducts (immutable) More...
 
const ObservationMatchesgetObservationMatches () const
 Return the registered observation matches (immutable) More...
 
const ObservationMatchGroupsgetObservationMatchGroups () const
 Return the registered groups of observation matches (immutable) More...
 
void addScore (ObservationMatchRef match_ref, ScoreTypeRef score_ref, double value)
 Add a score to an input match (e.g. PSM) More...
 
void setCurrentProcessingStep (ProcessingStepRef step_ref)
 Set a data processing step that will apply to all subsequent "register..." calls. More...
 
ProcessingStepRef getCurrentProcessingStep ()
 Return the current processing step (set via setCurrentProcessingStep()). More...
 
void clearCurrentProcessingStep ()
 Cancel the effect of setCurrentProcessingStep(). More...
 
std::vector< ObservationMatchRefgetBestMatchPerObservation (ScoreTypeRef score_ref, bool require_score=false) const
 Return the best match for each observation, according to a given score type. More...
 
std::pair< ObservationMatchRef, ObservationMatchRefgetMatchesForObservation (ObservationRef obs_ref) const
 Get range of matches (cf. equal_range) for a given observation. More...
 
ScoreTypeRef findScoreType (const String &score_name) const
 Look up a score type by name. More...
 
void calculateCoverages (bool check_molecule_length=false)
 Calculate sequence coverages of parent sequences. More...
 
void cleanup (bool require_observation_match=true, bool require_identified_sequence=true, bool require_parent_match=true, bool require_parent_group=false, bool require_match_group=false)
 Clean up the data structure after filtering parts of it. More...
 
bool empty () const
 Return whether the data structure is empty (no data) More...
 
RefTranslator merge (const IdentificationData &other)
 Merge in data from another instance. More...
 
void swap (IdentificationData &other)
 Swap contents with a second instance. More...
 
void clear ()
 Clear all contents. More...
 
template<class ScoredProcessingResults >
ScoreTypeRef pickScoreType (const ScoredProcessingResults &container, bool all_elements=false, bool any_score=false) const
 
void setMetaValue (const ObservationMatchRef ref, const String &key, const DataValue &value)
 Set a meta value on a stored input match. More...
 
void setMetaValue (const ObservationRef ref, const String &key, const DataValue &value)
 Set a meta value on a stored input item. More...
 
void setMetaValue (const IdentifiedMolecule &var, const String &key, const DataValue &value)
 Set a meta value on a stored identified molecule (variant) More...
 
void setMetaValue (const String &name, const DataValue &value)
 Sets the DataValue corresponding to a name. More...
 
void setMetaValue (UInt index, const DataValue &value)
 Sets the DataValue corresponding to an index. More...
 
- Public Member Functions inherited from MetaInfoInterface
 MetaInfoInterface ()
 Constructor. More...
 
 MetaInfoInterface (const MetaInfoInterface &rhs)
 Copy constructor. More...
 
 MetaInfoInterface (MetaInfoInterface &&) noexcept
 Move constructor. More...
 
 ~MetaInfoInterface ()
 Destructor. More...
 
MetaInfoInterfaceoperator= (const MetaInfoInterface &rhs)
 Assignment operator. More...
 
MetaInfoInterfaceoperator= (MetaInfoInterface &&) noexcept
 Move assignment operator. More...
 
void swap (MetaInfoInterface &rhs)
 Swap contents. More...
 
bool operator== (const MetaInfoInterface &rhs) const
 Equality operator. More...
 
bool operator!= (const MetaInfoInterface &rhs) const
 Equality operator. More...
 
const DataValuegetMetaValue (const String &name, const DataValue &default_value=DataValue::EMPTY) const
 Returns the value corresponding to a string, or a default value (default: DataValue::EMPTY) if not found. More...
 
const DataValuegetMetaValue (UInt index, const DataValue &default_value=DataValue::EMPTY) const
 Returns the value corresponding to an index, or a default value (default: DataValue::EMPTY) if not found. More...
 
bool metaValueExists (const String &name) const
 Returns whether an entry with the given name exists. More...
 
bool metaValueExists (UInt index) const
 Returns whether an entry with the given index exists. More...
 
void setMetaValue (const String &name, const DataValue &value)
 Sets the DataValue corresponding to a name. More...
 
void setMetaValue (UInt index, const DataValue &value)
 Sets the DataValue corresponding to an index. More...
 
void removeMetaValue (const String &name)
 Removes the DataValue corresponding to name if it exists. More...
 
void removeMetaValue (UInt index)
 Removes the DataValue corresponding to index if it exists. More...
 
void addMetaValues (const MetaInfoInterface &from)
 function to copy all meta values from one object to this one More...
 
void getKeys (std::vector< String > &keys) const
 Fills the given vector with a list of all keys for which a value is set. More...
 
void getKeys (std::vector< UInt > &keys) const
 Fills the given vector with a list of all keys for which a value is set. More...
 
bool isMetaEmpty () const
 Returns if the MetaInfo is empty. More...
 
void clearMetaInfo ()
 Removes all meta values. More...
 

Protected Member Functions

void checkScoreTypes_ (const std::map< ScoreTypeRef, double > &scores) const
 Helper function to check if all score types are valid. More...
 
void checkAppliedProcessingSteps_ (const AppliedProcessingSteps &steps_and_scores) const
 Helper function to check if all applied processing steps are valid. More...
 
void checkParentMatches_ (const ParentMatches &matches, MoleculeType expected_type) const
 Helper function to check if all parent matches are valid. More...
 
void mergeScoredProcessingResults_ (ScoredProcessingResult &result, const ScoredProcessingResult &other, const RefTranslator &trans)
 Helper function to merge scored processing results while updating references (to processing steps and score types) More...
 
template<typename ContainerType , typename ElementType >
ContainerType::iterator insertIntoMultiIndex_ (ContainerType &container, const ElementType &element)
 Helper function for adding entries (derived from ScoredProcessingResult) to a boost::multi_index_container structure. More...
 
template<typename ContainerType , typename ElementType >
ContainerType::iterator insertIntoMultiIndex_ (ContainerType &container, const ElementType &element, AddressLookup &lookup)
 Variant of insertIntoMultiIndex_() that also updates a look-up table of valid references (addresses) More...
 
template<typename RefType , typename ContainerType >
void setMetaValue_ (const RefType ref, const String &key, const DataValue &value, ContainerType &container, const AddressLookup &lookup=AddressLookup())
 Helper function to add a meta value to an element in a multi-index container. More...
 
- Protected Member Functions inherited from MetaInfoInterface
void createIfNotExists_ ()
 Creates the MetaInfo object if it does not exist. More...
 

Static Protected Member Functions

template<typename RefType , typename ContainerType >
static bool isValidReference_ (RefType ref, ContainerType &container)
 Check whether a reference points to an element in a container. More...
 
template<typename RefType >
static bool isValidHashedReference_ (RefType ref, const AddressLookup &lookup)
 Check validity of a reference based on a look-up table of addresses. More...
 
template<typename ContainerType , typename PredicateType >
static void removeFromSetIf_ (ContainerType &container, PredicateType predicate)
 Remove elements from a set (or ordered multi_index_container) if they fulfill a predicate. More...
 
template<typename ContainerType >
static void removeFromSetIfNotHashed_ (ContainerType &container, const AddressLookup &lookup)
 Remove elements from a set (or ordered multi_index_container) if they don't occur in a look-up table. More...
 
template<typename ContainerType >
static void updateAddressLookup_ (const ContainerType &container, AddressLookup &lookup)
 Recreate the address look-up table for a container. More...
 

Protected Attributes

InputFiles input_files_
 
ProcessingSoftwares processing_softwares_
 
ProcessingSteps processing_steps_
 
DBSearchParams db_search_params_
 
DBSearchSteps db_search_steps_
 
ScoreTypes score_types_
 
Observations observations_
 
ParentSequences parents_
 
ParentGroupSets parent_groups_
 
IdentifiedPeptides identified_peptides_
 
IdentifiedCompounds identified_compounds_
 
IdentifiedOligos identified_oligos_
 
Adducts adducts_
 
ObservationMatches observation_matches_
 
ObservationMatchGroups observation_match_groups_
 
ProcessingStepRef current_step_ref_
 Reference to the current data processing step (see setCurrentProcessingStep()) More...
 
bool no_checks_
 Suppress validity checks in register... calls? More...
 
AddressLookup observation_lookup_
 
AddressLookup parent_lookup_
 
AddressLookup identified_peptide_lookup_
 
AddressLookup identified_compound_lookup_
 
AddressLookup identified_oligo_lookup_
 
AddressLookup observation_match_lookup_
 
- Protected Attributes inherited from MetaInfoInterface
MetaInfometa_
 Pointer to the MetaInfo object (0 by default) More...
 

Friends

class IDFilter
 
class MapAlignmentTransformer
 

Additional Inherited Members

- Static Public Member Functions inherited from MetaInfoInterface
static MetaInfoRegistrymetaRegistry ()
 Returns a reference to the MetaInfoRegistry. More...
 

Detailed Description

Representation of spectrum identification results and associated data.

This class provides capabilities for storing spectrum identification results from different types of experiments/molecules (proteomics: peptides/proteins, metabolomics: small molecules, "nucleomics": RNA).

The class design has the following goals:

The following important subordinate classes are provided to represent different types of data:

Class Represents Key Proteomics example Corresponding legacy class
ProcessingStep Information about a data processing step that was applied (e.g. input files, software used, parameters) Combined information Mascot search ProteinIdentification
Observation A search query (with identifier, RT, m/z) from an input file, i.e. an MS2 spectrum or feature (for accurate mass search) File/Identifier MS2 spectrum PeptideIdentification
ParentSequence An entry in a FASTA file with associated information (sequence, coverage, etc.) Accession Protein ProteinHit
IdentifiedPeptide/-Oligo/-Compound An identified molecule of the respective type Sequence (or identifier for a compound) Peptide PeptideHit
ObservationMatch A match between a query (Observation), identified molecule (Identified...), and optionally adduct Combination of query/molecule/adduct references Peptide-spectrum match (PSM) PeptideIdentification/PeptideHit

To populate an IdentificationData instance with data, "register..." functions are used. These functions return "references" (implemented as iterators) that can be used to refer to stored data items and thus form connections. For example, a protein can be stored using registerParentSequence, which returns a corresponding reference. This reference can be used to build an IdentifiedPeptide object that references the protein. An identified peptide referencing a protein can only be registered if that protein has been registered already, to ensure data consistency. Given the identified peptide, information about the associated protein can be retrieved efficiently by simply dereferencing the reference.

To ensure non-redundancy, many data types have a "key" (see table above) to which a uniqueness constraint applies. This means only one item of such a type with a given key can be stored in an IdentificationData object. If items with an existing key are registered subsequently, attempts are made to merge new information (e.g. additional scores) into the existing entry. The details of this merging are handled in the merge function in each data class.

Warning
This class is not thread-safe while being modified.

Member Typedef Documentation

◆ AddressLookup

using AddressLookup = std::unordered_set<uintptr_t>

◆ AdductOpt

◆ AdductRef

◆ Adducts

◆ AppliedProcessingStep

◆ AppliedProcessingSteps

◆ DBSearchParam

◆ DBSearchParams

◆ DBSearchSteps

◆ IdentifiedCompound

◆ IdentifiedCompoundRef

◆ IdentifiedCompounds

◆ IdentifiedMolecule

◆ IdentifiedOligo

◆ IdentifiedOligoRef

◆ IdentifiedOligos

◆ IdentifiedPeptide

◆ IdentifiedPeptideRef

◆ IdentifiedPeptides

◆ InputFile

◆ InputFileRef

◆ InputFiles

◆ MassType

◆ MatchGroupRef

◆ MoleculeType

◆ Observation

◆ ObservationMatch

◆ ObservationMatches

◆ ObservationMatchGroup

◆ ObservationMatchGroups

◆ ObservationMatchRef

◆ ObservationRef

◆ Observations

◆ ParentGroup

◆ ParentGroupRef

◆ ParentGroups

◆ ParentGroupSet

◆ ParentGroupSets

◆ ParentMatch

◆ ParentMatches

◆ ParentSequence

◆ ParentSequenceRef

◆ ParentSequences

◆ PeakAnnotations

◆ ProcessingSoftware

◆ ProcessingSoftwareRef

◆ ProcessingSoftwares

◆ ProcessingStep

◆ ProcessingStepRef

◆ ProcessingSteps

◆ ScoredProcessingResult

◆ ScoreType

◆ ScoreTypeRef

◆ ScoreTypes

◆ SearchParamRef

Constructor & Destructor Documentation

◆ IdentificationData() [1/3]

IdentificationData ( )
inline

Default constructor.

◆ IdentificationData() [2/3]

Copy constructor.

Copy-constructing is expensive due to the necessary "rewiring" of references. Use the move constructor where possible.

◆ IdentificationData() [3/3]

IdentificationData ( IdentificationData &&  other)
inlinenoexcept

Move constructor.

Member Function Documentation

◆ addScore()

void addScore ( ObservationMatchRef  match_ref,
ScoreTypeRef  score_ref,
double  value 
)

Add a score to an input match (e.g. PSM)

◆ calculateCoverages()

void calculateCoverages ( bool  check_molecule_length = false)

Calculate sequence coverages of parent sequences.

Referenced by NucleicAcidSearchEngine::main_().

◆ checkAppliedProcessingSteps_()

void checkAppliedProcessingSteps_ ( const AppliedProcessingSteps steps_and_scores) const
protected

Helper function to check if all applied processing steps are valid.

◆ checkParentMatches_()

void checkParentMatches_ ( const ParentMatches matches,
MoleculeType  expected_type 
) const
protected

Helper function to check if all parent matches are valid.

◆ checkScoreTypes_()

void checkScoreTypes_ ( const std::map< ScoreTypeRef, double > &  scores) const
protected

Helper function to check if all score types are valid.

◆ cleanup()

void cleanup ( bool  require_observation_match = true,
bool  require_identified_sequence = true,
bool  require_parent_match = true,
bool  require_parent_group = false,
bool  require_match_group = false 
)

Clean up the data structure after filtering parts of it.

Make sure there are no invalid references or "orphan" data entries.

Parameters
require_observation_matchRemove identified molecules, observations and adducts that aren't part of observation matches?
require_identified_sequenceRemove parent sequences (proteins/RNAs) that aren't referenced by identified peptides/oligonucleotides?
require_parent_matchRemove identified peptides/oligonucleotides that don't reference a parent sequence (protein/RNA)?
require_parent_groupRemove parent sequences that aren't part of parent sequence groups?
require_match_groupRemove input matches that aren't part of match groups?

Referenced by IDFilter::filterObservationMatchesByFunctor(), and NucleicAcidSearchEngine::postProcessHits_().

◆ clear()

void clear ( )

Clear all contents.

◆ clearCurrentProcessingStep()

void clearCurrentProcessingStep ( )

Cancel the effect of setCurrentProcessingStep().

◆ empty()

bool empty ( ) const

Return whether the data structure is empty (no data)

◆ findScoreType()

ScoreTypeRef findScoreType ( const String score_name) const

Look up a score type by name.

Returns
Reference to the score type, if found; otherwise getScoreTypes().end()

Referenced by NucleicAcidSearchEngine::calculateAndFilterFDR_().

◆ getAdducts()

const Adducts& getAdducts ( ) const
inline

Return the registered adducts (immutable)

◆ getBestMatchPerObservation()

std::vector<ObservationMatchRef> getBestMatchPerObservation ( ScoreTypeRef  score_ref,
bool  require_score = false 
) const

Return the best match for each observation, according to a given score type.

Parameters
score_refScore type to use
require_scoreExclude matches without score of this type, even if they are the only matches for their observations?

◆ getCurrentProcessingStep()

ProcessingStepRef getCurrentProcessingStep ( )

Return the current processing step (set via setCurrentProcessingStep()).

If no current processing step has been set, processing_steps.end() is returned.

Referenced by NucleicAcidSearchEngine::main_(), and NucleicAcidSearchEngine::postProcessHits_().

◆ getDBSearchParams()

const DBSearchParams& getDBSearchParams ( ) const
inline

Return the registered database search parameters (immutable)

Referenced by NucleicAcidSearchEngine::main_().

◆ getDBSearchSteps()

const DBSearchSteps& getDBSearchSteps ( ) const
inline

Return the registered database search steps (immutable)

◆ getIdentifiedCompounds()

const IdentifiedCompounds& getIdentifiedCompounds ( ) const
inline

Return the registered compounds (immutable)

◆ getIdentifiedOligos()

const IdentifiedOligos& getIdentifiedOligos ( ) const
inline

Return the registered identified oligonucleotides (immutable)

Referenced by NucleicAcidSearchEngine::main_().

◆ getIdentifiedPeptides()

const IdentifiedPeptides& getIdentifiedPeptides ( ) const
inline

Return the registered identified peptides (immutable)

◆ getInputFiles()

const InputFiles& getInputFiles ( ) const
inline

Return the registered input files (immutable)

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ getMatchesForObservation()

std::pair<ObservationMatchRef, ObservationMatchRef> getMatchesForObservation ( ObservationRef  obs_ref) const

Get range of matches (cf. equal_range) for a given observation.

◆ getObservationMatches()

const ObservationMatches& getObservationMatches ( ) const
inline

Return the registered observation matches (immutable)

Referenced by NucleicAcidSearchEngine::calculateAndFilterFDR_(), and NucleicAcidSearchEngine::generateLFQInput_().

◆ getObservationMatchGroups()

const ObservationMatchGroups& getObservationMatchGroups ( ) const
inline

Return the registered groups of observation matches (immutable)

◆ getObservations()

const Observations& getObservations ( ) const
inline

Return the registered observations (immutable)

Referenced by NucleicAcidSearchEngine::calculateAndFilterFDR_(), and NucleicAcidSearchEngine::main_().

◆ getParentGroupSets()

const ParentGroupSets& getParentGroupSets ( ) const
inline

Return the registered parent sequence groupings (immutable)

◆ getParentSequences()

const ParentSequences& getParentSequences ( ) const
inline

Return the registered parent sequences (immutable)

Referenced by NucleicAcidSearchEngine::main_().

◆ getProcessingSoftwares()

const ProcessingSoftwares& getProcessingSoftwares ( ) const
inline

Return the registered data processing software (immutable)

◆ getProcessingSteps()

const ProcessingSteps& getProcessingSteps ( ) const
inline

Return the registered data processing steps (immutable)

◆ getScoreTypes()

const ScoreTypes& getScoreTypes ( ) const
inline

Return the registered score types (immutable)

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ insertIntoMultiIndex_() [1/2]

ContainerType::iterator insertIntoMultiIndex_ ( ContainerType &  container,
const ElementType &  element 
)
inlineprotected

Helper function for adding entries (derived from ScoredProcessingResult) to a boost::multi_index_container structure.

◆ insertIntoMultiIndex_() [2/2]

ContainerType::iterator insertIntoMultiIndex_ ( ContainerType &  container,
const ElementType &  element,
AddressLookup lookup 
)
inlineprotected

Variant of insertIntoMultiIndex_() that also updates a look-up table of valid references (addresses)

◆ isValidHashedReference_()

static bool isValidHashedReference_ ( RefType  ref,
const AddressLookup lookup 
)
inlinestaticprotected

Check validity of a reference based on a look-up table of addresses.

◆ isValidReference_()

static bool isValidReference_ ( RefType  ref,
ContainerType &  container 
)
inlinestaticprotected

Check whether a reference points to an element in a container.

◆ merge()

RefTranslator merge ( const IdentificationData other)

Merge in data from another instance.

Can be used to make a deep copy by calling merge() on an empty object. The returned translation table allows updating of references that are held externally.

Parameters
otherInstance to merge in.
Returns
Translation table for references (old -> new)

◆ mergeScoredProcessingResults_()

void mergeScoredProcessingResults_ ( ScoredProcessingResult result,
const ScoredProcessingResult other,
const RefTranslator trans 
)
protected

Helper function to merge scored processing results while updating references (to processing steps and score types)

Parameters
resultInstance that gets updated
otherInstance to merge into result
transMapping of corresponding references between other and result

◆ pickScoreType()

ScoreTypeRef pickScoreType ( const ScoredProcessingResults &  container,
bool  all_elements = false,
bool  any_score = false 
) const
inline

Pick a score type for operations (e.g. filtering) on a container of scored processing results (e.g. input matches, identified peptides, ...).

If all_elements is false, only the first element with a score will be considered (which is sufficient if all elements were processed in the same way). If all_elements is true, the score type supported by the highest number of elements will be chosen.

If any_score is false, only the primary score from the most recent processing step (that assigned a score) is taken into account. If any_score is true, all score types assigned across all elements are considered (this implies all_elements = true).

Parameters
containerContainer with elements derived from ScoredProcessingResult
all_elementsConsider all elements?
any_scoreConsider any score (or just primary/most recent ones)?
Returns
Reference to the chosen score type (or getScoreTypes().end() if there were no scores)

◆ registerAdduct()

AdductRef registerAdduct ( const AdductInfo adduct)

Register an adduct.

Returns
Reference to the registered adduct

Referenced by NucleicAcidSearchEngine::main_().

◆ registerDBSearchParam()

SearchParamRef registerDBSearchParam ( const DBSearchParam param)

Register database search parameters.

Returns
Reference to the registered search parameters

Referenced by NucleicAcidSearchEngine::main_().

◆ registerIdentifiedCompound()

IdentifiedCompoundRef registerIdentifiedCompound ( const IdentifiedCompound compound)

Register an identified compound (small molecule)

Returns
Reference to the registered compound

◆ registerIdentifiedOligo()

IdentifiedOligoRef registerIdentifiedOligo ( const IdentifiedOligo oligo)

Register an identified RNA oligonucleotide.

Returns
Reference to the registered oligonucleotide

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ registerIdentifiedPeptide()

IdentifiedPeptideRef registerIdentifiedPeptide ( const IdentifiedPeptide peptide)

Register an identified peptide.

Returns
Reference to the registered peptide

◆ registerInputFile()

InputFileRef registerInputFile ( const InputFile file)

Register an input file.

Returns
Reference to the registered file

Referenced by NucleicAcidSearchEngine::main_().

◆ registerObservation()

ObservationRef registerObservation ( const Observation obs)

Register an observation (e.g. MS2 spectrum or feature)

Returns
Reference to the registered observation

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ registerObservationMatch()

ObservationMatchRef registerObservationMatch ( const ObservationMatch match)

Register an observation match (e.g. peptide-spectrum match)

Returns
Reference to the registered observation match

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ registerObservationMatchGroup()

MatchGroupRef registerObservationMatchGroup ( const ObservationMatchGroup group)

Register a group of observation matches that belong together.

Returns
Reference to the registered group of observation matches

◆ registerParentGroupSet()

void registerParentGroupSet ( const ParentGroupSet groups)

Register a grouping of parent sequences (e.g. protein inference result)

◆ registerParentSequence()

ParentSequenceRef registerParentSequence ( const ParentSequence parent)

Register a parent sequence (e.g. protein or intact RNA)

Returns
Reference to the registered parent sequence

Referenced by NucleicAcidSearchEngine::main_().

◆ registerProcessingSoftware()

ProcessingSoftwareRef registerProcessingSoftware ( const ProcessingSoftware software)

Register data processing software.

Returns
Reference to the registered software

Referenced by NucleicAcidSearchEngine::main_().

◆ registerProcessingStep() [1/2]

ProcessingStepRef registerProcessingStep ( const ProcessingStep step)

Register a data processing step.

Returns
Reference to the registered processing step

Referenced by NucleicAcidSearchEngine::main_().

◆ registerProcessingStep() [2/2]

ProcessingStepRef registerProcessingStep ( const ProcessingStep step,
SearchParamRef  search_ref 
)

Register a database search step with associated parameters.

Returns
Reference to the registered processing step

◆ registerScoreType()

ScoreTypeRef registerScoreType ( const ScoreType score)

Register a score type.

Returns
Reference to the registered score type

Referenced by NucleicAcidSearchEngine::main_().

◆ removeFromSetIf_()

static void removeFromSetIf_ ( ContainerType &  container,
PredicateType  predicate 
)
inlinestaticprotected

Remove elements from a set (or ordered multi_index_container) if they fulfill a predicate.

Referenced by IDFilter::filterObservationMatchesByFunctor().

◆ removeFromSetIfNotHashed_()

static void removeFromSetIfNotHashed_ ( ContainerType &  container,
const AddressLookup lookup 
)
inlinestaticprotected

Remove elements from a set (or ordered multi_index_container) if they don't occur in a look-up table.

◆ setCurrentProcessingStep()

void setCurrentProcessingStep ( ProcessingStepRef  step_ref)

Set a data processing step that will apply to all subsequent "register..." calls.

This step will be appended to the list of processing steps for all relevant elements that are registered subsequently (unless it is already the last entry in the list). If a score type without a software reference is registered, the software reference of this processing step will be applied. Effective until clearCurrentProcessingStep() is called.

Referenced by NucleicAcidSearchEngine::main_().

◆ setMetaValue() [1/5]

void setMetaValue ( const IdentifiedMolecule var,
const String key,
const DataValue value 
)

Set a meta value on a stored identified molecule (variant)

◆ setMetaValue() [2/5]

void setMetaValue ( const ObservationMatchRef  ref,
const String key,
const DataValue value 
)

Set a meta value on a stored input match.

◆ setMetaValue() [3/5]

void setMetaValue ( const ObservationRef  ref,
const String key,
const DataValue value 
)

Set a meta value on a stored input item.

◆ setMetaValue() [4/5]

void setMetaValue

Sets the DataValue corresponding to a name.

◆ setMetaValue() [5/5]

void setMetaValue

Sets the DataValue corresponding to an index.

◆ setMetaValue_()

void setMetaValue_ ( const RefType  ref,
const String key,
const DataValue value,
ContainerType &  container,
const AddressLookup lookup = AddressLookup() 
)
inlineprotected

Helper function to add a meta value to an element in a multi-index container.

◆ swap()

void swap ( IdentificationData other)

Swap contents with a second instance.

◆ updateAddressLookup_()

static void updateAddressLookup_ ( const ContainerType &  container,
AddressLookup lookup 
)
inlinestaticprotected

Recreate the address look-up table for a container.

Friends And Related Function Documentation

◆ IDFilter

friend class IDFilter
friend

◆ MapAlignmentTransformer

friend class MapAlignmentTransformer
friend

Member Data Documentation

◆ adducts_

Adducts adducts_
protected

◆ current_step_ref_

ProcessingStepRef current_step_ref_
protected

Reference to the current data processing step (see setCurrentProcessingStep())

◆ db_search_params_

DBSearchParams db_search_params_
protected

◆ db_search_steps_

DBSearchSteps db_search_steps_
protected

◆ identified_compound_lookup_

AddressLookup identified_compound_lookup_
protected

◆ identified_compounds_

IdentifiedCompounds identified_compounds_
protected

◆ identified_oligo_lookup_

AddressLookup identified_oligo_lookup_
protected

◆ identified_oligos_

IdentifiedOligos identified_oligos_
protected

◆ identified_peptide_lookup_

AddressLookup identified_peptide_lookup_
protected

◆ identified_peptides_

IdentifiedPeptides identified_peptides_
protected

◆ input_files_

InputFiles input_files_
protected

◆ no_checks_

bool no_checks_
protected

Suppress validity checks in register... calls?

This is useful in situations where validity is already guaranteed (e.g. copying).

◆ observation_lookup_

AddressLookup observation_lookup_
protected

◆ observation_match_groups_

ObservationMatchGroups observation_match_groups_
protected

◆ observation_match_lookup_

AddressLookup observation_match_lookup_
protected

◆ observation_matches_

ObservationMatches observation_matches_
protected

◆ observations_

Observations observations_
protected

◆ parent_groups_

ParentGroupSets parent_groups_
protected

◆ parent_lookup_

AddressLookup parent_lookup_
protected

◆ parents_

ParentSequences parents_
protected

◆ processing_softwares_

ProcessingSoftwares processing_softwares_
protected

◆ processing_steps_

ProcessingSteps processing_steps_
protected

◆ score_types_

ScoreTypes score_types_
protected