bpp-phyl  2.2.0
AbstractNonHomogeneousTreeLikelihood.h
Go to the documentation of this file.
1 //
2 // File: AbstractNonHomogeneousTreeLikelihood.h
3 // Created by: Julien Dutheil
4 // Created on: Tue Oct 09 16:07 2007
5 // From file: AbstractHomogeneousTreeLikelihood.h
6 //
7 
8 /*
9 Copyright or © or Copr. Bio++ Development Team, (November 16, 2004)
10 
11 This software is a computer program whose purpose is to provide classes
12 for phylogenetic data analysis.
13 
14 This software is governed by the CeCILL license under French law and
15 abiding by the rules of distribution of free software. You can use,
16 modify and/ or redistribute the software under the terms of the CeCILL
17 license as circulated by CEA, CNRS and INRIA at the following URL
18 "http://www.cecill.info".
19 
20 As a counterpart to the access to the source code and rights to copy,
21 modify and redistribute granted by the license, users are provided only
22 with a limited warranty and the software's author, the holder of the
23 economic rights, and the successive licensors have only limited
24 liability.
25 
26 In this respect, the user's attention is drawn to the risks associated
27 with loading, using, modifying and/or developing or reproducing the
28 software by the user in light of its specific status of free software,
29 that may mean that it is complicated to manipulate, and that also
30 therefore means that it is reserved for developers and experienced
31 professionals having in-depth computer knowledge. Users are therefore
32 encouraged to load and test the software's suitability as regards their
33 requirements in conditions enabling the security of their systems and/or
34 data to be ensured and, more generally, to use and operate it in the
35 same conditions as regards security.
36 
37 The fact that you are presently reading this means that you have had
38 knowledge of the CeCILL license and that you accept its terms.
39 */
40 
41 #ifndef _ABSTRACTNONHOMOGENEOUSTREELIKELIHOOD_H_
42 #define _ABSTRACTNONHOMOGENEOUSTREELIKELIHOOD_H_
43 
46 
47 //From the STL:
48 #include <memory>
49 
50 namespace bpp
51 {
52 
59  public virtual NonHomogeneousTreeLikelihood,
61 {
62  public:
63 
66  {
67  private:
68  std::vector<ConstNoPartitionSiteModelDescription> siteModelDescriptions_;
69  size_t index_;
70  size_t nbModels_;
71 
72  public:
74  siteModelDescriptions_(), index_(0), nbModels_(modelSet->getNumberOfModels())
75  {
76  for (size_t i = 0; i < nbModels_; ++i)
78  }
79 
80  public:
81  ConstSiteModelDescription* next() throw (Exception)
82  {
83  if (!hasNext())
84  throw Exception("AbstractNonHomogeneousTreeLikelihood::ConstHomogeneousSiteModelIterator::next(). No more site in the set.");
85  return &siteModelDescriptions_[index_++];
86  }
87 
88  bool hasNext() const { return index_ < nbModels_; }
89  };
90 
91  protected:
93  ParameterList brLenParameters_;
94 
95  mutable std::map<int, VVVdouble> pxy_;
96 
97  mutable std::map<int, VVVdouble> dpxy_;
98 
99  mutable std::map<int, VVVdouble> d2pxy_;
100 
101  std::vector<double> rootFreqs_;
102 
109  std::vector<Node*> nodes_;
110 
114  mutable std::map<int, const Node*> idToNode_;
115 
116  //some values we'll need:
117  size_t nbSites_, //the number of sites in the container
118  nbDistinctSites_, //the number of distinct sites
119  nbClasses_, //the number of rate classes
120  nbStates_, //the number of states in the alphabet
121  nbNodes_; //the number of nodes in the tree
122 
123  bool verbose_;
124 
127  std::auto_ptr<Constraint> brLenConstraint_;
128 
131 
132  public:
134  const Tree& tree,
135  SubstitutionModelSet* modelSet,
136  DiscreteDistribution* rDist,
137  bool verbose = true,
138  bool reparametrizeRoot = true)
139  throw (Exception);
140 
147 
154 
156 
157  private:
158 
162  void init_(
163  const Tree& tree,
164  SubstitutionModelSet* modelSet,
165  DiscreteDistribution* rDist,
166  bool verbose) throw (Exception);
167 
168  public:
169 
177  size_t getNumberOfStates() const { return modelSet_->getNumberOfStates(); }
178 
179  const std::vector<int>& getAlphabetStates() const { return modelSet_->getAlphabetStates(); }
180 
181  int getAlphabetStateAsInt(size_t i) const { return modelSet_->getAlphabetStateAsInt(i); }
182 
183  std::string getAlphabetStateAsChar(size_t i) const { return modelSet_->getAlphabetStateAsChar(i); }
184 
185  void initialize() throw(Exception);
186 
187  ParameterList getBranchLengthsParameters() const;
188 
189  ParameterList getSubstitutionModelParameters() const;
190 
191  ParameterList getRateDistributionParameters() const
192  {
194  }
195 
197  {
198  return modelSet_->getModelForNode(nodeId);
199  }
200 
202  {
203  return modelSet_->getModelForNode(nodeId);
204  }
205 
206  const std::vector<double>& getRootFrequencies(size_t siteIndex) const { return rootFreqs_; }
207 
208  VVVdouble getTransitionProbabilitiesPerRateClass(int nodeId, size_t siteIndex) const { return pxy_[nodeId]; }
209 
211  {
213  }
214 
216  {
218  }
219 
230 
232 
233  void setSubstitutionModelSet(SubstitutionModelSet* modelSet) throw (Exception);
234 
235  ParameterList getRootFrequenciesParameters() const
236  {
238  }
241  public: //Specific methods:
242 
247  virtual void initParameters();
248 
254  virtual void applyParameters() throw (Exception);
255 
256  virtual void initBranchLengthsParameters();
257 
258  virtual void setMinimumBranchLength(double minimum) throw (Exception)
259  {
260  if (minimum > maximumBrLen_)
261  throw Exception("AbstractNonHomogeneousTreeLikelihood::setMinimumBranchLength. Minimum branch length sould be lower than the maximum one: " + TextTools::toString(maximumBrLen_));
262  minimumBrLen_ = minimum;
263  if (brLenConstraint_.get()) brLenConstraint_.release();
264  brLenConstraint_.reset(new IntervalConstraint(minimumBrLen_, maximumBrLen_, true, true));
266  }
267 
268  virtual void setMaximumBranchLength(double maximum) throw (Exception)
269  {
270  if (maximum < minimumBrLen_)
271  throw Exception("AbstractNonHomogeneousTreeLikelihood::setMaximumBranchLength. Maximum branch length sould be higher than the minimum one: " + TextTools::toString(minimumBrLen_));
272  maximumBrLen_ = maximum;
273  if (brLenConstraint_.get()) brLenConstraint_.release();
274  brLenConstraint_.reset(new IntervalConstraint(minimumBrLen_, maximumBrLen_, true, true));
276  }
277 
278  virtual double getMinimumBranchLength() const { return minimumBrLen_; }
279  virtual double getMaximumBranchLength() const { return maximumBrLen_; }
280 
281 
282  protected:
286  virtual void computeAllTransitionProbabilities();
290  virtual void computeTransitionProbabilitiesForNode(const Node * node);
291 
292 };
293 
294 } //end of namespace bpp.
295 
296 #endif //_ABSTRACTNONHOMOGENEOUSTREELIKELIHOOD_H_
297 
Substitution models manager for non-homogeneous / non-reversible models of evolution.
AbstractNonHomogeneousTreeLikelihood(const Tree &tree, SubstitutionModelSet *modelSet, DiscreteDistribution *rDist, bool verbose=true, bool reparametrizeRoot=true)
ParameterList getSubstitutionModelParameters() const
Get the parameters associated to substitution model(s).
Interface for all substitution models.
Partial implementation for branch non-homogeneous models of the TreeLikelihood interface.
const SubstitutionModel * getModelForNode(int nodeId) const
Get the model associated to a particular node id.
ParameterList getRootFrequenciesParameters() const
Get the parameters corresponding to the root frequencies.
std::vector< Node * > nodes_
Pointer toward all nodes in the tree.
const SubstitutionModel * getModel(size_t i) const
Get one model from the set knowing its index.
virtual const std::vector< int > & getAlphabetStates() const
Partial implementation of the DiscreteRatesAcrossSitesTreeLikelihood interface.
virtual void computeAllTransitionProbabilities()
Fill the pxy_, dpxy_ and d2pxy_ arrays for all nodes.
Iterates through all models used for all sites on a given branch.
ParameterList getRateDistributionParameters() const
Get the parameters associated to the rate distirbution.
virtual void applyParameters()
All parameters are stored in a parameter list. This function apply these parameters to the substituti...
ParameterList getRateDistributionParameters() const
Get the parameters associated to the rate distirbution.
virtual std::string getAlphabetStateAsChar(size_t index) const
Interface for phylogenetic tree objects.
Definition: Tree.h:148
Iterates through all models used for all branches on a given site.
virtual void initParameters()
This builds the parameters list from all parametrizable objects, i.e. substitution model...
ConstBranchModelIterator * getNewBranchModelIterator(int nodeId) const
const std::vector< double > & getRootFrequencies(size_t siteIndex) const
Get the values of the frequencies for each state in the alphabet at the root node.
VVVdouble getTransitionProbabilitiesPerRateClass(int nodeId, size_t siteIndex) const
Retrieves all Pij(t) for a particular branch, defined by the upper node.
virtual void computeTransitionProbabilitiesForNode(const Node *node)
Fill the pxy_, dpxy_ and d2pxy_ arrays for one node.
Exception thrown when something is wrong with a particular node.
std::map< int, const Node * > idToNode_
An index linking nodes to their id, for faster access than the getNode() method.
The phylogenetic node class.
Definition: Node.h:90
ParameterList getBranchLengthsParameters() const
Get the branch lengths parameters.
virtual int getAlphabetStateAsInt(size_t index) const
Specialization of the TreeLikelihood interface for the branch non-homogeneous and non-stationary mode...
const SubstitutionModel * getSubstitutionModelForNode(int nodeId) const
Get the substitution model associated to a given node.
size_t getNumberOfStates() const
Get the number of states associated to this model set.
const std::vector< int > & getNodesWithModel(size_t i) const
Get a list of nodes id for which the given model is associated.
SubstitutionModel * getSubstitutionModelForNode(int nodeId)
Get the substitution model associated to a given node.
A pair of SubstitutionModel / BranchIterator.
void init_(const Tree &tree, SubstitutionModelSet *modelSet, DiscreteDistribution *rDist, bool verbose)
Method called by constructor.
const SubstitutionModelSet * getSubstitutionModelSet() const
ConstSiteModelIterator * getNewSiteModelIterator(size_t siteIndex) const