OpenMS  2.7.0
Classes | Public Types | Public Member Functions | Static 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 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...
 

Public Types

using MoleculeType = IdentificationDataInternal::MoleculeType
 
using MassType = IdentificationDataInternal::MassType
 
using InputFiles = IdentificationDataInternal::InputFiles
 
using InputFileRef = IdentificationDataInternal::InputFileRef
 
using DataProcessingSoftware = IdentificationDataInternal::DataProcessingSoftware
 
using DataProcessingSoftwares = IdentificationDataInternal::DataProcessingSoftwares
 
using ProcessingSoftwareRef = IdentificationDataInternal::ProcessingSoftwareRef
 
using DataProcessingStep = IdentificationDataInternal::DataProcessingStep
 
using DataProcessingSteps = IdentificationDataInternal::DataProcessingSteps
 
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 AppliedProcessingStep = IdentificationDataInternal::AppliedProcessingStep
 
using AppliedProcessingSteps = IdentificationDataInternal::AppliedProcessingSteps
 
using DataQuery = IdentificationDataInternal::DataQuery
 
using DataQueries = IdentificationDataInternal::DataQueries
 
using DataQueryRef = IdentificationDataInternal::DataQueryRef
 
using ParentMolecule = IdentificationDataInternal::ParentMolecule
 
using ParentMolecules = IdentificationDataInternal::ParentMolecules
 
using ParentMoleculeRef = IdentificationDataInternal::ParentMoleculeRef
 
using MoleculeParentMatch = IdentificationDataInternal::MoleculeParentMatch
 
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 PeakAnnotations = IdentificationDataInternal::PeakAnnotations
 
using IdentifiedMoleculeRef = IdentificationDataInternal::IdentifiedMoleculeRef
 
using MoleculeQueryMatch = IdentificationDataInternal::MoleculeQueryMatch
 
using MoleculeQueryMatches = IdentificationDataInternal::MoleculeQueryMatches
 
using QueryMatchRef = IdentificationDataInternal::QueryMatchRef
 
using QueryMatchGroup = IdentificationDataInternal::QueryMatchGroup
 
using QueryMatchGroups = IdentificationDataInternal::QueryMatchGroups
 
using MatchGroupRef = IdentificationDataInternal::MatchGroupRef
 
using ParentMoleculeGroup = IdentificationDataInternal::ParentMoleculeGroup
 
using ParentMoleculeGroups = IdentificationDataInternal::ParentMoleculeGroups
 
using ParentGroupRef = IdentificationDataInternal::ParentGroupRef
 
using ParentMoleculeGrouping = IdentificationDataInternal::ParentMoleculeGrouping
 
using ParentMoleculeGroupings = IdentificationDataInternal::ParentMoleculeGroupings
 
using AddressLookup = boost::unordered_set< uintptr_t >
 

Public Member Functions

 IdentificationData ()
 Default constructor. More...
 
 IdentificationData (const IdentificationData &other)=delete
 
 IdentificationData (IdentificationData &&other)
 Move constructor. More...
 
InputFileRef registerInputFile (const String &file)
 Register an input file. More...
 
ProcessingSoftwareRef registerDataProcessingSoftware (const DataProcessingSoftware &software)
 Register data processing software. More...
 
SearchParamRef registerDBSearchParam (const DBSearchParam &param)
 Register database search parameters. More...
 
ProcessingStepRef registerDataProcessingStep (const DataProcessingStep &step)
 Register a data processing step. More...
 
ProcessingStepRef registerDataProcessingStep (const DataProcessingStep &step, SearchParamRef search_ref)
 Register a database search step with associated parameters. More...
 
ScoreTypeRef registerScoreType (const ScoreType &score)
 Register a score type. More...
 
DataQueryRef registerDataQuery (const DataQuery &query)
 Register a data query (e.g. MS2 spectrum or feature) More...
 
ParentMoleculeRef registerParentMolecule (const ParentMolecule &parent)
 Register a parent molecule (e.g. protein or intact RNA) More...
 
void registerParentMoleculeGrouping (const ParentMoleculeGrouping &grouping)
 Register a grouping of parent molecules (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...
 
QueryMatchRef registerMoleculeQueryMatch (const MoleculeQueryMatch &match)
 Register a molecule-query match (e.g. peptide-spectrum match) More...
 
MatchGroupRef registerQueryMatchGroup (const QueryMatchGroup &group)
 Register a group of associated molecule-query matches. More...
 
const InputFilesgetInputFiles () const
 Return the registered input files (immutable) More...
 
const DataProcessingSoftwaresgetDataProcessingSoftwares () const
 Return the registered data processing software (immutable) More...
 
const DataProcessingStepsgetDataProcessingSteps () 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 DataQueriesgetDataQueries () const
 Return the registered data queries (immutable) More...
 
const ParentMoleculesgetParentMolecules () const
 Return the registered parent molecules (immutable) More...
 
const ParentMoleculeGroupingsgetParentMoleculeGroupings () const
 Return the registered parent molecule 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 MoleculeQueryMatchesgetMoleculeQueryMatches () const
 Return the registered molecule-query matches (immutable) More...
 
const QueryMatchGroupsgetQueryMatchGroups () const
 Return the registered groups of molecule-query matches (immutable) More...
 
void addScore (QueryMatchRef match_ref, ScoreTypeRef score_ref, double value)
 Add a score to a molecule-query 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 ()
 
void clearCurrentProcessingStep ()
 Cancel the effect of setCurrentProcessingStep(). More...
 
std::vector< QueryMatchRefgetBestMatchPerQuery (ScoreTypeRef score_ref) const
 Return the best match for each data query, according to a given score type. More...
 
std::pair< ScoreTypeRef, bool > findScoreType (const String &score_name) const
 Look up a score type by name. More...
 
void calculateCoverages (bool check_molecule_length=false)
 Calculate sequence coverage of parent molecules. More...
 
void cleanup (bool require_query_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...
 
- 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...
 
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...
 

Static Public Member Functions

static bool isBetterScore (double first, double second, bool higher_better)
 Helper function to compare two scores. More...
 
- Static Public Member Functions inherited from MetaInfoInterface
static MetaInfoRegistrymetaRegistry ()
 Returns a reference to the MetaInfoRegistry. 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...
 
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...
 
- 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_
 
DataProcessingSoftwares processing_softwares_
 
DataProcessingSteps processing_steps_
 
DBSearchParams db_search_params_
 
DBSearchSteps db_search_steps_
 
ScoreTypes score_types_
 
DataQueries data_queries_
 
ParentMolecules parent_molecules_
 
ParentMoleculeGroupings parent_molecule_groupings_
 
IdentifiedPeptides identified_peptides_
 
IdentifiedCompounds identified_compounds_
 
IdentifiedOligos identified_oligos_
 
MoleculeQueryMatches query_matches_
 
QueryMatchGroups query_match_groups_
 
ProcessingStepRef current_step_ref_
 Reference to the current data processing step (see setCurrentProcessingStep()) More...
 
AddressLookup data_query_lookup_
 
AddressLookup parent_molecule_lookup_
 
AddressLookup identified_peptide_lookup_
 
AddressLookup identified_compound_lookup_
 
AddressLookup identified_oligo_lookup_
 
AddressLookup query_match_lookup_
 
- Protected Attributes inherited from MetaInfoInterface
MetaInfometa_
 Pointer to the MetaInfo object (0 by default) More...
 

Friends

class IDFilter
 

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
DataProcessingStep Information about a data processing step that was applied (e.g. input files, software used, parameters) Combined information Mascot search ProteinIdentification
DataQuery A search query (with identifier, RT, m/z), i.e. an MS2 spectrum or feature (for accurate mass search) Identifier MS2 spectrum PeptideIdentification
ParentMolecule 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
MoleculeQueryMatch A match between a query (DataQuery) and identified molecule (Identified...) Combination of query and molecule 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 registerParentMolecule, 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.

Member Typedef Documentation

◆ AddressLookup

using AddressLookup = boost::unordered_set<uintptr_t>

◆ AppliedProcessingStep

◆ AppliedProcessingSteps

◆ DataProcessingSoftware

◆ DataProcessingSoftwares

◆ DataProcessingStep

◆ DataProcessingSteps

◆ DataQueries

◆ DataQuery

◆ DataQueryRef

◆ DBSearchParam

◆ DBSearchParams

◆ DBSearchSteps

◆ IdentifiedCompound

◆ IdentifiedCompoundRef

◆ IdentifiedCompounds

◆ IdentifiedMoleculeRef

◆ IdentifiedOligo

◆ IdentifiedOligoRef

◆ IdentifiedOligos

◆ IdentifiedPeptide

◆ IdentifiedPeptideRef

◆ IdentifiedPeptides

◆ InputFileRef

◆ InputFiles

◆ MassType

◆ MatchGroupRef

◆ MoleculeParentMatch

◆ MoleculeQueryMatch

◆ MoleculeQueryMatches

◆ MoleculeType

◆ ParentGroupRef

◆ ParentMatches

◆ ParentMolecule

◆ ParentMoleculeGroup

◆ ParentMoleculeGrouping

◆ ParentMoleculeGroupings

◆ ParentMoleculeGroups

◆ ParentMoleculeRef

◆ ParentMolecules

◆ PeakAnnotations

◆ ProcessingSoftwareRef

◆ ProcessingStepRef

◆ QueryMatchGroup

◆ QueryMatchGroups

◆ QueryMatchRef

◆ ScoreType

◆ ScoreTypeRef

◆ ScoreTypes

◆ SearchParamRef

Constructor & Destructor Documentation

◆ IdentificationData() [1/3]

IdentificationData ( )
inline

Default constructor.

◆ IdentificationData() [2/3]

IdentificationData ( const IdentificationData other)
delete

◆ IdentificationData() [3/3]

IdentificationData ( IdentificationData &&  other)
inline

Move constructor.

Member Function Documentation

◆ addScore()

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

Add a score to a molecule-query match (e.g. PSM)

◆ calculateCoverages()

void calculateCoverages ( bool  check_molecule_length = false)

Calculate sequence coverage of parent molecules.

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_query_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_query_matchRemove identified molecules and data queries that aren't part of molecule-query matches?
require_identified_sequenceRemove parent molecules (proteins/RNAs) that aren't referenced by identified peptides/oligonucleotides?
require_parent_matchRemove identified peptides/oligonucleotides that don't reference a parent molecule (protein/RNA)?
require_parent_groupRemove parent molecules that aren't part of parent molecule groups?
require_match_groupRemove molecule-query matches that aren't part of match groups?

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ clearCurrentProcessingStep()

void clearCurrentProcessingStep ( )

Cancel the effect of setCurrentProcessingStep().

◆ findScoreType()

std::pair<ScoreTypeRef, bool> findScoreType ( const String score_name) const

Look up a score type by name.

Returns
A pair: 1. Reference to the score type, if found; 2. Boolean indicating success or failure

Referenced by NucleicAcidSearchEngine::calculateAndFilterFDR_().

◆ getBestMatchPerQuery()

std::vector<QueryMatchRef> getBestMatchPerQuery ( ScoreTypeRef  score_ref) const

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

◆ 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::postProcessHits_().

◆ getDataProcessingSoftwares()

const DataProcessingSoftwares& getDataProcessingSoftwares ( ) const
inline

Return the registered data processing software (immutable)

◆ getDataProcessingSteps()

const DataProcessingSteps& getDataProcessingSteps ( ) const
inline

Return the registered data processing steps (immutable)

◆ getDataQueries()

const DataQueries& getDataQueries ( ) const
inline

Return the registered data queries (immutable)

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

◆ getDBSearchParams()

const DBSearchParams& getDBSearchParams ( ) const
inline

Return the registered database search parameters (immutable)

◆ 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_().

◆ getMoleculeQueryMatches()

const MoleculeQueryMatches& getMoleculeQueryMatches ( ) const
inline

Return the registered molecule-query matches (immutable)

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

◆ getParentMoleculeGroupings()

const ParentMoleculeGroupings& getParentMoleculeGroupings ( ) const
inline

Return the registered parent molecule groupings (immutable)

◆ getParentMolecules()

const ParentMolecules& getParentMolecules ( ) const
inline

Return the registered parent molecules (immutable)

◆ getQueryMatchGroups()

const QueryMatchGroups& getQueryMatchGroups ( ) const
inline

Return the registered groups of molecule-query matches (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)

◆ isBetterScore()

static bool isBetterScore ( double  first,
double  second,
bool  higher_better 
)
inlinestatic

Helper function to compare two scores.

◆ 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.

◆ registerDataProcessingSoftware()

ProcessingSoftwareRef registerDataProcessingSoftware ( const DataProcessingSoftware software)

Register data processing software.

Returns
Reference to the registered software

Referenced by NucleicAcidSearchEngine::registerIDMetaData_().

◆ registerDataProcessingStep() [1/2]

ProcessingStepRef registerDataProcessingStep ( const DataProcessingStep step)

Register a data processing step.

Returns
Reference to the registered processing step

Referenced by NucleicAcidSearchEngine::registerIDMetaData_().

◆ registerDataProcessingStep() [2/2]

ProcessingStepRef registerDataProcessingStep ( const DataProcessingStep step,
SearchParamRef  search_ref 
)

Register a database search step with associated parameters.

Returns
Reference to the registered processing step

◆ registerDataQuery()

DataQueryRef registerDataQuery ( const DataQuery query)

Register a data query (e.g. MS2 spectrum or feature)

Returns
Reference to the registered data query

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ registerDBSearchParam()

SearchParamRef registerDBSearchParam ( const DBSearchParam param)

Register database search parameters.

Returns
Reference to the registered search parameters

Referenced by NucleicAcidSearchEngine::registerIDMetaData_().

◆ 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 String file)

Register an input file.

Returns
Reference to the registered file

Referenced by NucleicAcidSearchEngine::registerIDMetaData_().

◆ registerMoleculeQueryMatch()

QueryMatchRef registerMoleculeQueryMatch ( const MoleculeQueryMatch match)

Register a molecule-query match (e.g. peptide-spectrum match)

Returns
Reference to the registered molecule-query match

Referenced by NucleicAcidSearchEngine::postProcessHits_().

◆ registerParentMolecule()

ParentMoleculeRef registerParentMolecule ( const ParentMolecule parent)

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

Returns
Reference to the registered parent molecule

◆ registerParentMoleculeGrouping()

void registerParentMoleculeGrouping ( const ParentMoleculeGrouping grouping)

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

◆ registerQueryMatchGroup()

MatchGroupRef registerQueryMatchGroup ( const QueryMatchGroup group)

Register a group of associated molecule-query matches.

Returns
Reference to the registered group of matches

◆ registerScoreType()

ScoreTypeRef registerScoreType ( const ScoreType score)

Register a score type.

Returns
Reference to the registered score type

Referenced by NucleicAcidSearchEngine::registerIDMetaData_().

◆ removeFromSetIf_()

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

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

◆ 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::registerIDMetaData_().

◆ 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

Member Data Documentation

◆ current_step_ref_

ProcessingStepRef current_step_ref_
protected

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

◆ data_queries_

DataQueries data_queries_
protected

◆ data_query_lookup_

AddressLookup data_query_lookup_
protected

◆ 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

◆ parent_molecule_groupings_

ParentMoleculeGroupings parent_molecule_groupings_
protected

◆ parent_molecule_lookup_

AddressLookup parent_molecule_lookup_
protected

◆ parent_molecules_

ParentMolecules parent_molecules_
protected

◆ processing_softwares_

DataProcessingSoftwares processing_softwares_
protected

◆ processing_steps_

DataProcessingSteps processing_steps_
protected

◆ query_match_groups_

QueryMatchGroups query_match_groups_
protected

◆ query_match_lookup_

AddressLookup query_match_lookup_
protected

◆ query_matches_

MoleculeQueryMatches query_matches_
protected

◆ score_types_

ScoreTypes score_types_
protected