bpp-seq-omics  2.4.1
VcfOutputMafIterator.cpp
Go to the documentation of this file.
1 //
2 // File: VcfOutputMafIterator.cpp
3 // Authors: Julien Dutheil
4 // Created: Tue Jan 05 2013
5 //
6 
7 /*
8 Copyright or © or Copr. Bio++ Development Team, (2010)
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 "VcfOutputMafIterator.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 VcfOutputMafIterator::writeHeader_(std::ostream& out) const
59 {
60  time_t t = time(0); // get current time
61  struct tm* ct = localtime(&t);
62  out << "##fileformat=VCFv4.0" << endl;
63  out << "##fileDate=" << (ct->tm_year + 1900) << (ct->tm_mon + 1) << ct->tm_mday << endl;
64  out << "##source=Bio++" << endl;
65  out << "##FILTER=<ID=PASS,Description=\"All filters passed\">" << endl;
66  out << "##FILTER=<ID=gap,Description=\"At least one sequence contains a gap\">" << endl;
67  out << "##FILTER=<ID=unk,Description=\"At least one sequence contains an unresolved character\">" << endl;
68  if (genotypes_.size() > 0)
69  out << "##FORMAT=<ID=GT,Number=1,Type=String,Description=\"Genotype\">" << endl;
70  out << "##INFO=<ID=AC,Number=A,Type=Integer,Description=\"Total number of alternate alleles in called genotypes\">" << endl;
71  //There are more options in the header that we may want to support...
72 
73  //Now write the header line:
74  out << "#CHROM\tPOS\tID\tREF\tALT\tQUAL\tFILTER\tINFO";
75  if (genotypes_.size() > 0) {
76  out << "\tFORMAT";
77  for (size_t i = 0; i < genotypes_.size(); ++i) {
78  out << "\t" << genotypes_[i][0];
79  if (genotypes_[i].size() > 1) {
80  //Polyploid case
81  for (size_t j = 1; j < genotypes_[i].size(); ++j) {
82  out << "-" << genotypes_[i][j];
83  }
84  }
85  }
86  }
87  out << endl;
88 }
89 
90 void VcfOutputMafIterator::writeBlock_(std::ostream& out, const MafBlock& block) const
91 {
92  const MafSequence& refSeq = block.getSequenceForSpecies(refSpecies_);
93  string chr = refSeq.getChromosome();
94  SequenceWalker walker(refSeq);
95  size_t offset = refSeq.start();
96  int gap = refSeq.getAlphabet()->getGapCharacterCode();
97  map<int, string> chars;
98  for (int i = 0; i < static_cast<int>(AlphabetTools::DNA_ALPHABET.getNumberOfTypes()); ++i)
99  chars[i] = AlphabetTools::DNA_ALPHABET.intToChar(i);
100  VectorSiteContainer sites(block.getAlignment());
101  //Where to store genotype information, if any:
102  vector<int> gt(genotypes_.size());
103  //Now we look all sites for SNPs:
104  for (size_t i = 0; i < sites.getNumberOfSites(); i++) {
105  if (refSeq[i] == gap)
106  continue;
107  string filter = "";
108  if (SiteTools::hasGap(sites.getSite(i))) {
109  filter = "gap";
110  }
111  if (SymbolListTools::hasUnresolved(sites.getSite(i))) {
112  if (filter != "")
113  filter += ";";
114  filter += "unk";
115  }
116  if (filter == "")
117  filter = "PASS";
118 
119  map<int, size_t> counts;
120  SiteTools::getCounts(sites.getSite(i), counts);
121  int ref = refSeq[i];
122  string alt = "";
123  string ac = "";
124 
125  map<int, int> snps;
126  int c = 0;
127  for (int x = 0; x < 4; ++x) {
128  if (x != ref) {
129  size_t f = counts[x];
130  if (f > 0) {
131  if (alt != "") {
132  alt += ",";
133  ac += ",";
134  }
135  alt += chars[x];
136  ac += TextTools::toString(f);
137  snps[x] = ++c;
138  }
139  } else {
140  snps[x] = 0;
141  }
142  }
143  if (ac == "" && outputAll_) {
144  ac = TextTools::toString(counts[ref]);
145  }
146  if (ac != "") {
147  out << chr << "\t" << (offset + walker.getSequencePosition(i) + 1) << "\t.\t" << chars[refSeq[i]] << "\t" << alt << "\t.\t" << filter << "\tAC=" << ac;
148  //Write genotpyes:
149  if (genotypes_.size() > 0) {
150  out << "\tGT";
151  for (size_t g = 0; g < genotypes_.size(); ++g) {
152  string geno = "";
153  for (auto x: genotypes_[g]) {
154  if (geno != "") geno += "|"; //Polyploid
155  vector<const MafSequence*> sequences = block.getSequencesForSpecies(x);
156  if (sequences.size() == 0)
157  geno += (generateDiploids_ ? ".|." : ".");
158  else if (sequences.size() > 1)
159  throw Exception("VcfOutputMafIterator::writeBlock(). Duplicated sequence for species '" + x + "'.");
160  else {
161  int state = (*sequences[0])[i];
163  geno += (generateDiploids_ ? ".|." : ".");
164  else {
165  geno += TextTools::toString(snps[state]);
166  if (generateDiploids_)
167  geno += "|" + TextTools::toString(snps[state]);
168  }
169  }
170  }
171  out << "\t" << geno;
172  }
173  }
174  out << endl;
175  }
176  }
177 }
178 
static const DNA DNA_ALPHABET
virtual int getGapCharacterCode() 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
std::vector< const MafSequence * > getSequencesForSpecies(const std::string &species) const
Definition: MafBlock.h:149
AlignedSequenceContainer & getAlignment()
Definition: MafBlock.h:108
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
bool isUnresolved(int state) const
size_t getSequencePosition(size_t alnPos)
static bool hasGap(const IntCoreSymbolList &site)
static void getCounts(const IntCoreSymbolList &list, std::map< int, size_t > &counts)
void writeBlock_(std::ostream &out, const MafBlock &block) const
void writeHeader_(std::ostream &out) const
const Site & getSite(size_t siteIndex) const
size_t getNumberOfSites() const
std::string toString(T t)