bpp-seq-omics  2.4.1
MsmcOutputMafIterator.cpp
Go to the documentation of this file.
1 //
2 // File: MsmcOutputMafIterator.cpp
3 // Authors: Julien Dutheil
4 // Created: Tue Jan 06 2015
5 //
6 
7 /*
8 Copyright or © or Copr. Bio++ Development Team, (2015)
9 
10 This software is a computer program whose purpose is to provide classes
11 for sequences analysis.
12 
13 This software is governed by the CeCILL license under French law and
14 abiding by the rules of distribution of free software. You can use,
15 modify and/ or redistribute the software under the terms of the CeCILL
16 license as circulated by CEA, CNRS and INRIA at the following URL
17 "http://www.cecill.info".
18 
19 As a counterpart to the access to the source code and rights to copy,
20 modify and redistribute granted by the license, users are provided only
21 with a limited warranty and the software's author, the holder of the
22 economic rights, and the successive licensors have only limited
23 liability.
24 
25 In this respect, the user's attention is drawn to the risks associated
26 with loading, using, modifying and/or developing or reproducing the
27 software by the user in light of its specific status of free software,
28 that may mean that it is complicated to manipulate, and that also
29 therefore means that it is reserved for developers and experienced
30 professionals having in-depth computer knowledge. Users are therefore
31 encouraged to load and test the software's suitability as regards their
32 requirements in conditions enabling the security of their systems and/or
33 data to be ensured and, more generally, to use and operate it in the
34 same conditions as regards security.
35 
36 The fact that you are presently reading this means that you have had
37 knowledge of the CeCILL license and that you accept its terms.
38 */
39 
40 #include "MsmcOutputMafIterator.h"
41 
42 //From bpp-seq:
46 #include <Bpp/Seq/SiteTools.h>
47 #include <Bpp/Seq/SequenceWalker.h>
48 
49 using namespace bpp;
50 
51 //From the STL:
52 #include <string>
53 #include <numeric>
54 #include <ctime>
55 
56 using namespace std;
57 
58 void MsmcOutputMafIterator::writeBlock_(std::ostream& out, const MafBlock& block)
59 {
60  //Preliminary stuff...
61 
63  for (size_t i = 0; i < species_.size(); ++i) {
64  if (block.hasSequenceForSpecies(species_[i])) {
65  sites.addSequence(block.getSequenceForSpecies(species_[i]));
66  //Note: in case of duplicates, this takes the first sequence.
67  } else {
68  //Block with missing species are ignored.
69  return;
70  }
71  }
72  //Get the reference species for coordinates:
73  if (! block.hasSequenceForSpecies(refSpecies_))
74  return;
75  const MafSequence& refSeq = block.getSequenceForSpecies(refSpecies_);
76  string chr = refSeq.getChromosome();
77  if (chr != currentChr_) {
78  currentChr_ = chr;
79  nbOfCalledSites_ = 0; //Reset count of called sites.
80  lastPosition_ = 0;
81  } else {
82  //Check that block are ordered according to reference sequence:
83  if (refSeq.start() < lastPosition_)
84  throw Exception("MsmcOutputMafIterator: blocks are not projected according to reference sequence: " + refSeq.getDescription() + "<!>" + TextTools::toString(lastPosition_) + ".");
85  lastPosition_ = refSeq.stop();
86  }
87 
88  SequenceWalker walker(refSeq);
89  size_t offset = refSeq.start();
90  int gap = refSeq.getAlphabet()->getGapCharacterCode();
91 
92  //Now we shall scan all sites for SNPs:
93  for (size_t i = 0; i < sites.getNumberOfSites(); i++) {
94  if (refSeq[i] == gap)
95  continue;
96 
97  //We call SNPs only at position without gap or unresolved characters:
98  if (SiteTools::isComplete(sites.getSite(i))) {
99  nbOfCalledSites_++;
100 
101  if (!SiteTools::isConstant(sites.getSite(i))) {
102  string pos = "NA";
103  if (refSeq[i] != gap) {
104  pos = TextTools::toString(offset + walker.getSequencePosition(i) + 1);
105  }
106  out << chr << "\t" << pos << "\t" << nbOfCalledSites_ << "\t" << sites.getSite(i).toString() << endl;
107  //Reset number of called sites
108  nbOfCalledSites_ = 0;
109  }
110  }
111  }
112 }
113 
static const DNA DNA_ALPHABET
virtual int getGapCharacterCode() const=0
virtual std::string toString() const=0
virtual const Alphabet * getAlphabet() const=0
A synteny block data structure, the basic unit of a MAF alignement file.
Definition: MafBlock.h:57
bool hasSequenceForSpecies(const std::string &species) const
Definition: MafBlock.h:129
const MafSequence & getSequenceForSpecies(const std::string &species) const
Definition: MafBlock.h:139
A sequence class which is used to store data from MAF files.
Definition: MafSequence.h:64
std::string getDescription() const
Definition: MafSequence.h:188
const std::string & getChromosome() const
Definition: MafSequence.h:164
size_t start() const
Definition: MafSequence.h:116
size_t stop() const
Definition: MafSequence.h:121
void writeBlock_(std::ostream &out, const MafBlock &block)
size_t getSequencePosition(size_t alnPos)
static bool isComplete(const IntCoreSymbolList &site)
static bool isConstant(const IntCoreSymbolList &site, bool ignoreUnknown=false, bool unresolvedRaisesException=true)
void addSequence(const Sequence &sequence, bool checkName=true)
const Site & getSite(size_t siteIndex) const
size_t getNumberOfSites() const
std::string toString(T t)