|
OpenMS
2.5.0
|
Go to the documentation of this file.
147 inline ExitCodes run(std::vector<FASTAFile::FASTAEntry>& proteins, std::vector<ProteinIdentification>& prot_ids, std::vector<PeptideIdentification>& pep_ids)
150 return run<TFI_Vector>(protein_container, prot_ids, pep_ids);
192 if (decoy_string_.empty())
200 OPENMS_LOG_WARN <<
"Unable to determine decoy string automatically (not enough decoys were detected)! Using default " << (r.is_prefix ?
"prefix" :
"suffix") <<
" decoy string '" << r.name <<
"'\n"
201 <<
"If you think that this is incorrect, please provide a decoy_string and its position manually!" << std::endl;
203 prefix_ = r.is_prefix;
204 decoy_string_ = r.name;
206 OPENMS_LOG_INFO <<
"Using " << (prefix_ ?
"prefix" :
"suffix") <<
" decoy string '" << decoy_string_ <<
"'" << std::endl;
216 bool xtandem_fix_parameters =
true;
217 bool msgfplus_fix_parameters =
true;
220 for (
const auto& prot_id : prot_ids)
222 String search_engine = prot_id.getSearchEngine();
224 if (search_engine !=
"XTANDEM") { xtandem_fix_parameters =
false; }
225 if (!(search_engine ==
"MSGFPLUS" || search_engine ==
"MS-GF+")) { msgfplus_fix_parameters =
false; }
229 if (msgfplus_fix_parameters && enzyme.
getEnzymeName() ==
"Trypsin")
231 OPENMS_LOG_WARN <<
"MSGFPlus detected but enzyme cutting rules were set to Trypsin. Correcting to Trypsin/P to copy with special cutting rule in MSGFPlus." << std::endl;
239 const size_t PROTEIN_CACHE_SIZE = 4e5;
241 this->startProgress(0, 1,
"Load first DB chunk");
242 proteins.cacheChunk(PROTEIN_CACHE_SIZE);
245 if (proteins.empty())
247 OPENMS_LOG_ERROR <<
"Error: An empty database was provided. Mapping makes no sense. Aborting..." << std::endl;
248 return DATABASE_EMPTY;
253 OPENMS_LOG_WARN <<
"Warning: An empty set of peptide identifications was provided. Output will be empty as well." << std::endl;
254 if (!keep_unreferenced_proteins_)
257 for (std::vector<ProteinIdentification>::iterator it = prot_ids.begin();
258 it != prot_ids.end(); ++it)
260 it->getHits().clear();
263 return PEPTIDE_IDS_EMPTY;
268 std::vector<bool> protein_is_decoy;
269 std::vector<std::string> protein_accessions;
271 bool invalid_protein_sequence =
false;
278 bool has_illegal_AAs(
false);
280 for (std::vector<PeptideIdentification>::const_iterator it1 = pep_ids.begin(); it1 != pep_ids.end(); ++it1)
283 const std::vector<PeptideHit>& hits = it1->getHits();
284 for (std::vector<PeptideHit>::const_iterator it2 = hits.begin(); it2 != hits.end(); ++it2)
290 String seq = it2->getSequence().toUnmodifiedString().
remove(
'*');
293 OPENMS_LOG_ERROR <<
"Peptide sequence '" << it2->getSequence() <<
"' contains one or more ambiguous amino acids (B|J|Z|X).\n";
294 has_illegal_AAs =
true;
300 appendValue(pep_DB, seq.c_str());
305 OPENMS_LOG_ERROR <<
"One or more peptides contained illegal amino acids. This is not allowed!"
306 <<
"\nPlease either remove the peptide or replace it with one of the unambiguous ones (while allowing for ambiguous AA's to match the protein)." << std::endl;;
309 OPENMS_LOG_INFO <<
"Mapping " << length(pep_DB) <<
" peptides to " << (proteins.size() == PROTEIN_CACHE_SIZE ?
"? (unknown number of)" :
String(proteins.size())) <<
" proteins." << std::endl;
311 if (length(pep_DB) == 0)
313 OPENMS_LOG_WARN <<
"Warning: Peptide identifications have no hits inside! Output will be empty as well." << std::endl;
314 return PEPTIDE_IDS_EMPTY;
320 OPENMS_LOG_INFO <<
"Searching with up to " << aaa_max_ <<
" ambiguous amino acid(s) and " << mm_max_ <<
" mismatch(es)!" << std::endl;
331 uint16_t count_j_proteins(0);
332 bool has_active_data =
true;
333 const std::string jumpX(aaa_max_ + mm_max_ + 1,
'X');
335 this->startProgress(0, proteins.size() == PROTEIN_CACHE_SIZE ? std::numeric_limits<SignedSize>::max() : proteins.size(),
"Aho-Corasick");
336 std::atomic<int> progress_prots(0);
348 #pragma omp barrier // all threads need to be here, since we are about to swap protein data
351 DEBUG_ONLY std::cerr <<
" activating cache ...\n";
352 has_active_data = proteins.activateCache();
353 protein_accessions.resize(proteins.getChunkOffset() + proteins.chunkSize());
356 if (!has_active_data)
break;
361 DEBUG_ONLY std::cerr <<
"Filling Protein Cache ...";
362 proteins.cacheChunk(PROTEIN_CACHE_SIZE);
363 protein_is_decoy.resize(proteins.getChunkOffset() + prot_count);
366 const String& seq = proteins.chunkAt(i).identifier;
367 protein_is_decoy[i + proteins.getChunkOffset()] = (prefix_ ? seq.
hasPrefix(decoy_string_) : seq.
hasSuffix(decoy_string_));
371 DEBUG_ONLY std::cerr <<
" starting for loop \n";
373 #pragma omp for schedule(dynamic, 100) nowait
377 if (omp_get_thread_num() == 0)
379 this->setProgress(progress_prots);
382 prot = proteins.chunkAt(i).sequence;
386 if (prot.
has(
'[') || prot.
has(
'('))
388 invalid_protein_sequence =
true;
407 Size prot_idx = i + proteins.getChunkOffset();
416 size_t offset = -1, start = 0;
417 while ((offset = prot.find(jumpX, offset + 1)) != std::string::npos)
420 addHits_(fuzzyAC, pattern, pep_DB, prot.
substr(start, offset + jumpX.size() - start), prot, prot_idx, (
int)start, func_threads);
422 while (offset + jumpX.size() < prot.size() && prot[offset + jumpX.size()] ==
'X') ++offset;
427 if (start < prot.size())
429 addHits_(fuzzyAC, pattern, pep_DB, prot.
substr(start), prot, prot_idx, (
int)start, func_threads);
434 addHits_(fuzzyAC, pattern, pep_DB, prot, prot, prot_idx, 0, func_threads);
439 protein_accessions[prot_idx] = proteins.chunkAt(i).identifier;
440 acc_to_prot_thread[protein_accessions[prot_idx]] = prot_idx;
447 #pragma omp critical(PeptideIndexer_joinAC)
452 func.
merge(func_threads);
454 acc_to_prot.insert(acc_to_prot_thread.begin(), acc_to_prot_thread.end());
455 acc_to_prot_thread.clear();
461 std::cout <<
"Merge took: " << s.
toString() <<
"\n";
463 std::cout << mu.
delta(
"Aho-Corasick") <<
"\n\n";
469 <<
" ... rejected by enzyme filter: " << func.
filter_rejected << std::endl;
471 if (count_j_proteins)
473 OPENMS_LOG_WARN <<
"PeptideIndexer found " << count_j_proteins <<
" protein sequences in your database containing the amino acid 'J'."
474 <<
"To match 'J' in a protein, an ambiguous amino acid placeholder for I/L will be used.\n"
475 <<
"This costs runtime and eats into the 'aaa_max' limit, leaving less opportunity for B/Z/X matches.\n"
476 <<
"If you want 'J' to be treated as unambiguous, enable '-IL_equivalent'!" << std::endl;
486 for (
Size run_idx = 0; run_idx < prot_ids.size(); ++run_idx)
488 runid_to_runidx[prot_ids[run_idx].getIdentifier()] = run_idx;
492 Size stats_matched_unique(0);
493 Size stats_matched_multi(0);
494 Size stats_unmatched(0);
495 Size stats_count_m_t(0);
496 Size stats_count_m_d(0);
497 Size stats_count_m_td(0);
502 for (std::vector<PeptideIdentification>::iterator it1 = pep_ids.begin(); it1 != pep_ids.end(); ++it1)
505 Size run_idx = runid_to_runidx[it1->getIdentifier()];
507 std::vector<PeptideHit>& hits = it1->getHits();
509 for (std::vector<PeptideHit>::iterator it2 = hits.begin(); it2 != hits.end(); ++it2)
512 it2->setPeptideEvidences(std::vector<PeptideEvidence>());
517 bool matches_target(
false);
518 bool matches_decoy(
false);
520 std::set<Size> prot_indices;
522 for (std::set<PeptideProteinMatchInformation>::const_iterator it_i = func.
pep_to_prot[pep_idx].begin();
525 prot_indices.insert(it_i->protein_index);
526 const String& accession = protein_accessions[it_i->protein_index];
527 PeptideEvidence pe(accession, it_i->position, it_i->position + (
int)it2->getSequence().size() - 1, it_i->AABefore, it_i->AAAfter);
528 it2->addPeptideEvidence(pe);
530 runidx_to_protidx[run_idx].insert(it_i->protein_index);
532 if (protein_is_decoy[it_i->protein_index])
534 matches_decoy =
true;
538 matches_target =
true;
542 if (matches_decoy && matches_target)
544 it2->setMetaValue(
"target_decoy",
"target+decoy");
547 else if (matches_target)
549 it2->setMetaValue(
"target_decoy",
"target");
552 else if (matches_decoy)
554 it2->setMetaValue(
"target_decoy",
"decoy");
559 if (prot_indices.size() == 1)
561 it2->setMetaValue(
"protein_references",
"unique");
562 ++stats_matched_unique;
564 else if (prot_indices.size() > 1)
566 it2->setMetaValue(
"protein_references",
"non-unique");
567 ++stats_matched_multi;
571 it2->setMetaValue(
"protein_references",
"unmatched");
573 if (stats_unmatched < 15)
OPENMS_LOG_INFO <<
"Unmatched peptide: " << it2->getSequence() <<
"\n";
574 else if (stats_unmatched == 15)
OPENMS_LOG_INFO <<
"Unmatched peptide: ...\n";
582 Size total_peptides = stats_count_m_t + stats_count_m_d + stats_count_m_td + stats_unmatched;
586 OPENMS_LOG_INFO <<
" unmatched : " << stats_unmatched <<
" (" << stats_unmatched * 100 / total_peptides <<
" %)\n";
588 OPENMS_LOG_INFO <<
" match to target DB only: " << stats_count_m_t <<
" (" << stats_count_m_t * 100 / total_peptides <<
" %)\n";
589 OPENMS_LOG_INFO <<
" match to decoy DB only : " << stats_count_m_d <<
" (" << stats_count_m_d * 100 / total_peptides <<
" %)\n";
590 OPENMS_LOG_INFO <<
" match to both : " << stats_count_m_td <<
" (" << stats_count_m_td * 100 / total_peptides <<
" %)\n";
593 OPENMS_LOG_INFO <<
" no match (to 0 protein) : " << stats_unmatched <<
"\n";
594 OPENMS_LOG_INFO <<
" unique match (to 1 protein) : " << stats_matched_unique <<
"\n";
595 OPENMS_LOG_INFO <<
" non-unique match (to >1 protein): " << stats_matched_multi << std::endl;
598 Size stats_matched_proteins(0), stats_matched_new_proteins(0), stats_orphaned_proteins(0), stats_proteins_target(0), stats_proteins_decoy(0);
601 for (
Size run_idx = 0; run_idx < prot_ids.size(); ++run_idx)
603 std::set<Size> masterset = runidx_to_protidx[run_idx];
605 std::vector<ProteinHit>& phits = prot_ids[run_idx].getHits();
608 std::vector<ProteinHit> orphaned_hits;
609 for (std::vector<ProteinHit>::iterator p_hit = phits.begin(); p_hit != phits.end(); ++p_hit)
611 const String& acc = p_hit->getAccession();
612 if (!acc_to_prot.
has(acc))
614 ++stats_orphaned_proteins;
615 if (keep_unreferenced_proteins_)
617 p_hit->setMetaValue(
"target_decoy",
"");
618 orphaned_hits.push_back(*p_hit);
623 phits = orphaned_hits;
628 phits.reserve(phits.size() + masterset.size());
629 for (std::set<Size>::const_iterator it = masterset.begin(); it != masterset.end(); ++it)
634 if (write_protein_sequence_ || write_protein_description_)
636 proteins.readAt(fe, *it);
637 if (write_protein_sequence_)
641 if (write_protein_description_)
646 if (protein_is_decoy[*it])
649 ++stats_proteins_decoy;
654 ++stats_proteins_target;
656 phits.push_back(hit);
657 ++stats_matched_new_proteins;
659 stats_matched_proteins += phits.size();
666 OPENMS_LOG_INFO <<
" total proteins searched: " << proteins.size() <<
"\n";
667 OPENMS_LOG_INFO <<
" matched proteins : " << stats_matched_proteins <<
" (" << stats_matched_new_proteins <<
" new)\n";
668 if (stats_matched_proteins)
670 OPENMS_LOG_INFO <<
" matched target proteins: " << stats_proteins_target <<
" (" << stats_proteins_target * 100 / stats_matched_proteins <<
" %)\n";
671 OPENMS_LOG_INFO <<
" matched decoy proteins : " << stats_proteins_decoy <<
" (" << stats_proteins_decoy * 100 / stats_matched_proteins <<
" %)\n";
673 OPENMS_LOG_INFO <<
" orphaned proteins : " << stats_orphaned_proteins << (keep_unreferenced_proteins_ ?
" (all kept)" :
" (all removed)\n");
678 bool has_error =
false;
680 if (invalid_protein_sequence)
682 OPENMS_LOG_ERROR <<
"Error: One or more protein sequences contained the characters '[' or '(', which are illegal in protein sequences."
683 <<
"\nPeptide hits might be masked by these characters (which usually indicate presence of modifications).\n";
687 if ((stats_count_m_d + stats_count_m_td) == 0)
689 String msg(
"No peptides were matched to the decoy portion of the database! Did you provide the correct concatenated database? Are your 'decoy_string' (=" +
String(decoy_string_) +
") and 'decoy_string_position' (=" +
String(param_.getValue(
"decoy_string_position")) +
") settings correct?");
690 if (missing_decoy_action_ ==
"error")
692 OPENMS_LOG_ERROR <<
"Error: " << msg <<
"\nSet 'missing_decoy_action' to 'warn' if you are sure this is ok!\nAborting ..." << std::endl;
695 else if (missing_decoy_action_ ==
"warn")
697 OPENMS_LOG_WARN <<
"Warn: " << msg <<
"\nSet 'missing_decoy_action' to 'error' if you want to elevate this to an error!" << std::endl;
704 if ((!allow_unmatched_) && (stats_unmatched > 0))
706 OPENMS_LOG_ERROR <<
"PeptideIndexer found unmatched peptides, which could not be associated to a protein.\n"
707 <<
"Potential solutions:\n"
708 <<
" - check your FASTA database for completeness\n"
709 <<
" - set 'enzyme:specificity' to match the identification parameters of the search engine\n"
710 <<
" - some engines (e.g. X! Tandem) employ loose cutting rules generating non-tryptic peptides;\n"
711 <<
" if you trust them, disable enzyme specificity\n"
712 <<
" - increase 'aaa_max' to allow more ambiguous amino acids\n"
713 <<
" - as a last resort: use the 'allow_unmatched' option to accept unmatched peptides\n"
714 <<
" (note that unmatched peptides cannot be used for FDR calculation or quantification)\n";
720 OPENMS_LOG_ERROR <<
"Result files will be written, but PeptideIndexer will exit with an error code." << std::endl;
721 return UNEXPECTED_RESULT;
726 const String& getDecoyString()
const;
728 bool isPrefix()
const;
739 const std::tuple<const Size&, const Int&, const char&, const char&>
tie()
const
741 return std::tie(protein_index,
position, AABefore, AAAfter);
745 return tie() < other.
tie();
749 return tie() == other.
tie();
756 typedef std::map<OpenMS::Size, std::set<PeptideProteinMatchInformation> >
MapType;
767 pep_to_prot(), filter_passed(0), filter_rejected(0), enzyme_(enzyme), xtandem_(xtandem)
773 if (pep_to_prot.empty())
779 for (FoundProteinFunctor::MapType::const_iterator it = other.
pep_to_prot.begin(); it != other.
pep_to_prot.end(); ++it)
781 this->pep_to_prot[it->first].insert(other.
pep_to_prot[it->first].begin(), other.
pep_to_prot[it->first].end());
805 pep_to_prot[idx_pep].insert(match);
823 const seqan::Peptide& tmp_pep = pep_DB[fuzzyAC.
getHitDBIndex()];
828 void updateMembers_()
override;
Base class for TOPP applications.
Definition: TOPPBase.h:144
Param copy(const String &prefix, bool remove_prefix=false) const
Returns a new Param object containing all entries that start with prefix.
void addHits_(AhoCorasickAmbiguous &fuzzyAC, const AhoCorasickAmbiguous::FuzzyACPattern &pattern, const AhoCorasickAmbiguous::PeptideDB &pep_DB, const String &prot, const String &full_prot, SignedSize idx_prot, Int offset, FoundProteinFunctor &func_threads) const
Definition: PeptideIndexing.h:818
Definition: PeptideIndexing.h:134
String enzyme_name_
Definition: PeptideIndexing.h:833
A convenience class to report either absolute or delta (between two timepoints) RAM usage.
Definition: SysInfo.h:83
void setDescription(const String &description)
sets the description of the protein
void store(const String &filename, const std::vector< ProteinIdentification > &protein_ids, const std::vector< PeptideIdentification > &peptide_ids, const String &document_id="")
Stores the data in an idXML file.
void setSequence(const String &sequence)
sets the protein sequence
void setEnzyme(const String &name)
Sets the enzyme for the digestion (by name)
OpenMS::Size filter_passed
Definition: PeptideIndexing.h:758
bool findNext(const FuzzyACPattern &pattern)
Enumerate hits.
Definition: AhoCorasickAmbiguous.h:1037
void setProtein(const String &protein_sequence)
Reset to new protein sequence. All previous data is forgotten.
Definition: AhoCorasickAmbiguous.h:1024
ProteaseDigestion enzyme_
Definition: PeptideIndexing.h:762
String & substitute(char from, char to)
Replaces all occurrences of the character from by the character to.
bool has(const Key &key) const
Test whether the map contains the given key.
Definition: Map.h:108
String sequence
Definition: FASTAFile.h:80
String decoy_string_
Definition: PeptideIndexing.h:830
A more convenient string class.
Definition: String.h:58
void addHit(const OpenMS::Size idx_pep, const OpenMS::Size idx_prot, const OpenMS::Size len_pep, const OpenMS::String &seq_prot, OpenMS::Int position)
Definition: PeptideIndexing.h:792
String description
Definition: FASTAFile.h:79
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
Representation of a protein hit.
Definition: ProteinHit.h:57
Int getHitProteinPosition()
Offset into protein sequence where hit was found.
Definition: AhoCorasickAmbiguous.h:1057
String< AAcid, Alloc< void > > AAString
Definition: AhoCorasickAmbiguous.h:206
Extended Aho-Corasick algorithm capable of matching ambiguous amino acids in the pattern (i....
Definition: AhoCorasickAmbiguous.h:970
const DataValue & getValue(const String &key) const
Returns a value of a parameter.
ExitCodes run(std::vector< FASTAFile::FASTAEntry > &proteins, std::vector< ProteinIdentification > &prot_ids, std::vector< PeptideIdentification > &pep_ids)
forward for old interface and pyOpenMS; use run<T>() for more control
Definition: PeptideIndexing.h:147
String delta(const String &event="delta")
#define DEBUG_ONLY
Definition: AhoCorasickAmbiguous.h:46
static void initPattern(const PeptideDB &pep_db, const int aaa_max, const int mm_max, FuzzyACPattern &pattern)
Construct a trie from a set of peptide sequences (which are to be found in a protein).
Definition: AhoCorasickAmbiguous.h:991
#define OPENMS_LOG_WARN
Macro if a warning, a piece of information which should be read by the user, should be logged.
Definition: LogStream.h:460
static Specificity getSpecificityByName(const String &name)
OpenMS::Size filter_rejected
Definition: PeptideIndexing.h:759
FASTA entry type (identifier, description and sequence)
Definition: FASTAFile.h:76
Definition: PeptideIndexing.h:136
bool hasPrefix(const String &string) const
true if String begins with string, false otherwise
MapType pep_to_prot
Definition: PeptideIndexing.h:757
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:91
void load(const String &filename, std::vector< ProteinIdentification > &protein_ids, std::vector< PeptideIdentification > &peptide_ids)
Loads the identifications of an idXML file without identifier.
void setAccession(const String &accession)
sets the accession of the protein
static const char N_TERMINAL_AA
Definition: PeptideEvidence.h:60
ExitCodes run(FASTAContainer< T > &proteins, std::vector< ProteinIdentification > &prot_ids, std::vector< PeptideIdentification > &pep_ids)
Re-index peptide identifications honoring enzyme cutting rules, ambiguous amino acids and target/deco...
Definition: PeptideIndexing.h:189
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
bool toBool() const
Conversion to bool.
Size getHitDBIndex()
Get index of hit into peptide database of the pattern.
Definition: AhoCorasickAmbiguous.h:1047
void merge(FoundProteinFunctor &other)
Definition: PeptideIndexing.h:771
FoundProteinFunctor(const ProteaseDigestion &enzyme, bool xtandem)
Definition: PeptideIndexing.h:766
Base class for all classes that want to report their progress.
Definition: ProgressLogger.h:54
FASTAContainer<TFI_Vector> simply takes an existing vector of FASTAEntries and provides the same inte...
Definition: FASTAContainer.h:243
double getClockTime() const
Int aaa_max_
Definition: PeptideIndexing.h:842
bool has(Byte byte) const
true if String contains the byte, false otherwise
Representation of a peptide evidence.
Definition: PeptideEvidence.h:50
bool write_protein_sequence_
Definition: PeptideIndexing.h:836
::seqan::Pattern< PeptideDB, ::seqan::FuzzyAC > FuzzyACPattern
Definition: AhoCorasickAmbiguous.h:974
void setParameters(const Param ¶m)
Sets the parameters.
bool hasSuffix(const String &string) const
true if String ends with string, false otherwise
String missing_decoy_action_
Definition: PeptideIndexing.h:832
const Param & getParameters() const
Non-mutable access to the parameters.
::seqan::StringSet<::seqan::AAString > PeptideDB
Definition: AhoCorasickAmbiguous.h:973
String toString() const
get a compact representation of the current time status.
ExitCodes
Exit codes.
Definition: PeptideIndexing.h:130
void after()
record data for the second timepoint
StopWatch Class.
Definition: StopWatch.h:59
bool prefix_
Definition: PeptideIndexing.h:831
bool isValidProduct(const String &protein, int pep_pos, int pep_length, bool ignore_missed_cleavages=true, bool allow_nterm_protein_cleavage=false, bool allow_random_asp_pro_cleavage=false) const
Variant of EnzymaticDigestion::isValidProduct() with support for n-term protein cleavage and random D...
int main(int argc, const char **argv)
Definition: INIFileEditor.cpp:73
std::map< OpenMS::Size, std::set< PeptideProteinMatchInformation > > MapType
Definition: PeptideIndexing.h:756
bool write_protein_description_
Definition: PeptideIndexing.h:837
bool update(const Param &p_outdated, const bool add_unknown=false)
Rescue parameter values from p_outdated to current param.
ptrdiff_t SignedSize
Signed Size type e.g. used as pointer difference.
Definition: Types.h:134
String getEnzymeName() const
Returns the enzyme for the digestion.
#define OPENMS_LOG_ERROR
Macro to be used if non-fatal error are reported (processing continues)
Definition: LogStream.h:455
Logger::LogStream OpenMS_Log_debug
Global static instance of a LogStream to capture messages classified as debug output....
String substr(size_t pos=0, size_t n=npos) const
Wrapper for the STL substr() method. Returns a String object with its contents initialized to a subst...
Int mm_max_
Definition: PeptideIndexing.h:843
Class for the enzymatic digestion of proteins.
Definition: ProteaseDigestion.h:60
static String findDatabase(const String &db_name)
Definition: PeptideIndexing.h:133
Refreshes the protein references for all peptide hits in a vector of PeptideIdentifications and adds ...
Definition: PeptideIndexing.h:124
static bool readable(const String &file)
Return true if the file exists and is readable.
static const char C_TERMINAL_AA
Definition: PeptideEvidence.h:61
bool allow_unmatched_
Definition: PeptideIndexing.h:839
bool xtandem_
Definition: PeptideIndexing.h:763
String enzyme_specificity_
Definition: PeptideIndexing.h:834
Management and storage of parameters / INI files.
Definition: Param.h:73
Definition: PeptideIndexing.h:753
String & remove(char what)
Remove all occurrences of the character what.
bool isAmbiguous(AAcid c)
Definition: AhoCorasickAmbiguous.h:578
bool IL_equivalent_
Definition: PeptideIndexing.h:840
Map class based on the STL map (containing several convenience functions)
Definition: Map.h:50
#define OPENMS_LOG_INFO
Macro if a information, e.g. a status should be reported.
Definition: LogStream.h:465
void setSpecificity(Specificity spec)
Sets the specificity for the digestion (default is SPEC_FULL).
static String & toUpper(String &this_s)
Definition: StringUtils.h:832
Definition: PeptideIndexing.h:132
Definition: PeptideIndexing.h:135
Size< TNeedle >::Type position(const PatternAuxData< TNeedle > &dh)
Definition: AhoCorasickAmbiguous.h:561
static Result findDecoyString(FASTAContainer< T > &proteins)
Heuristic to determine the decoy string given a set of protein names.
Definition: FASTAContainer.h:359
void setLogType(LogType type) const
Sets the progress log that should be used. The default type is NONE!
template parameter for vector-based FASTA access
Definition: FASTAContainer.h:82
Used to load and store idXML files.
Definition: IdXMLFile.h:63
FASTAContainer<TFI_File> will make FASTA entries available chunk-wise from start to end by loading it...
Definition: FASTAContainer.h:93
bool keep_unreferenced_proteins_
Definition: PeptideIndexing.h:838