OpenMS  2.4.0
PeptideIndexing.h
Go to the documentation of this file.
1 // --------------------------------------------------------------------------
2 // OpenMS -- Open-Source Mass Spectrometry
3 // --------------------------------------------------------------------------
4 // Copyright The OpenMS Team -- Eberhard Karls University Tuebingen,
5 // ETH Zurich, and Freie Universitaet Berlin 2002-2018.
6 //
7 // This software is released under a three-clause BSD license:
8 // * Redistributions of source code must retain the above copyright
9 // notice, this list of conditions and the following disclaimer.
10 // * Redistributions in binary form must reproduce the above copyright
11 // notice, this list of conditions and the following disclaimer in the
12 // documentation and/or other materials provided with the distribution.
13 // * Neither the name of any author or any participating institution
14 // may be used to endorse or promote products derived from this software
15 // without specific prior written permission.
16 // For a full list of authors, refer to the file AUTHORS.
17 // --------------------------------------------------------------------------
18 // THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
19 // AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
20 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
21 // ARE DISCLAIMED. IN NO EVENT SHALL ANY OF THE AUTHORS OR THE CONTRIBUTING
22 // INSTITUTIONS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
23 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
24 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
25 // OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
26 // WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
27 // OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
28 // ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
29 //
30 // --------------------------------------------------------------------------
31 // $Maintainer: Chris Bielow $
32 // $Authors: Andreas Bertsch, Chris Bielow $
33 // --------------------------------------------------------------------------
34 
35 #pragma once
36 
37 
54 #include <OpenMS/SYSTEM/SysInfo.h>
55 
56 #include <atomic>
57 #include <algorithm>
58 #include <fstream>
59 
60 
61 namespace OpenMS
62 {
63 
124  class OPENMS_DLLAPI PeptideIndexing :
125  public DefaultParamHandler, public ProgressLogger
126  {
127 public:
128 
131  {
138  };
139 
141  PeptideIndexing();
142 
144  ~PeptideIndexing() override;
145 
146 
148  inline ExitCodes run(std::vector<FASTAFile::FASTAEntry>& proteins, std::vector<ProteinIdentification>& prot_ids, std::vector<PeptideIdentification>& pep_ids)
149  {
150  FASTAContainer<TFI_Vector> protein_container(proteins);
151  return run<TFI_Vector>(protein_container, prot_ids, pep_ids);
152  }
153 
189  template<typename T>
190  ExitCodes run(FASTAContainer<T>& proteins, std::vector<ProteinIdentification>& prot_ids, std::vector<PeptideIdentification>& pep_ids)
191  {
192  // no decoy string provided? try to deduce from data
193  if (decoy_string_.empty())
194  {
195  bool is_decoy_string_auto_successful = findDecoyString_(proteins);
196 
197  if (!is_decoy_string_auto_successful && contains_decoys_)
198  {
199  return DECOYSTRING_EMPTY;
200  }
201  else if (!is_decoy_string_auto_successful && !contains_decoys_)
202  {
203  LOG_WARN << "Unable to determine decoy string automatically, not enough decoys were detected! Using default " << (prefix_ ? "prefix" : "suffix") << " decoy string '" << decoy_string_ << "\n"
204  << "If you think that this is false, please provide a decoy_string and its position manually!" << std::endl;
205  }
206  else
207  {
208  // decoy string and position was extracted successfully
209  LOG_INFO << "Using " << (prefix_ ? "prefix" : "suffix") << " decoy string '" << decoy_string_ << "'" << std::endl;
210  }
211 
212  proteins.reset();
213  }
214 
215  //---------------------------------------------------------------
216  // parsing parameters, correcting xtandem and MSGFPlus parameters
217  //---------------------------------------------------------------
218  ProteaseDigestion enzyme;
219  enzyme.setEnzyme(enzyme_name_);
220  enzyme.setSpecificity(enzyme.getSpecificityByName(enzyme_specificity_));
221 
222  bool xtandem_fix_parameters = true, msgfplus_fix_parameters = true;
223 
224  // specificity is none or semi? don't automate xtandem
227  {
228  xtandem_fix_parameters = false;
229  }
230 
231  // enzyme is already Trypsin/P? don't automate MSGFPlus
232  if (enzyme.getEnzymeName() == "Trypsin/P") { msgfplus_fix_parameters = false; }
233 
234  // determine if search engine is solely xtandem or MSGFPlus
235  for (const auto& prot_id : prot_ids)
236  {
237  if (!msgfplus_fix_parameters && !xtandem_fix_parameters) { break; }
238  String se = prot_id.getSearchEngine();
239  std::string search_engine = StringUtils::toUpper(se);
240  if (search_engine != "XTANDEM") { xtandem_fix_parameters = false; }
241  if (search_engine != "MSGFPLUS" || "MS-GF+") { msgfplus_fix_parameters = false; }
242  }
243 
244  // solely MSGFPlus -> Trypsin P as enzyme
245  if (msgfplus_fix_parameters && enzyme.getEnzymeName() == "Trypsin")
246  {
247  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;
248  enzyme.setEnzyme("Trypsin/P");
249  }
250 
251  //-------------------------------------------------------------
252  // calculations
253  //-------------------------------------------------------------
254  // cache the first proteins
255  const size_t PROTEIN_CACHE_SIZE = 4e5; // 400k should be enough for most DB's and is not too hard on memory either (~200 MB FASTA)
256 
257  proteins.cacheChunk(PROTEIN_CACHE_SIZE);
258 
259  if (proteins.empty()) // we do not allow an empty database
260  {
261  LOG_ERROR << "Error: An empty database was provided. Mapping makes no sense. Aborting..." << std::endl;
262  return DATABASE_EMPTY;
263  }
264 
265  if (pep_ids.empty()) // Aho-Corasick requires non-empty input; but we allow this case, since the TOPP tool should not crash when encountering a bad raw file (with no PSMs)
266  {
267  LOG_WARN << "Warning: An empty set of peptide identifications was provided. Output will be empty as well." << std::endl;
268  if (!keep_unreferenced_proteins_)
269  {
270  // delete only protein hits, not whole ID runs incl. meta data:
271  for (std::vector<ProteinIdentification>::iterator it = prot_ids.begin();
272  it != prot_ids.end(); ++it)
273  {
274  it->getHits().clear();
275  }
276  }
277  return PEPTIDE_IDS_EMPTY;
278  }
279 
280  FoundProteinFunctor func(enzyme, xtandem_fix_parameters); // store the matches
281  Map<String, Size> acc_to_prot; // map: accessions --> FASTA protein index
282  std::vector<bool> protein_is_decoy; // protein index -> is decoy?
283  std::vector<std::string> protein_accessions; // protein index -> accession
284 
285  bool invalid_protein_sequence = false; // check for proteins with modifications, i.e. '[' or '(', and throw an exception
286 
287  { // new scope - forget data after search
288 
289  /*
290  BUILD Peptide DB
291  */
292  bool has_illegal_AAs(false);
294  for (std::vector<PeptideIdentification>::const_iterator it1 = pep_ids.begin(); it1 != pep_ids.end(); ++it1)
295  {
296  //String run_id = it1->getIdentifier();
297  const std::vector<PeptideHit>& hits = it1->getHits();
298  for (std::vector<PeptideHit>::const_iterator it2 = hits.begin(); it2 != hits.end(); ++it2)
299  {
300  //
301  // Warning:
302  // do not skip over peptides here, since the results are iterated in the same way
303  //
304  String seq = it2->getSequence().toUnmodifiedString().remove('*'); // make a copy, i.e. do NOT change the peptide sequence!
305  if (seqan::isAmbiguous(seqan::AAString(seq.c_str())))
306  { // do not quit here, to show the user all sequences .. only quit after loop
307  LOG_ERROR << "Peptide sequence '" << it2->getSequence() << "' contains one or more ambiguous amino acids (B|J|Z|X).\n";
308  has_illegal_AAs = true;
309  }
310  if (IL_equivalent_) // convert L to I;
311  {
312  seq.substitute('L', 'I');
313  }
314  appendValue(pep_DB, seq.c_str());
315  }
316  }
317  if (has_illegal_AAs)
318  {
319  LOG_ERROR << "One or more peptides contained illegal amino acids. This is not allowed!"
320  << "\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;;
321  }
322 
323  LOG_INFO << "Mapping " << length(pep_DB) << " peptides to " << (proteins.size() == PROTEIN_CACHE_SIZE ? "? (unknown number of)" : String(proteins.size())) << " proteins." << std::endl;
324 
325  if (length(pep_DB) == 0)
326  { // Aho-Corasick will crash if given empty needles as input
327  LOG_WARN << "Warning: Peptide identifications have no hits inside! Output will be empty as well." << std::endl;
328  return PEPTIDE_IDS_EMPTY;
329  }
330 
331  /*
332  Aho Corasick (fast)
333  */
334  LOG_INFO << "Searching with up to " << aaa_max_ << " ambiguous amino acid(s) and " << mm_max_ << " mismatch(es)!" << std::endl;
336  LOG_INFO << "Building trie ...";
337  StopWatch s;
338  s.start();
340  AhoCorasickAmbiguous::initPattern(pep_DB, aaa_max_, mm_max_, pattern);
341  s.stop();
342  LOG_INFO << " done (" << int(s.getClockTime()) << "s)" << std::endl;
343  s.reset();
344 
345  uint16_t count_j_proteins(0);
346  bool has_active_data = true; // becomes false if end of FASTA file is reached
347  const std::string jumpX(aaa_max_ + mm_max_ + 1, 'X'); // jump over stretches of 'X' which cost a lot of time; +1 because AXXA is a valid hit for aaa_max == 2 (cannot split it)
348  this->startProgress(0, proteins.size(), "Aho-Corasick");
349  std::atomic<int> progress_prots(0);
350 #ifdef _OPENMP
351 #pragma omp parallel
352 #endif
353  {
354  FoundProteinFunctor func_threads(enzyme, xtandem_fix_parameters);
355  Map<String, Size> acc_to_prot_thread; // map: accessions --> FASTA protein index
356  AhoCorasickAmbiguous fuzzyAC;
357  String prot;
358 
359  while (true)
360  {
361  #pragma omp barrier // all threads need to be here, since we are about to swap protein data
362  #pragma omp single
363  {
364  DEBUG_ONLY std::cerr << " activating cache ...\n";
365  has_active_data = proteins.activateCache(); // swap in last cache
366  protein_accessions.resize(proteins.getChunkOffset() + proteins.chunkSize());
367  } // implicit barrier here
368 
369  if (!has_active_data) break; // leave while-loop
370  SignedSize prot_count = (SignedSize)proteins.chunkSize();
371 
372  #pragma omp master
373  {
374  DEBUG_ONLY std::cerr << "Filling Protein Cache ...";
375  proteins.cacheChunk(PROTEIN_CACHE_SIZE);
376  protein_is_decoy.resize(proteins.getChunkOffset() + prot_count);
377  for (SignedSize i = 0; i < prot_count; ++i)
378  { // do this in master only, to avoid false sharing
379  const String& seq = proteins.chunkAt(i).identifier;
380  protein_is_decoy[i + proteins.getChunkOffset()] = (prefix_ ? seq.hasPrefix(decoy_string_) : seq.hasSuffix(decoy_string_));
381  }
382  DEBUG_ONLY std::cerr << " done" << std::endl;
383  }
384  DEBUG_ONLY std::cerr << " starting for loop \n";
385  // search all peptides in each protein
386  #pragma omp for schedule(dynamic, 100) nowait
387  for (SignedSize i = 0; i < prot_count; ++i)
388  {
389  ++progress_prots; // atomic
390  if (omp_get_thread_num() == 0)
391  {
392  this->setProgress(progress_prots);
393  }
394 
395  prot = proteins.chunkAt(i).sequence;
396  prot.remove('*');
397 
398  // check for invalid sequences with modifications
399  if (prot.has('[') || prot.has('('))
400  {
401  invalid_protein_sequence = true; // not omp-critical because its write-only
402  // we cannot throw an exception here, since we'd need to catch it within the parallel region
403  }
404 
405  // convert L/J to I; also replace 'J' in proteins
406  if (IL_equivalent_)
407  {
408  prot.substitute('L', 'I');
409  prot.substitute('J', 'I');
410  }
411  else
412  { // warn if 'J' is found (it eats into aaa_max)
413  if (prot.has('J'))
414  {
415  #pragma omp atomic
416  ++count_j_proteins;
417  }
418  }
419 
420  Size prot_idx = i + proteins.getChunkOffset();
421 
422  // test if protein was a hit
423  Size hits_total = func_threads.filter_passed + func_threads.filter_rejected;
424 
425  // check if there are stretches of 'X'
426  if (prot.has('X'))
427  {
428  // create chunks of the protein (splitting it at stretches of 'X..X') and feed them to AC one by one
429  size_t offset = -1, start = 0;
430  while ((offset = prot.find(jumpX, offset + 1)) != std::string::npos)
431  {
432  //std::cout << "found X..X at " << offset << " in protein " << proteins[i].identifier << "\n";
433  addHits_(fuzzyAC, pattern, pep_DB, prot.substr(start, offset + jumpX.size() - start), prot, prot_idx, (int)start, func_threads);
434  // skip ahead while we encounter more X...
435  while (offset + jumpX.size() < prot.size() && prot[offset + jumpX.size()] == 'X') ++offset;
436  start = offset;
437  //std::cout << " new start: " << start << "\n";
438  }
439  // last chunk
440  if (start < prot.size())
441  {
442  addHits_(fuzzyAC, pattern, pep_DB, prot.substr(start), prot, prot_idx, (int)start, func_threads);
443  }
444  }
445  else
446  {
447  addHits_(fuzzyAC, pattern, pep_DB, prot, prot, prot_idx, 0, func_threads);
448  }
449  // was protein found?
450  if (hits_total < func_threads.filter_passed + func_threads.filter_rejected)
451  {
452  protein_accessions[prot_idx] = proteins.chunkAt(i).identifier;
453  acc_to_prot_thread[protein_accessions[prot_idx]] = prot_idx;
454  }
455  } // end parallel FOR
456 
457  // join results again
458  DEBUG_ONLY std::cerr << " critical now \n";
459  #ifdef _OPENMP
460  #pragma omp critical(PeptideIndexer_joinAC)
461  #endif
462  {
463  s.start();
464  // hits
465  func.merge(func_threads);
466  // accession -> index
467  acc_to_prot.insert(acc_to_prot_thread.begin(), acc_to_prot_thread.end());
468  acc_to_prot_thread.clear();
469  s.stop();
470  } // OMP end critical
471  } // end readChunk
472  } // OMP end parallel
473  this->endProgress();
474  std::cout << "Merge took: " << s.toString() << "\n";
475  mu.after();
476  std::cout << mu.delta("Aho-Corasick") << "\n\n";
477 
478  LOG_INFO << "\nAho-Corasick done:\n found " << func.filter_passed << " hits for " << func.pep_to_prot.size() << " of " << length(pep_DB) << " peptides.\n";
479 
480  // write some stats
481  LOG_INFO << "Peptide hits passing enzyme filter: " << func.filter_passed << "\n"
482  << " ... rejected by enzyme filter: " << func.filter_rejected << std::endl;
483 
484  if (count_j_proteins)
485  {
486  LOG_WARN << "PeptideIndexer found " << count_j_proteins << " protein sequences in your database containing the amino acid 'J'."
487  << "To match 'J' in a protein, an ambiguous amino acid placeholder for I/L will be used.\n"
488  << "This costs runtime and eats into the 'aaa_max' limit, leaving less opportunity for B/Z/X matches.\n"
489  << "If you want 'J' to be treated as unambiguous, enable '-IL_equivalent'!" << std::endl;
490  }
491 
492  } // end local scope
493 
494  //
495  // do mapping
496  //
497  // index existing proteins
498  Map<String, Size> runid_to_runidx; // identifier to index
499  for (Size run_idx = 0; run_idx < prot_ids.size(); ++run_idx)
500  {
501  runid_to_runidx[prot_ids[run_idx].getIdentifier()] = run_idx;
502  }
503 
504  // for peptides --> proteins
505  Size stats_matched_unique(0);
506  Size stats_matched_multi(0);
507  Size stats_unmatched(0); // no match to DB
508  Size stats_count_m_t(0); // match to Target DB
509  Size stats_count_m_d(0); // match to Decoy DB
510  Size stats_count_m_td(0); // match to T+D DB
511 
512  Map<Size, std::set<Size> > runidx_to_protidx; // in which protID do appear which proteins (according to mapped peptides)
513 
514  Size pep_idx(0);
515  for (std::vector<PeptideIdentification>::iterator it1 = pep_ids.begin(); it1 != pep_ids.end(); ++it1)
516  {
517  // which ProteinIdentification does the peptide belong to?
518  Size run_idx = runid_to_runidx[it1->getIdentifier()];
519 
520  std::vector<PeptideHit>& hits = it1->getHits();
521 
522  for (std::vector<PeptideHit>::iterator it2 = hits.begin(); it2 != hits.end(); ++it2)
523  {
524  // clear protein accessions
525  it2->setPeptideEvidences(std::vector<PeptideEvidence>());
526 
527  //
528  // is this a decoy hit?
529  //
530  bool matches_target(false);
531  bool matches_decoy(false);
532 
533  std::set<Size> prot_indices;
534  // add new protein references
535  for (std::set<PeptideProteinMatchInformation>::const_iterator it_i = func.pep_to_prot[pep_idx].begin();
536  it_i != func.pep_to_prot[pep_idx].end(); ++it_i)
537  {
538  prot_indices.insert(it_i->protein_index);
539  const String& accession = protein_accessions[it_i->protein_index];
540  PeptideEvidence pe(accession, it_i->position, it_i->position + (int)it2->getSequence().size() - 1, it_i->AABefore, it_i->AAAfter);
541  it2->addPeptideEvidence(pe);
542 
543  runidx_to_protidx[run_idx].insert(it_i->protein_index); // fill protein hits
544 
545  if (protein_is_decoy[it_i->protein_index])
546  {
547  matches_decoy = true;
548  }
549  else
550  {
551  matches_target = true;
552  }
553  }
554 
555  if (matches_decoy && matches_target)
556  {
557  it2->setMetaValue("target_decoy", "target+decoy");
558  ++stats_count_m_td;
559  }
560  else if (matches_target)
561  {
562  it2->setMetaValue("target_decoy", "target");
563  ++stats_count_m_t;
564  }
565  else if (matches_decoy)
566  {
567  it2->setMetaValue("target_decoy", "decoy");
568  ++stats_count_m_d;
569  } // else: could match to no protein (i.e. both are false)
570  //else ... // not required (handled below; see stats_unmatched);
571 
572  if (prot_indices.size() == 1)
573  {
574  it2->setMetaValue("protein_references", "unique");
575  ++stats_matched_unique;
576  }
577  else if (prot_indices.size() > 1)
578  {
579  it2->setMetaValue("protein_references", "non-unique");
580  ++stats_matched_multi;
581  }
582  else
583  {
584  it2->setMetaValue("protein_references", "unmatched");
585  ++stats_unmatched;
586  if (stats_unmatched < 15) LOG_INFO << "Unmatched peptide: " << it2->getSequence() << "\n";
587  else if (stats_unmatched == 15) LOG_INFO << "Unmatched peptide: ...\n";
588  }
589 
590  ++pep_idx; // next hit
591  }
592 
593  }
594 
595  Size total_peptides = stats_count_m_t + stats_count_m_d + stats_count_m_td + stats_unmatched;
596  LOG_INFO << "-----------------------------------\n";
597  LOG_INFO << "Peptide statistics\n";
598  LOG_INFO << "\n";
599  LOG_INFO << " unmatched : " << stats_unmatched << " (" << stats_unmatched * 100 / total_peptides << " %)\n";
600  LOG_INFO << " target/decoy:\n";
601  LOG_INFO << " match to target DB only: " << stats_count_m_t << " (" << stats_count_m_t * 100 / total_peptides << " %)\n";
602  LOG_INFO << " match to decoy DB only : " << stats_count_m_d << " (" << stats_count_m_d * 100 / total_peptides << " %)\n";
603  LOG_INFO << " match to both : " << stats_count_m_td << " (" << stats_count_m_td * 100 / total_peptides << " %)\n";
604  LOG_INFO << "\n";
605  LOG_INFO << " mapping to proteins:\n";
606  LOG_INFO << " no match (to 0 protein) : " << stats_unmatched << "\n";
607  LOG_INFO << " unique match (to 1 protein) : " << stats_matched_unique << "\n";
608  LOG_INFO << " non-unique match (to >1 protein): " << stats_matched_multi << std::endl;
609 
611  Size stats_matched_proteins(0), stats_matched_new_proteins(0), stats_orphaned_proteins(0), stats_proteins_target(0), stats_proteins_decoy(0);
612 
613  // all peptides contain the correct protein hit references, now update the protein hits
614  for (Size run_idx = 0; run_idx < prot_ids.size(); ++run_idx)
615  {
616  std::set<Size> masterset = runidx_to_protidx[run_idx]; // all protein matches from above
617 
618  std::vector<ProteinHit>& phits = prot_ids[run_idx].getHits();
619  {
620  // go through existing protein hits and count orphaned proteins (with no peptide hits)
621  std::vector<ProteinHit> orphaned_hits;
622  for (std::vector<ProteinHit>::iterator p_hit = phits.begin(); p_hit != phits.end(); ++p_hit)
623  {
624  const String& acc = p_hit->getAccession();
625  if (!acc_to_prot.has(acc)) // acc_to_prot only contains found proteins from current run
626  { // old hit is orphaned
627  ++stats_orphaned_proteins;
628  if (keep_unreferenced_proteins_)
629  {
630  p_hit->setMetaValue("target_decoy", "");
631  orphaned_hits.push_back(*p_hit);
632  }
633  }
634  }
635  // only keep orphaned hits (if any)
636  phits = orphaned_hits;
637  }
638 
639  // add new protein hits
641  phits.reserve(phits.size() + masterset.size());
642  for (std::set<Size>::const_iterator it = masterset.begin(); it != masterset.end(); ++it)
643  {
644  ProteinHit hit;
645  hit.setAccession(protein_accessions[*it]);
646 
647  if (write_protein_sequence_ || write_protein_description_)
648  {
649  proteins.readAt(fe, *it);
650  if (write_protein_sequence_)
651  {
652  hit.setSequence(fe.sequence);
653  } // no else, since sequence is empty by default
654  if (write_protein_description_)
655  {
656  hit.setDescription(fe.description);
657  } // no else, since description is empty by default
658  }
659  if (protein_is_decoy[*it])
660  {
661  hit.setMetaValue("target_decoy", "decoy");
662  ++stats_proteins_decoy;
663  }
664  else
665  {
666  hit.setMetaValue("target_decoy", "target");
667  ++stats_proteins_target;
668  }
669  phits.push_back(hit);
670  ++stats_matched_new_proteins;
671  }
672  stats_matched_proteins += phits.size();
673  }
674 
675 
676  LOG_INFO << "-----------------------------------\n";
677  LOG_INFO << "Protein statistics\n";
678  LOG_INFO << "\n";
679  LOG_INFO << " total proteins searched: " << proteins.size() << "\n";
680  LOG_INFO << " matched proteins : " << stats_matched_proteins << " (" << stats_matched_new_proteins << " new)\n";
681  if (stats_matched_proteins)
682  { // prevent Division-by-0 Exception
683  LOG_INFO << " matched target proteins: " << stats_proteins_target << " (" << stats_proteins_target * 100 / stats_matched_proteins << " %)\n";
684  LOG_INFO << " matched decoy proteins : " << stats_proteins_decoy << " (" << stats_proteins_decoy * 100 / stats_matched_proteins << " %)\n";
685  }
686  LOG_INFO << " orphaned proteins : " << stats_orphaned_proteins << (keep_unreferenced_proteins_ ? " (all kept)" : " (all removed)\n");
687  LOG_INFO << "-----------------------------------" << std::endl;
688 
689 
691  bool has_error = false;
692 
693  if (invalid_protein_sequence)
694  {
695  LOG_ERROR << "Error: One or more protein sequences contained the characters '[' or '(', which are illegal in protein sequences."
696  << "\nPeptide hits might be masked by these characters (which usually indicate presence of modifications).\n";
697  has_error = true;
698  }
699 
700  if ((stats_count_m_d + stats_count_m_td) == 0)
701  {
702  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?");
703  if (missing_decoy_action_ == "error")
704  {
705  LOG_ERROR << "Error: " << msg << "\nSet 'missing_decoy_action' to 'warn' if you are sure this is ok!\nAborting ..." << std::endl;
706  has_error = true;
707  }
708  else if (missing_decoy_action_ == "warn")
709  {
710  LOG_WARN << "Warn: " << msg << "\nSet 'missing_decoy_action' to 'error' if you want to elevate this to an error!" << std::endl;
711  }
712  else // silent
713  {
714  }
715  }
716 
717  if ((!allow_unmatched_) && (stats_unmatched > 0))
718  {
719  LOG_ERROR << "PeptideIndexer found unmatched peptides, which could not be associated to a protein.\n"
720  << "Potential solutions:\n"
721  << " - check your FASTA database for completeness\n"
722  << " - set 'enzyme:specificity' to match the identification parameters of the search engine\n"
723  << " - some engines (e.g. X! Tandem) employ loose cutting rules generating non-tryptic peptides;\n"
724  << " if you trust them, disable enzyme specificity\n"
725  << " - increase 'aaa_max' to allow more ambiguous amino acids\n"
726  << " - as a last resort: use the 'allow_unmatched' option to accept unmatched peptides\n"
727  << " (note that unmatched peptides cannot be used for FDR calculation or quantification)\n";
728  has_error = true;
729  }
730 
731  if (has_error)
732  {
733  LOG_ERROR << "Result files will be written, but PeptideIndexer will exit with an error code." << std::endl;
734  return UNEXPECTED_RESULT;
735  }
736  return EXECUTION_OK;
737  }
738 
739  const String& getDecoyString() const;
740 
741  bool isPrefix() const;
742 
743  protected:
744  using DecoyStringToAffixCount = std::map<std::string, std::pair<int, int>>;
745  using CaseInsensitiveToCaseSensitiveDecoy = std::map<std::string, std::string>;
747 
748  template<typename T>
750  {
751  // common decoy strings in FASTA files
752  // note: decoy prefixes/suffices must be provided in lower case
753  std::vector<std::string> affixes = {"decoy", "dec", "reverse", "rev", "__id_decoy", "xxx", "shuffled", "shuffle", "pseudo", "random"};
754 
755  // map decoys to counts of occurrences as prefix/suffix
756  DecoyStringToAffixCount decoy_count;
757  // map case insensitive strings back to original case (as used in fasta)
758  CaseInsensitiveToCaseSensitiveDecoy decoy_case_sensitive;
759 
760  // assume that it contains decoys for now
761  contains_decoys_ = true;
762 
763  // setup prefix- and suffix regex strings
764  const std::string regexstr_prefix = std::string("^(") + ListUtils::concatenate<std::string>(affixes, "_*|") + "_*)";
765  const std::string regexstr_suffix = std::string("(") + ListUtils::concatenate<std::string>(affixes, "_*|") + "_*)$";
766 
767  // setup regexes
768  const boost::regex pattern_prefix(regexstr_prefix);
769  const boost::regex pattern_suffix(regexstr_suffix);
770 
771  int all_prefix_occur(0), all_suffix_occur(0), all_proteins_count(0);
772 
773  const size_t PROTEIN_CACHE_SIZE = 4e5;
774 
775  while (true)
776  {
777  proteins.cacheChunk(PROTEIN_CACHE_SIZE);
778  if (!proteins.activateCache()) break;
779 
780  auto prot_count = (SignedSize) proteins.chunkSize();
781  all_proteins_count += prot_count;
782 
783  {
784  for (SignedSize i = 0; i < prot_count; ++i)
785  {
786  String seq = proteins.chunkAt(i).identifier;
787 
788  String seq_lower = seq;
789  seq_lower.toLower();
790 
791  boost::smatch sm;
792  // search for prefix
793  bool found_prefix = boost::regex_search(seq_lower, sm, pattern_prefix);
794  if (found_prefix)
795  {
796  std::string match = sm[0];
797  all_prefix_occur++;
798 
799  // increase count of observed prefix
800  decoy_count[match].first++;
801 
802  // store observed (case sensitive and with special characters)
803  std::string seq_decoy = StringUtils::prefix(seq, match.length());
804  decoy_case_sensitive[match] = seq_decoy;
805  }
806 
807  // search for suffix
808  bool found_suffix = boost::regex_search(seq_lower, sm, pattern_suffix);
809  if (found_suffix)
810  {
811  std::string match = sm[0];
812  all_suffix_occur++;
813 
814  // increase count of observed suffix
815  decoy_count[match].second++;
816 
817  // store observed (case sensitive and with special characters)
818  std::string seq_decoy = StringUtils::suffix(seq, match.length());
819  decoy_case_sensitive[match] = seq_decoy;
820  }
821  }
822  }
823  }
824 
825  // DEBUG ONLY: print counts of found decoys
826  for (auto &a : decoy_count) LOG_DEBUG << a.first << "\t" << a.second.first << "\t" << a.second.second << std::endl;
827 
828  // less than 40% of proteins are decoys -> won't be able to determine a decoy string and its position
829  // return default values
830  if (all_prefix_occur + all_suffix_occur < 0.4 * all_proteins_count) {
831  decoy_string_ = "DECOY_";
832  prefix_ = true;
833 
834  contains_decoys_ = false;
835  return false;
836  }
837 
838  if (all_prefix_occur == all_suffix_occur)
839  {
840  LOG_ERROR << "Unable to determine decoy string!" << std::endl;
841  return false;
842  }
843 
844  // Decoy prefix occurred at least 80% of all prefixes + observed in at least 40% of all proteins -> set it as prefix decoy
845  for (const auto& pair : decoy_count)
846  {
847  const std::string & case_insensitive_decoy_string = pair.first;
848  const std::pair<int, int>& prefix_suffix_counts = pair.second;
849  double freq_prefix = static_cast<double>(prefix_suffix_counts.first) / static_cast<double>(all_prefix_occur);
850  double freq_prefix_in_proteins = static_cast<double>(prefix_suffix_counts.first) / static_cast<double>(all_proteins_count);
851 
852  if (freq_prefix >= 0.8 && freq_prefix_in_proteins >= 0.4)
853  {
854  prefix_ = true;
855  decoy_string_ = decoy_case_sensitive[case_insensitive_decoy_string];
856 
857  if (prefix_suffix_counts.first != all_prefix_occur)
858  {
859  LOG_WARN << "More than one decoy prefix observed!" << std::endl;
860  LOG_WARN << "Using most frequent decoy prefix (" << (int) (freq_prefix * 100) <<"%)" << std::endl;
861  }
862 
863  return true;
864  }
865  }
866 
867  // Decoy suffix occurred at least 80% of all suffixes + observed in at least 40% of all proteins -> set it as suffix decoy
868  for (const auto& pair : decoy_count)
869  {
870  const std::string& case_insensitive_decoy_string = pair.first;
871  const std::pair<int, int>& prefix_suffix_counts = pair.second;
872  double freq_suffix = static_cast<double>(prefix_suffix_counts.second) / static_cast<double>(all_suffix_occur);
873  double freq_suffix_in_proteins = static_cast<double>(prefix_suffix_counts.second) / static_cast<double>(all_proteins_count);
874 
875  if (freq_suffix >= 0.8 && freq_suffix_in_proteins >= 0.4)
876  {
877  prefix_ = false;
878  decoy_string_ = decoy_case_sensitive[case_insensitive_decoy_string];
879 
880  if (prefix_suffix_counts.second != all_suffix_occur)
881  {
882  LOG_WARN << "More than one decoy suffix observed!" << std::endl;
883  LOG_WARN << "Using most frequent decoy suffix (" << (int) (freq_suffix * 100) <<"%)" << std::endl;
884  }
885 
886  return true;
887  }
888  }
889 
890  LOG_ERROR << "Unable to determine decoy string and its position. Please provide a decoy string and its position as parameters." << std::endl;
891  return false;
892  }
893 
894 
896  {
899 
902 
904  char AABefore;
905 
907  char AAAfter;
908 
909  bool operator<(const PeptideProteinMatchInformation& other) const
910  {
911  if (protein_index != other.protein_index)
912  {
913  return protein_index < other.protein_index;
914  }
915  else if (position != other.position)
916  {
917  return position < other.position;
918  }
919  else if (AABefore != other.AABefore)
920  {
921  return AABefore < other.AABefore;
922  }
923  else if (AAAfter != other.AAAfter)
924  {
925  return AAAfter < other.AAAfter;
926  }
927  return false;
928  }
929 
931  {
932  return protein_index == other.protein_index &&
933  position == other.position &&
934  AABefore == other.AABefore &&
935  AAAfter == other.AAAfter;
936  }
937 
938  };
940  {
941  public:
942  typedef std::map<OpenMS::Size, std::set<PeptideProteinMatchInformation> > MapType;
943 
946 
949 
952 
953  private:
955 
957  bool xtandem_;
958 
959  public:
960  explicit FoundProteinFunctor(const ProteaseDigestion& enzyme, bool xtandem) :
961  pep_to_prot(), filter_passed(0), filter_rejected(0), enzyme_(enzyme), xtandem_(xtandem)
962  {
963  }
964 
966  {
967  if (pep_to_prot.empty())
968  { // first merge is easy
969  pep_to_prot.swap(other.pep_to_prot);
970  }
971  else
972  {
973  for (FoundProteinFunctor::MapType::const_iterator it = other.pep_to_prot.begin(); it != other.pep_to_prot.end(); ++it)
974  { // augment set
975  this->pep_to_prot[it->first].insert(other.pep_to_prot[it->first].begin(), other.pep_to_prot[it->first].end());
976  }
977  other.pep_to_prot.clear();
978  }
979  // cheap members
980  this->filter_passed += other.filter_passed;
981  other.filter_passed = 0;
982  this->filter_rejected += other.filter_rejected;
983  other.filter_rejected = 0;
984  }
985 
986  void addHit(const OpenMS::Size idx_pep,
987  const OpenMS::Size idx_prot,
988  const OpenMS::Size len_pep,
989  const OpenMS::String& seq_prot,
991  {
992  if (enzyme_.isValidProduct(seq_prot, position, len_pep, true, true, xtandem_))
993  {
995  match.protein_index = idx_prot;
996  match.position = position;
997  match.AABefore = (position == 0) ? PeptideEvidence::N_TERMINAL_AA : seq_prot[position - 1];
998  match.AAAfter = (position + len_pep >= seq_prot.size()) ? PeptideEvidence::C_TERMINAL_AA : seq_prot[position + len_pep];
999  pep_to_prot[idx_pep].insert(match);
1000  ++filter_passed;
1001  }
1002  else
1003  {
1004  //std::cerr << "REJECTED Peptide " << seq_pep << " with hit to protein "
1005  // << seq_prot << " at position " << position << std::endl;
1006  ++filter_rejected;
1007  }
1008  }
1009 
1010  };
1011 
1012  inline 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
1013  {
1014  fuzzyAC.setProtein(prot);
1015  while (fuzzyAC.findNext(pattern))
1016  {
1017  const seqan::Peptide& tmp_pep = pep_DB[fuzzyAC.getHitDBIndex()];
1018  func_threads.addHit(fuzzyAC.getHitDBIndex(), idx_prot, length(tmp_pep), full_prot, fuzzyAC.getHitProteinPosition() + offset);
1019  }
1020 
1021  }
1022 
1023  void updateMembers_() override;
1024 
1026  bool prefix_;
1030 
1036 
1039 
1040  };
1041 }
1042 
void setMetaValue(const String &name, const DataValue &value)
Sets the DataValue corresponding to a name.
std::map< std::string, std::string > CaseInsensitiveToCaseSensitiveDecoy
Definition: PeptideIndexing.h:745
bool has(Byte byte) const
true if String contains the byte, false otherwise
char AABefore
the amino acid after the peptide in the protein
Definition: PeptideIndexing.h:904
A more convenient string class.
Definition: String.h:58
void after()
record data for the second timepoint
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:190
bool contains_decoys_
Definition: PeptideIndexing.h:746
#define LOG_INFO
Macro if a information, e.g. a status should be reported.
Definition: LogStream.h:456
OpenMS::Int position
the position of the peptide in the protein
Definition: PeptideIndexing.h:901
void setSequence(const String &sequence)
sets the protein sequence
Size< TNeedle >::Type position(const PatternAuxData< TNeedle > &dh)
Definition: AhoCorasickAmbiguous.h:561
bool operator==(const PeptideProteinMatchInformation &other) const
Definition: PeptideIndexing.h:930
StopWatch Class.
Definition: StopWatch.h:59
bool operator<(const PeptideProteinMatchInformation &other) const
Definition: PeptideIndexing.h:909
OpenMS::Size protein_index
index of the protein the peptide is contained in
Definition: PeptideIndexing.h:898
void setProtein(const String &protein_sequence)
Reset to new protein sequence. All previous data is forgotten.
Definition: AhoCorasickAmbiguous.h:1024
bool findDecoyString_(FASTAContainer< T > &proteins)
Definition: PeptideIndexing.h:749
A convenience class to report either absolute or delta (between two timepoints) RAM usage...
Definition: SysInfo.h:83
::seqan::Pattern< PeptideDB, ::seqan::FuzzyAC > FuzzyACPattern
Definition: AhoCorasickAmbiguous.h:974
Definition: PeptideIndexing.h:939
ptrdiff_t SignedSize
Signed Size type e.g. used as pointer difference.
Definition: Types.h:134
Refreshes the protein references for all peptide hits in a vector of PeptideIdentifications and adds ...
Definition: PeptideIndexing.h:124
bool keep_unreferenced_proteins_
Definition: PeptideIndexing.h:1033
no requirements on start / end
Definition: EnzymaticDigestion.h:70
String decoy_string_
Definition: PeptideIndexing.h:1025
String getEnzymeName() const
Returns the enzyme for the digestion.
Main OpenMS namespace.
Definition: FeatureDeconvolution.h:46
#define DEBUG_ONLY
Definition: AhoCorasickAmbiguous.h:46
std::map< std::string, std::pair< int, int > > DecoyStringToAffixCount
Definition: PeptideIndexing.h:744
Definition: PeptideIndexing.h:133
String missing_decoy_action_
Definition: PeptideIndexing.h:1027
void setSpecificity(Specificity spec)
Sets the specificity for the digestion (default is SPEC_FULL).
#define LOG_DEBUG
Macro for general debugging information.
Definition: LogStream.h:460
String toString() const
get a compact representation of the current time status.
String & remove(char what)
Remove all occurrences of the character what.
ExitCodes
Exit codes.
Definition: PeptideIndexing.h:130
bool prefix_
Definition: PeptideIndexing.h:1026
String enzyme_name_
Definition: PeptideIndexing.h:1028
#define LOG_ERROR
Macro to be used if non-fatal error are reported (processing continues)
Definition: LogStream.h:448
Size getHitDBIndex()
Get index of hit into peptide database of the pattern.
Definition: AhoCorasickAmbiguous.h:1047
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...
#define LOG_WARN
Macro if a warning, a piece of information which should be read by the user, should be logged...
Definition: LogStream.h:452
std::map< OpenMS::Size, std::set< PeptideProteinMatchInformation > > MapType
Definition: PeptideIndexing.h:942
bool findNext(const FuzzyACPattern &pattern)
Enumerate hits.
Definition: AhoCorasickAmbiguous.h:1037
Int aaa_max_
Definition: PeptideIndexing.h:1037
Int getHitProteinPosition()
Offset into protein sequence where hit was found.
Definition: AhoCorasickAmbiguous.h:1057
template parameter for vector-based FASTA access
Definition: FASTAContainer.h:76
Extended Aho-Corasick algorithm capable of matching ambiguous amino acids in the pattern (i...
Definition: AhoCorasickAmbiguous.h:970
Specificity getSpecificity() const
Returns the specificity for the digestion.
void merge(FoundProteinFunctor &other)
Definition: PeptideIndexing.h:965
static const char N_TERMINAL_AA
Definition: PeptideEvidence.h:60
bool write_protein_sequence_
Definition: PeptideIndexing.h:1031
Class for the enzymatic digestion of proteins.
Definition: ProteaseDigestion.h:60
Definition: PeptideIndexing.h:134
String & toLower()
Converts the string to lowercase.
bool isAmbiguous(AAcid c)
Definition: AhoCorasickAmbiguous.h:578
semi specific, i.e., one of the two cleavage sites must fulfill requirements
Definition: EnzymaticDigestion.h:69
Definition: PeptideIndexing.h:132
static String suffix(const String &this_s, size_t length)
Definition: StringUtils.h:269
MapType pep_to_prot
peptide index –> protein indices
Definition: PeptideIndexing.h:945
Definition: PeptideIndexing.h:135
void setDescription(const String &description)
sets the description of the protein
static String & toUpper(String &this_s)
Definition: StringUtils.h:732
bool IL_equivalent_
Definition: PeptideIndexing.h:1035
Representation of a peptide evidence.
Definition: PeptideEvidence.h:50
bool has(const Key &key) const
Test whether the map contains the given key.
Definition: Map.h:108
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
void setEnzyme(const String &name)
Sets the enzyme for the digestion (by name)
Definition: PeptideIndexing.h:137
void setAccession(const String &accession)
sets the accession of the protein
Representation of a protein hit.
Definition: ProteinHit.h:57
String enzyme_specificity_
Definition: PeptideIndexing.h:1029
FoundProteinFunctor(const ProteaseDigestion &enzyme, bool xtandem)
Definition: PeptideIndexing.h:960
static Specificity getSpecificityByName(const String &name)
String delta(const String &event="delta")
bool hasPrefix(const String &string) const
true if String begins with string, false otherwise
size_t Size
Size type e.g. used as variable which can hold result of size()
Definition: Types.h:127
String & substitute(char from, char to)
Replaces all occurrences of the character from by the character to.
OpenMS::Size filter_rejected
number of rejected hits (not passing addHit())
Definition: PeptideIndexing.h:951
bool write_protein_description_
Definition: PeptideIndexing.h:1032
Base class for all classes that want to report their progress.
Definition: ProgressLogger.h:54
Int mm_max_
Definition: PeptideIndexing.h:1038
String sequence
Definition: FASTAFile.h:80
static String prefix(const String &this_s, size_t length)
Definition: StringUtils.h:260
FASTA entry type (identifier, description and sequence)
Definition: FASTAFile.h:76
A base class for all classes handling default parameters.
Definition: DefaultParamHandler.h:91
Definition: PeptideIndexing.h:136
bool allow_unmatched_
Definition: PeptideIndexing.h:1034
static const char C_TERMINAL_AA
Definition: PeptideEvidence.h:61
bool xtandem_
are we checking xtandem cleavage rules?
Definition: PeptideIndexing.h:957
String< AAcid, Alloc< void > > AAString
Definition: AhoCorasickAmbiguous.h:206
int Int
Signed integer type.
Definition: Types.h:102
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...
char AAAfter
the amino acid before the peptide in the protein
Definition: PeptideIndexing.h:907
Map class based on the STL map (containing several convenience functions)
Definition: Map.h:50
FASTAContainer<TFI_Vector> simply takes an existing vector of FASTAEntries and provides the same inte...
Definition: FASTAContainer.h:237
::seqan::StringSet<::seqan::AAString > PeptideDB
Definition: AhoCorasickAmbiguous.h:973
String description
Definition: FASTAFile.h:79
bool hasSuffix(const String &string) const
true if String ends with string, false otherwise
OpenMS::Size filter_passed
number of accepted hits (passing addHit() constraints)
Definition: PeptideIndexing.h:948
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:148
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:986
double getClockTime() const
ProteaseDigestion enzyme_
Definition: PeptideIndexing.h:954
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:1012