bpp-phyl  2.2.0
AbstractTreeLikelihood.h
Go to the documentation of this file.
1 //
2 // File: AbstractTreeLikelihood.h
3 // Created by: Julien Dutheil
4 // Created on: Fri Oct 17 17:57:21 2003
5 //
6 
7 /*
8 Copyright or © or Copr. Bio++ Development Team, (November 16, 2004)
9 
10 This software is a computer program whose purpose is to provide classes
11 for phylogenetic data 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 #ifndef _ABSTRACTTREELIKELIHOOD_H_
41 #define _ABSTRACTTREELIKELIHOOD_H_
42 
43 #include "TreeLikelihood.h"
44 #include "../Tree.h"
45 #include "../TreeTemplate.h"
46 
47 #include <Bpp/Numeric/AbstractParametrizable.h>
48 
49 //From SeqLib:
50 #include <Bpp/Seq/Container/SiteContainer.h>
51 
52 namespace bpp
53 {
54 
70  public virtual TreeLikelihood,
71  public AbstractParametrizable
72 {
73  public:
80  public BranchIterator
81  {
82  private:
83  std::vector<int> nodesId_;
84  size_t index_;
85 
86  public:
87  SimpleBranchIterator(const std::vector<int>& nodesId) :
88  nodesId_(nodesId), index_(0) {}
89 
90  public:
91  int next() throw (Exception)
92  {
93  if (!hasNext())
94  throw Exception("AbstractTreeLikelihood::SimpleBranchIterator::next(). No more branch in the set.");
95  return nodesId_[index_++];
96  }
97 
98  bool hasNext() const { return index_ < nodesId_.size(); }
99 
100  };
101 
110  public SiteIterator
111  {
112  private:
113  size_t maxIndex_;
114  size_t index_;
115  size_t offset_;
116 
117  public:
118  SimpleSiteIterator(size_t nbSites, size_t offset = 0) :
119  maxIndex_(nbSites), index_(0), offset_(offset) {}
120 
121  public:
122  size_t next() throw (Exception)
123  {
124  if (!hasNext())
125  throw Exception("AbstractTreeLikelihood::SimpleSiteIterator::next(). No more site in the set.");
126  return offset_ + index_++;
127  }
128 
129  bool hasNext() const { return index_ < maxIndex_; }
130 
131  };
132 
140  {
141  private:
143  size_t nbSites_;
144 
145  public:
147  model_(model), nbSites_(nbSites) {}
148 
150  model_(bmd.model_),
151  nbSites_(bmd.nbSites_)
152  {}
153 
155  {
156  model_ = bmd.model_;
157  nbSites_ = bmd.nbSites_;
158  return *this;
159  }
160 
161  public:
162  const SubstitutionModel* getModel() const { return model_; }
163 
165  };
166 
169  {
170  private:
172  size_t index_;
173 
174  public:
176  branchModelDescription_(model, nbSites), index_(0) {}
177 
178  public:
180  {
181  if (!hasNext())
182  throw Exception("AbstractHomogeneousTreeLikelihood::ConstHomogeneousBranchModelIterator::next(). No more branch in the set.");
183  index_++;
184  return &branchModelDescription_;
185  }
186 
187  bool hasNext() const { return index_ == 0; }
188  };
189 
192  {
193  private:
195  std::vector<int> nodesId_;
196 
197  public:
198  ConstNoPartitionSiteModelDescription(const SubstitutionModel* model, const std::vector<int> nodesId) :
199  model_(model), nodesId_(nodesId) {}
200 
202  model_(smd.model_),
203  nodesId_(smd.nodesId_)
204  {}
205 
207  {
208  model_ = smd.model_;
209  nodesId_ = smd.nodesId_;
210  return *this;
211  }
212 
213  public:
214  const SubstitutionModel* getModel() const { return model_; }
215 
217  };
218 
224  protected:
225  const SiteContainer* data_;
230 
231  public:
233  AbstractParametrizable(""),
234  data_(0),
235  tree_(0),
238  initialized_(false) {}
239 
241  AbstractParametrizable(lik),
242  data_(0),
243  tree_(0),
247  {
248  if (lik.data_) data_ = dynamic_cast<SiteContainer*>(lik.data_->clone());
249  if (lik.tree_) tree_ = lik.tree_->clone();
250  }
251 
253  {
254  AbstractParametrizable::operator=(lik);
255  if (data_) delete data_;
256  if (lik.data_) data_ = dynamic_cast<SiteContainer*>(lik.data_->clone());
257  else data_ = 0;
258  if (tree_) delete tree_;
259  if (lik.tree_) tree_ = lik.tree_->clone();
260  else tree_ = 0;
264  return *this;
265  }
266 
273  {
274  if (data_) delete data_;
275  if (tree_) delete tree_;
276  }
277 
278  public:
284  const SiteContainer* getData() const { return data_; }
285  const Alphabet* getAlphabet() const { return data_->getAlphabet(); }
286  Vdouble getLikelihoodForEachSite() const;
287  Vdouble getLogLikelihoodForEachSite() const;
288  VVdouble getLikelihoodForEachSiteForEachState() const;
290  size_t getNumberOfSites() const { return data_->getNumberOfSites(); }
291  const Tree& getTree() const { return *tree_; }
297  bool isInitialized() const { return initialized_; }
298  void initialize() throw (Exception) { initialized_ = true; }
301 // protected:
302 //
303 // /**
304 // * @brief Recompute pxy_, dpxy_ and d2pxy_ arrays, and derivatives if needed.
305 // *
306 // * This method is called when some parameter has changed.
307 // *
308 // * @param params The parameters that changed.
309 // */
310 // virtual void fireParameterChanged(const ParameterList & params) = 0;
311 
312 };
313 
314 } //end of namespace bpp.
315 
316 #endif //_ABSTRACTTREELIKELIHOOD_H_
317 
A pair of SubstitutionModel / SiteIterator.
ConstNoPartitionSiteModelDescription & operator=(const ConstNoPartitionSiteModelDescription &smd)
ConstNoPartitionSiteModelDescription(const ConstNoPartitionSiteModelDescription &smd)
ConstNoPartitionSiteModelDescription(const SubstitutionModel *model, const std::vector< int > nodesId)
Interface for all substitution models.
virtual ~AbstractTreeLikelihood()
Abstract class destructor.
An iterator over a set of sites, speicfied by their position.
Iterates through all models used for all sites on a given branch.
VVdouble getLogLikelihoodForEachSiteForEachState() const
Get the logarithm of the likelihood for each site and for each state.
const SiteContainer * getData() const
Get the dataset for which the likelihood must be evaluated.
The TreeLikelihood interface.
void enableDerivatives(bool yn)
Tell if derivatives must be computed.
The phylogenetic tree class.
Vdouble getLogLikelihoodForEachSite() const
Get the logarithm of the likelihood for each site.
Interface for phylogenetic tree objects.
Definition: Tree.h:148
ConstNoPartitionBranchModelIterator(const SubstitutionModel *model, size_t nbSites)
ConstNoPartitionBranchModelDescription(const SubstitutionModel *model, size_t nbSites)
AbstractTreeLikelihood & operator=(const AbstractTreeLikelihood &lik)
ConstNoPartitionBranchModelDescription(const ConstNoPartitionBranchModelDescription &bmd)
VVdouble getLikelihoodForEachSiteForEachState() const
Get the likelihood for each site and for each state.
Partial implementation of the TreeLikelihood interface.
AbstractTreeLikelihood(const AbstractTreeLikelihood &lik)
Vdouble getLikelihoodForEachSite() const
Get the likelihood for each site.
size_t getNumberOfSites() const
Get the number of sites in the dataset.
const Alphabet * getAlphabet() const
Get the alphabet associated to the dataset.
const Tree & getTree() const
Get the tree (topology and branch lengths).
An iterator over a set of branches, specified by their node ids.
SimpleBranchIterator(const std::vector< int > &nodesId)
A pair of SubstitutionModel / BranchIterator.
ConstNoPartitionBranchModelDescription & operator=(const ConstNoPartitionBranchModelDescription &bmd)
void initialize()
Init the likelihood object.