bpp-seq-omics  2.4.1
FullGapFilterMafIterator.cpp
Go to the documentation of this file.
1 //
2 // File: FullGapFilterMafIterator.cpp
3 // Authors: Julien Dutheil
4 // Created: Tue Sep 07 2010
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 
41 
42 //From bpp-seq
44 #include <Bpp/Seq/SiteTools.h>
45 
46 using namespace bpp;
47 
48 //From the STL:
49 #include <string>
50 #include <numeric>
51 
52 using namespace std;
53 
55 {
56  MafBlock* block = iterator_->nextBlock();
57  if (!block) return 0;
58 
59  //We create a copy of the ingroup alignement for better efficiency:
61  for (size_t i = 0; i < species_.size(); ++i) {
62  if (block->hasSequenceForSpecies(species_[i])) {
63  vsc.addSequence(block->getSequenceForSpecies(species_[i]));
64  }
65  }
66  if (vsc.getNumberOfSequences() == 0) return block; //Block ignored as it does not contain any of the focus species.
67 
68  //Now check the positions that are only made of gaps:
69  if (verbose_) {
70  ApplicationTools::message->endLine();
71  ApplicationTools::displayTask("Cleaning block for gap sites", true);
72  }
73  size_t n = block->getNumberOfSites();
74  vector <size_t> start;
75  vector <unsigned int> count;
76  bool test = false;
77  for (size_t i = 0; i < n; ++i) {
78  const Site* site = &vsc.getSite(i);
79  if (SiteTools::isGapOnly(*site)) {
80  if (test) {
81  count[count.size() - 1]++;
82  } else {
83  start.push_back(i);
84  count.push_back(1);
85  test = true;
86  }
87  } else {
88  test = false;
89  }
90  }
91  //Now remove blocks:
92  size_t totalRemoved = 0;
93  for(size_t i = start.size(); i > 0; --i) {
94  if (verbose_)
95  ApplicationTools::displayGauge(start.size() - i, start.size() - 1, '=');
96  block->getAlignment().deleteSites(start[i - 1], count[i - 1]);
97  totalRemoved += count[i - 1];
98  }
99  if (verbose_)
101 
102  //Correct coordinates:
103  if (totalRemoved > 0) {
104  for (size_t i = 0; i < block->getNumberOfSequences(); ++i) {
105  const MafSequence* seq = &block->getSequence(i);
106  if (!VectorTools::contains(species_, seq->getSpecies())) {
108  }
109  }
110  }
111  if (logstream_) {
112  (*logstream_ << "FULL GAP CLEANER: " << totalRemoved << " positions have been removed.").endLine();
113  }
114  return block;
115 }
116 
void deleteSites(size_t siteIndex, size_t length)
static const DNA DNA_ALPHABET
static void displayTask(const std::string &text, bool eof=false)
static std::shared_ptr< OutputStream > message
static void displayTaskDone()
static void displayGauge(size_t iter, size_t total, char symbol='>', const std::string &mes="")
A synteny block data structure, the basic unit of a MAF alignement file.
Definition: MafBlock.h:57
size_t getNumberOfSequences() const
Definition: MafBlock.h:111
size_t getNumberOfSites() const
Definition: MafBlock.h:113
const MafSequence & getSequence(const std::string &name) const
Definition: MafBlock.h:121
AlignedSequenceContainer & getAlignment()
Definition: MafBlock.h:108
void removeCoordinatesFromSequence(size_t i)
Definition: MafBlock.h:170
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 & getSpecies() const
Definition: MafSequence.h:162
static bool isGapOnly(const IntCoreSymbolList &site)
void addSequence(const Sequence &sequence, bool checkName=true)
const Site & getSite(size_t siteIndex) const
size_t getNumberOfSequences() const
static bool contains(const std::vector< T > &vec, T el)
std::size_t count(const std::string &s, const std::string &pattern)