bpp-seq-omics  2.4.1
SequenceLDhotOutputMafIterator.cpp
Go to the documentation of this file.
1 //
2 // File: OutputSequenceLDHotMafIterator.cpp
3 // Authors: Julien Dutheil
4 // Created: Thr May 12 2016
5 //
6 
7 /*
8 Copyright or © or Copr. Bio++ Development Team, (2016)
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 
41 
42 //From bpp-seq:
45 #include <Bpp/Seq/SiteTools.h>
46 
47 using namespace bpp;
48 
49 //From the STL:
50 #include <string>
51 #include <numeric>
52 
53 using namespace std;
54 
56 {
57  MafBlock* block = iterator_->nextBlock();
58  if (block) {
59  string chr = "ChrNA";
60  string start = "StartNA";
61  string stop = "StopNA";
62  if (block->hasSequenceForSpecies(refSpecies_)) {
63  const MafSequence& refseq = block->getSequenceForSpecies(refSpecies_);
64  chr = refseq.getChromosome();
65  start = TextTools::toString(refseq.start());
66  stop = TextTools::toString(refseq.stop());
67  }
68  string file = file_;
69  TextTools::replaceAll(file, "%i", TextTools::toString(++currentBlockIndex_));
70  TextTools::replaceAll(file, "%c", chr);
71  TextTools::replaceAll(file, "%b", start);
72  TextTools::replaceAll(file, "%e", stop);
73  std::ofstream output(file.c_str(), ios::out);
74  writeBlock(output, *block);
75  }
76  return block;
77 }
78 
79 void SequenceLDhotOutputMafIterator::writeBlock(std::ostream& out, const MafBlock& block) const {
80  //First get alignment:
81  const SiteContainer& aln = block.getAlignment();
82  unique_ptr<VectorSiteContainer> variableSites(new VectorSiteContainer(aln.getSequencesNames(), &AlphabetTools::DNA_ALPHABET));
83 
84  //We first preparse the data:
85  //We assume all sequences are distinct:
86  size_t nbDistinct = aln.getNumberOfSequences();
87  size_t nbGenes = aln.getNumberOfSequences();
88  size_t nbLoci = 0;
89 
90  string positions = "";
91  for (size_t i = 0; i < aln.getNumberOfSites(); ++i) {
92  const Site& s = aln.getSite(i);
93  if (completeOnly_ && !SiteTools::isComplete(s)) {
94  continue;
95  }
96  unsigned int count = 0;
97  int x = -1;
98  for (size_t j = 0; j < s.size() && count < 2; ++j) {
99  if (!AlphabetTools::DNA_ALPHABET.isGap(s[j]) && !AlphabetTools::DNA_ALPHABET.isUnresolved(s[j])) {
100  if (count == 0) {
101  //First state found
102  count++;
103  //We record the state
104  x = s[j];
105  } else {
106  if (s[j] != x) {
107  //New state found
108  count++;
109  }
110  //Otherwise, same state as before.
111  }
112  }
113  }
114  if (count == 2) {
115  //At least two alleles (non-gap, non-unresolved) found in this position, so we record it
116  positions += " " + TextTools::toString(i + 1);
117  nbLoci++;
118  variableSites->addSite(aln.getSite(i));
119  }
120  }
121 
122  //Write header:
123  out << "Distinct = " << nbDistinct << endl;
124  out << "Genes = " << nbGenes << endl;
125  out << "Loci = " << nbLoci << endl;
126  out << "K = -4 %4-allele model with Haplotype Alleles specified by A,C,G,T" << endl;
127 
128  out << "Positions of loci:" << endl;
129  out << positions << endl;
130 
131  out << "Haplotypes" << endl;
132 
133  for (size_t i = 0; i < aln.getNumberOfSequences(); ++i) {
134  out << variableSites->getSequence(i).toString() << " 1" << endl;
135  }
136 
137  out << "#" << endl;
138 }
139 
virtual size_t getNumberOfSites() const=0
static const DNA DNA_ALPHABET
virtual size_t size() const=0
A synteny block data structure, the basic unit of a MAF alignement file.
Definition: MafBlock.h:57
AlignedSequenceContainer & getAlignment()
Definition: MafBlock.h:108
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
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
virtual size_t getNumberOfSequences() const=0
virtual std::vector< std::string > getSequencesNames() const=0
void writeBlock(std::ostream &out, const MafBlock &block) const
virtual const Site & getSite(size_t siteIndex) const=0
static bool isComplete(const IntCoreSymbolList &site)
void replaceAll(std::string &target, const std::string &query, const std::string &replacement)
std::string toString(T t)
std::size_t count(const std::string &s, const std::string &pattern)