bpp-phyl3  3.0.0
bpp::OnABranchPhyloLikelihood Class Referenceabstract

Wraps a dataflow graph as a function: resultNode = f(variableNodes). More...

#include <Bpp/Phyl/Likelihood/PhyloLikelihoods/OnABranchPhyloLikelihood.h>

+ Inheritance diagram for bpp::OnABranchPhyloLikelihood:
+ Collaboration diagram for bpp::OnABranchPhyloLikelihood:

Public Member Functions

 OnABranchPhyloLikelihood (Context &context, std::shared_ptr< LikelihoodCalculationSingleProcess > likCal, uint edgeId, const ParameterList &variableNodes)
 
 OnABranchPhyloLikelihood (Context &context, std::shared_ptr< LikelihoodCalculationSingleProcess > likCal, uint edgeId)
 
 OnABranchPhyloLikelihood (Context &context, std::shared_ptr< LikelihoodCalculationOnABranch > likCal, const ParameterList &variableNodes)
 
 OnABranchPhyloLikelihood (Context &context, std::shared_ptr< LikelihoodCalculationOnABranch > likCal)
 : the parameters the independent parameters of the LikelihoodCalculation More...
 
OnABranchPhyloLikelihoodclone () const override
 
size_t getNumberOfSites () const override
 Get the number of sites in the dataset. More...
 
size_t getNumberOfDistinctSites () const
 
bool isInitialized () const override
 Get the number of model classes. More...
 
void setModel (std::shared_ptr< ConfiguredModel > model)
 
LikelihoodCalculationlikelihoodCalculation () const override
 
std::shared_ptr< LikelihoodCalculationgetLikelihoodCalculation () const override
 
AlignedLikelihoodCalculationalignedLikelihoodCalculation () const override
 
std::shared_ptr< AlignedLikelihoodCalculationgetAlignedLikelihoodCalculation () const override
 
LikelihoodCalculationOnABranchlikelihoodCalculationOnABranch () const
 
std::shared_ptr< LikelihoodCalculationOnABranchgetLikelihoodCalculationOnABranch () const
 
ValueRef< RowLikgetFirstOrderDerivativeVector (const std::string &variable) const
 
ValueRef< RowLikfirstOrderDerivativeVector (const std::string &variable) const
 
ValueRef< RowLikgetSecondOrderDerivativeVector (const std::string &variable) const
 
ValueRef< RowLikgetSecondOrderDerivativeVector (const std::string &variable1, const std::string &variable2) const
 
ValueRef< RowLiksecondOrderDerivativeVector (const std::string &variable1, const std::string &variable2) const
 
VVDataLik getLikelihoodPerSitePerClass () const
 Get the posterior probabilities of each class, for each site. More...
 
std::vector< size_t > getClassWithMaxPostProbPerSite () const
 
DataLik getLikelihoodForASite (size_t site) const
 Get the likelihood for a site (on uncompressed data) More...
 
double getLogLikelihoodForASite (size_t site) const
 Get the log likelihood for a site, and its derivatives. More...
 
VDataLik getLikelihoodPerSite () const
 Get the likelihood for each site. More...
 
const Contextcontext () const override
 
Contextcontext () override
 
ValueRef< DataLikgetLikelihoodNode () const override
 
virtual void enableSecondOrderDerivatives (bool yn)=0
 
virtual bool enableSecondOrderDerivatives () const=0
 
virtual void enableSecondOrderDerivatives (bool yn) override
 
bool enableSecondOrderDerivatives () const override
 
virtual double getSecondOrderDerivative (const std::string &variable) const=0
 
virtual double getSecondOrderDerivative (const std::string &variable1, const std::string &variable2) const=0
 
double getSecondOrderDerivative (const std::string &variable) const override
 
double getSecondOrderDerivative (const std::string &variable1, const std::string &variable2) const override
 
virtual double d2f (const std::string &variable, const ParameterList &parameters)
 
virtual double d2f (const std::string &variable1, const std::string &variable2, const ParameterList &parameters)
 
virtual void enableFirstOrderDerivatives (bool yn)=0
 
virtual bool enableFirstOrderDerivatives () const=0
 
virtual void enableFirstOrderDerivatives (bool yn) override
 Tell if derivatives must be computed: for Function inheritance. More...
 
bool enableFirstOrderDerivatives () const override
 
virtual double getFirstOrderDerivative (const std::string &variable) const=0
 
double getFirstOrderDerivative (const std::string &variable) const override
 
virtual double df (const std::string &variable, const ParameterList &parameters)
 
virtual void setParameters (const ParameterList &parameters)=0
 
void setParameters (const ParameterList &parameters) override
 
virtual double getValue () const=0
 
double getValue () const override
 
virtual double f (const ParameterList &parameters)
 
virtual bool hasParameter (const std::string &name) const=0
 
virtual const ParameterListgetParameters () const=0
 
virtual const Parameterparameter (const std::string &name) const=0
 
virtual double getParameterValue (const std::string &name) const=0
 
virtual void setAllParametersValues (const ParameterList &parameters)=0
 
virtual void setParameterValue (const std::string &name, double value)=0
 
virtual void setParametersValues (const ParameterList &parameters)=0
 
virtual bool matchParametersValues (const ParameterList &parameters)=0
 
virtual void removeConstraint (const std::string &name)=0
 
virtual void setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint)=0
 
virtual size_t getNumberOfParameters () const=0
 
virtual void setNamespace (const std::string &prefix)=0
 
virtual std::string getNamespace () const=0
 
virtual std::string getParameterNameWithoutNamespace (const std::string &name) const=0
 
void shareParameters (const ParameterList &variableNodes)
 Share Parameters, that are DF_parameters. More...
 
ValueRef< DataLikfirstOrderDerivativeNode (const std::string &variable) const
 
ValueRef< DataLiksecondOrderDerivativeNode (const std::string &variable1, const std::string &variable2) const
 
bool hasParameter (const std::string &name) const override
 
const ParameterListgetParameters () const override
 
const Parameterparameter (const std::string &name) const override
 
const std::shared_ptr< Parameter > & getParameter (const std::string &name) const
 
double getParameterValue (const std::string &name) const override
 
void setAllParametersValues (const ParameterList &parameters) override
 
void setParameterValue (const std::string &name, double value) override
 
void setParametersValues (const ParameterList &parameters) override
 
bool matchParametersValues (const ParameterList &parameters) override
 
void removeConstraint (const std::string &name) override
 
void setConstraint (const std::string &name, std::shared_ptr< ConstraintInterface > constraint) override
 
size_t getNumberOfParameters () const override
 
void setNamespace (const std::string &prefix) override
 
std::string getNamespace () const override
 
std::string getParameterNameWithoutNamespace (const std::string &name) const override
 
virtual void fireParameterChanged (const ParameterList &parameters)
 
Retrieve some particular independent parameters subsets.
ParameterList getNonDerivableParameters () const override
 
ParameterList getDerivableParameters () const override
 
ParameterList getBranchLengthParameters () const override
 Get the independent branch lengths parameters. More...
 
ParameterList getSubstitutionModelParameters () const override
 Get the independent parameters associated to substitution model(s). More...
 
ParameterList getRateDistributionParameters () const override
 Get the independent parameters associated to the rate distribution(s). More...
 
ParameterList getRootFrequenciesParameters () const override
 Get the independent parameters associated to the root frequencies(s). More...
 
The data functions
virtual const Contextcontext () const =0
 
virtual Contextcontext ()=0
 
The likelihood functions.
virtual ValueRef< DataLikgetLikelihoodNode () const =0
 
double getLogLikelihood () const
 Get the logarithm of the likelihood for the whole dataset. More...
 

Protected Member Functions

void setNumberOfSites (size_t nbSites)
 
virtual ParameterListgetParameters_ ()=0
 
Node_DFaccessVariableNode (const std::string &name) const
 
const std::shared_ptr< Parameter > & getParameter (size_t i) const
 
std::shared_ptr< Parameter > & getParameter (size_t i)
 
virtual void addParameter_ (Parameter *parameter)
 
virtual void addParameters_ (const ParameterList &parameters)
 
virtual void shareParameter_ (const std::shared_ptr< Parameter > &parameter)
 
virtual void shareParameters_ (const ParameterList &parameters)
 
virtual void includeParameters_ (const ParameterList &parameters)
 
virtual void deleteParameter_ (size_t index)
 
virtual void deleteParameter_ (std::string &name)
 
virtual void deleteParameters_ (const std::vector< std::string > &names)
 
void resetParameters_ ()
 
ParametergetParameter_ (const std::string &name)
 
ParametergetParameter_ (size_t index)
 
const ParametergetParameter_ (size_t index) const
 
ParametergetParameterWithNamespace_ (const std::string &name)
 
const ParametergetParameterWithNamespace_ (const std::string &name) const
 
ParameterListgetParameters_ () override
 

Static Protected Member Functions

static Node_DFaccessVariableNode (const Parameter &param)
 

Protected Attributes

std::shared_ptr< LikelihoodCalculationOnABranchlikCal_
 
std::unordered_map< std::string, ValueRef< RowLik > > firstOrderDerivativeVectors_
 For Dataflow computing. More...
 
std::unordered_map< std::pair< std::string, std::string >, ValueRef< RowLik >, StringPairHashsecondOrderDerivativeVectors_
 
size_t nbSites_
 
Contextcontext_
 
DataLik minusLogLik_
 the value More...
 
std::unordered_map< std::string, ValueRef< DataLik > > firstOrderDerivativeNodes_
 For Dataflow computing. More...
 
std::unordered_map< std::pair< std::string, std::string >, ValueRef< DataLik >, StringPairHashsecondOrderDerivativeNodes_
 

Private Attributes

ParameterList parameters_
 
std::string prefix_
 

Detailed Description

Wraps a dataflow graph as a function: resultNode = f(variableNodes).

Definition at line 29 of file OnABranchPhyloLikelihood.h.

Constructor & Destructor Documentation

◆ OnABranchPhyloLikelihood() [1/4]

bpp::OnABranchPhyloLikelihood::OnABranchPhyloLikelihood ( Context context,
std::shared_ptr< LikelihoodCalculationSingleProcess likCal,
uint  edgeId,
const ParameterList variableNodes 
)
inline

Definition at line 47 of file OnABranchPhyloLikelihood.h.

References bpp::AbstractParametrizable::shareParameters_().

Referenced by clone().

◆ OnABranchPhyloLikelihood() [2/4]

bpp::OnABranchPhyloLikelihood::OnABranchPhyloLikelihood ( Context context,
std::shared_ptr< LikelihoodCalculationSingleProcess likCal,
uint  edgeId 
)
inline

Definition at line 59 of file OnABranchPhyloLikelihood.h.

◆ OnABranchPhyloLikelihood() [3/4]

bpp::OnABranchPhyloLikelihood::OnABranchPhyloLikelihood ( Context context,
std::shared_ptr< LikelihoodCalculationOnABranch likCal,
const ParameterList variableNodes 
)
inline

◆ OnABranchPhyloLikelihood() [4/4]

bpp::OnABranchPhyloLikelihood::OnABranchPhyloLikelihood ( Context context,
std::shared_ptr< LikelihoodCalculationOnABranch likCal 
)
inline

: the parameters the independent parameters of the LikelihoodCalculation

Definition at line 83 of file OnABranchPhyloLikelihood.h.

References likCal_, and bpp::AbstractParametrizable::shareParameters_().

Member Function Documentation

◆ accessVariableNode() [1/2]

◆ accessVariableNode() [2/2]

Node_DF& bpp::AbstractPhyloLikelihood::accessVariableNode ( const std::string &  name) const
inlineprotectedinherited

◆ alignedLikelihoodCalculation()

AlignedLikelihoodCalculation& bpp::OnABranchPhyloLikelihood::alignedLikelihoodCalculation ( ) const
inlineoverridevirtual
Returns
The LikelihoodCalculation.

Implements bpp::AlignedPhyloLikelihoodInterface.

Definition at line 212 of file OnABranchPhyloLikelihood.h.

References likCal_.

◆ clone()

OnABranchPhyloLikelihood* bpp::OnABranchPhyloLikelihood::clone ( ) const
inlineoverridevirtual

Implements bpp::AbstractParametrizable.

Definition at line 100 of file OnABranchPhyloLikelihood.h.

References OnABranchPhyloLikelihood().

◆ context() [1/4]

virtual const Context& bpp::PhyloLikelihoodInterface::context ( ) const
pure virtualinherited

◆ context() [2/4]

◆ context() [3/4]

Context& bpp::AbstractPhyloLikelihood::context ( )
inlineoverridevirtualinherited

◆ context() [4/4]

virtual Context& bpp::PhyloLikelihoodInterface::context ( )
pure virtualinherited

◆ enableFirstOrderDerivatives() [1/2]

bool bpp::AbstractPhyloLikelihood::enableFirstOrderDerivatives ( ) const
inlineoverridevirtualinherited

Implements bpp::SecondOrderDerivable.

Definition at line 115 of file AbstractPhyloLikelihood.h.

◆ enableFirstOrderDerivatives() [2/2]

virtual void bpp::AbstractPhyloLikelihood::enableFirstOrderDerivatives ( bool  yn)
inlineoverridevirtualinherited

Tell if derivatives must be computed: for Function inheritance.

Implements bpp::SecondOrderDerivable.

Definition at line 113 of file AbstractPhyloLikelihood.h.

◆ enableSecondOrderDerivatives() [1/2]

bool bpp::AbstractPhyloLikelihood::enableSecondOrderDerivatives ( ) const
inlineoverridevirtualinherited

Implements bpp::SecondOrderDerivable.

Definition at line 116 of file AbstractPhyloLikelihood.h.

◆ enableSecondOrderDerivatives() [2/2]

virtual void bpp::AbstractPhyloLikelihood::enableSecondOrderDerivatives ( bool  yn)
inlineoverridevirtualinherited

Implements bpp::SecondOrderDerivable.

Definition at line 114 of file AbstractPhyloLikelihood.h.

◆ firstOrderDerivativeNode()

◆ firstOrderDerivativeVector()

ValueRef<RowLik> bpp::OnABranchPhyloLikelihood::firstOrderDerivativeVector ( const std::string &  variable) const
inline

◆ getAlignedLikelihoodCalculation()

std::shared_ptr<AlignedLikelihoodCalculation> bpp::OnABranchPhyloLikelihood::getAlignedLikelihoodCalculation ( ) const
inlineoverridevirtual
Returns
A shared pointer toward the LikelihoodCalculation.

Implements bpp::AlignedPhyloLikelihoodInterface.

Definition at line 217 of file OnABranchPhyloLikelihood.h.

References likCal_.

◆ getBranchLengthParameters()

ParameterList bpp::OnABranchPhyloLikelihood::getBranchLengthParameters ( ) const
inlineoverridevirtual

Get the independent branch lengths parameters.

Returns
A ParameterList with all branch lengths.

Implements bpp::PhyloLikelihoodInterface.

Definition at line 162 of file OnABranchPhyloLikelihood.h.

◆ getClassWithMaxPostProbPerSite()

vector< size_t > OnABranchPhyloLikelihood::getClassWithMaxPostProbPerSite ( ) const

Definition at line 82 of file OnABranchPhyloLikelihood.cpp.

◆ getDerivableParameters()

ParameterList bpp::OnABranchPhyloLikelihood::getDerivableParameters ( ) const
inlineoverridevirtual

Implements bpp::PhyloLikelihoodInterface.

Definition at line 152 of file OnABranchPhyloLikelihood.h.

◆ getFirstOrderDerivative()

double bpp::AbstractPhyloLikelihood::getFirstOrderDerivative ( const std::string &  variable) const
inlineoverridevirtualinherited

◆ getFirstOrderDerivativeVector()

ValueRef<RowLik> bpp::OnABranchPhyloLikelihood::getFirstOrderDerivativeVector ( const std::string &  variable) const
inline

Definition at line 234 of file OnABranchPhyloLikelihood.h.

References firstOrderDerivativeVector().

◆ getLikelihoodCalculation()

std::shared_ptr<LikelihoodCalculation> bpp::OnABranchPhyloLikelihood::getLikelihoodCalculation ( ) const
inlineoverridevirtual
Returns
A shared pointer toward the LikelihoodCalculation.

Implements bpp::PhyloLikelihoodInterface.

Definition at line 207 of file OnABranchPhyloLikelihood.h.

References likCal_.

◆ getLikelihoodCalculationOnABranch()

std::shared_ptr<LikelihoodCalculationOnABranch> bpp::OnABranchPhyloLikelihood::getLikelihoodCalculationOnABranch ( ) const
inline

Definition at line 227 of file OnABranchPhyloLikelihood.h.

References likCal_.

Referenced by firstOrderDerivativeVector().

◆ getLikelihoodForASite()

DataLik bpp::AbstractAlignedPhyloLikelihood::getLikelihoodForASite ( size_t  site) const
inlinevirtualinherited

Get the likelihood for a site (on uncompressed data)

Parameters
siteThe site index to analyse.
Returns
The likelihood for site site.

Implements bpp::AlignedPhyloLikelihoodInterface.

Definition at line 143 of file AlignedPhyloLikelihood.h.

References bpp::AlignedPhyloLikelihoodInterface::alignedLikelihoodCalculation(), and bpp::AlignedLikelihoodCalculation::getLikelihoodForASite().

◆ getLikelihoodNode() [1/2]

virtual ValueRef<DataLik> bpp::PhyloLikelihoodInterface::getLikelihoodNode ( ) const
pure virtualinherited
Returns
the LikDF node where the Likelihood is computed.

Implemented in bpp::AbstractPhyloLikelihood.

◆ getLikelihoodNode() [2/2]

ValueRef<DataLik> bpp::AbstractPhyloLikelihood::getLikelihoodNode ( ) const
inlineoverridevirtualinherited

◆ getLikelihoodPerSite()

VDataLik bpp::AbstractAlignedPhyloLikelihood::getLikelihoodPerSite ( ) const
inlinevirtualinherited

◆ getLikelihoodPerSitePerClass()

VVDataLik OnABranchPhyloLikelihood::getLikelihoodPerSitePerClass ( ) const

Get the posterior probabilities of each class, for each site.

Returns
A 2D-vector (Sites X Class) of all posterior probabilities.

Definition at line 72 of file OnABranchPhyloLikelihood.cpp.

References bpp::copyEigenToBpp().

◆ getLogLikelihood()

double bpp::PhyloLikelihoodInterface::getLogLikelihood ( ) const
inlineinherited

Get the logarithm of the likelihood for the whole dataset.

Returns
The logarithm of the likelihood of the dataset.

Definition at line 66 of file PhyloLikelihood.h.

References bpp::SecondOrderDerivable::getValue().

◆ getLogLikelihoodForASite()

double bpp::AbstractAlignedPhyloLikelihood::getLogLikelihoodForASite ( size_t  site) const
inlinevirtualinherited

Get the log likelihood for a site, and its derivatives.

Parameters
siteThe site index to analyse.
Returns
The (D)log likelihood for site site.

Implements bpp::AlignedPhyloLikelihoodInterface.

Definition at line 154 of file AlignedPhyloLikelihood.h.

References bpp::convert(), and bpp::AlignedPhyloLikelihoodInterface::getAlignedLikelihoodCalculation().

◆ getNonDerivableParameters()

ParameterList bpp::OnABranchPhyloLikelihood::getNonDerivableParameters ( ) const
inlineoverridevirtual

Implements bpp::PhyloLikelihoodInterface.

Definition at line 147 of file OnABranchPhyloLikelihood.h.

◆ getNumberOfDistinctSites()

size_t bpp::OnABranchPhyloLikelihood::getNumberOfDistinctSites ( ) const
inline

◆ getNumberOfSites()

size_t bpp::OnABranchPhyloLikelihood::getNumberOfSites ( ) const
inlineoverridevirtual

Get the number of sites in the dataset.

Returns
the number of sites in the dataset.

Reimplemented from bpp::AbstractAlignedPhyloLikelihood.

Definition at line 107 of file OnABranchPhyloLikelihood.h.

References bpp::LikelihoodCalculationOnABranch::getNumberOfSites(), and likelihoodCalculationOnABranch().

◆ getRateDistributionParameters()

ParameterList bpp::OnABranchPhyloLikelihood::getRateDistributionParameters ( ) const
inlineoverridevirtual

Get the independent parameters associated to the rate distribution(s).

Returns
A ParameterList.

Implements bpp::PhyloLikelihoodInterface.

Definition at line 182 of file OnABranchPhyloLikelihood.h.

◆ getRootFrequenciesParameters()

ParameterList bpp::OnABranchPhyloLikelihood::getRootFrequenciesParameters ( ) const
inlineoverridevirtual

Get the independent parameters associated to the root frequencies(s).

Returns
A ParameterList.

Implements bpp::PhyloLikelihoodInterface.

Definition at line 193 of file OnABranchPhyloLikelihood.h.

◆ getSecondOrderDerivative() [1/2]

double bpp::AbstractPhyloLikelihood::getSecondOrderDerivative ( const std::string &  variable) const
inlineoverridevirtualinherited

Implements bpp::SecondOrderDerivable.

Definition at line 155 of file AbstractPhyloLikelihood.h.

◆ getSecondOrderDerivative() [2/2]

double bpp::AbstractPhyloLikelihood::getSecondOrderDerivative ( const std::string &  variable1,
const std::string &  variable2 
) const
inlineoverridevirtualinherited

◆ getSecondOrderDerivativeVector() [1/2]

ValueRef<RowLik> bpp::OnABranchPhyloLikelihood::getSecondOrderDerivativeVector ( const std::string &  variable) const
inline

Definition at line 254 of file OnABranchPhyloLikelihood.h.

◆ getSecondOrderDerivativeVector() [2/2]

ValueRef<RowLik> bpp::OnABranchPhyloLikelihood::getSecondOrderDerivativeVector ( const std::string &  variable1,
const std::string &  variable2 
) const
inline

Definition at line 259 of file OnABranchPhyloLikelihood.h.

References secondOrderDerivativeVector().

◆ getSubstitutionModelParameters()

ParameterList bpp::OnABranchPhyloLikelihood::getSubstitutionModelParameters ( ) const
inlineoverridevirtual

Get the independent parameters associated to substitution model(s).

Returns
A ParameterList.

Implements bpp::PhyloLikelihoodInterface.

Definition at line 172 of file OnABranchPhyloLikelihood.h.

◆ getValue()

double bpp::AbstractPhyloLikelihood::getValue ( ) const
inlineoverridevirtualinherited

◆ isInitialized()

bool bpp::OnABranchPhyloLikelihood::isInitialized ( ) const
inlineoverridevirtual

Get the number of model classes.

Returns
'true' is the likelihood function has been initialized.

Reimplemented from bpp::AbstractPhyloLikelihood.

Definition at line 132 of file OnABranchPhyloLikelihood.h.

References bpp::LikelihoodCalculation::isInitialized(), and likelihoodCalculationOnABranch().

◆ likelihoodCalculation()

LikelihoodCalculation& bpp::OnABranchPhyloLikelihood::likelihoodCalculation ( ) const
inlineoverridevirtual
Returns
The LikelihoodCalculation.

Implements bpp::PhyloLikelihoodInterface.

Definition at line 202 of file OnABranchPhyloLikelihood.h.

References likCal_.

◆ likelihoodCalculationOnABranch()

LikelihoodCalculationOnABranch& bpp::OnABranchPhyloLikelihood::likelihoodCalculationOnABranch ( ) const
inline

Definition at line 222 of file OnABranchPhyloLikelihood.h.

References likCal_.

Referenced by getNumberOfDistinctSites(), getNumberOfSites(), and isInitialized().

◆ secondOrderDerivativeNode()

ValueRef<DataLik> bpp::AbstractPhyloLikelihood::secondOrderDerivativeNode ( const std::string &  variable1,
const std::string &  variable2 
) const
inlineinherited

◆ secondOrderDerivativeVector()

ValueRef<RowLik> bpp::OnABranchPhyloLikelihood::secondOrderDerivativeVector ( const std::string &  variable1,
const std::string &  variable2 
) const
inline

◆ setModel()

void bpp::OnABranchPhyloLikelihood::setModel ( std::shared_ptr< ConfiguredModel model)
inline

Definition at line 137 of file OnABranchPhyloLikelihood.h.

References likCal_.

◆ setNumberOfSites()

void bpp::AbstractAlignedPhyloLikelihood::setNumberOfSites ( size_t  nbSites)
inlineprotectedinherited

◆ setParameters()

void bpp::AbstractPhyloLikelihood::setParameters ( const ParameterList parameters)
inlineoverridevirtualinherited

◆ shareParameters()

void bpp::AbstractPhyloLikelihood::shareParameters ( const ParameterList variableNodes)
inlineinherited

Member Data Documentation

◆ context_

◆ firstOrderDerivativeNodes_

std::unordered_map<std::string, ValueRef<DataLik> > bpp::AbstractPhyloLikelihood::firstOrderDerivativeNodes_
mutableprotectedinherited

For Dataflow computing.

Definition at line 44 of file AbstractPhyloLikelihood.h.

Referenced by bpp::AbstractPhyloLikelihood::firstOrderDerivativeNode().

◆ firstOrderDerivativeVectors_

std::unordered_map<std::string, ValueRef<RowLik> > bpp::OnABranchPhyloLikelihood::firstOrderDerivativeVectors_
mutableprotected

For Dataflow computing.

Definition at line 40 of file OnABranchPhyloLikelihood.h.

Referenced by firstOrderDerivativeVector().

◆ likCal_

◆ minusLogLik_

DataLik bpp::AbstractPhyloLikelihood::minusLogLik_
mutableprotectedinherited

◆ nbSites_

◆ secondOrderDerivativeNodes_

std::unordered_map<std::pair<std::string, std::string>, ValueRef<DataLik>, StringPairHash> bpp::AbstractPhyloLikelihood::secondOrderDerivativeNodes_
mutableprotectedinherited

◆ secondOrderDerivativeVectors_

std::unordered_map<std::pair<std::string, std::string>, ValueRef<RowLik>, StringPairHash> bpp::OnABranchPhyloLikelihood::secondOrderDerivativeVectors_
mutableprotected

Definition at line 44 of file OnABranchPhyloLikelihood.h.

Referenced by secondOrderDerivativeVector().


The documentation for this class was generated from the following files: