void AutocorrelatedBranchMatrixDistribution::recursiveSimulate(const TopologyNode& node, RbVector< RateMatrix > *values, const std::vector< double > &scaledParent) { // get the index size_t nodeIndex = node.getIndex(); // first we simulate our value RandomNumberGenerator* rng = GLOBAL_RNG; // do we keep our parents values? double u = rng->uniform01(); if ( u < changeProbability->getValue() ) { // change // draw a new value for the base frequencies std::vector<double> newParent = RbStatistics::Dirichlet::rv(scaledParent, *rng); std::vector<double> newScaledParent = newParent; // compute the new scaled parent std::vector<double>::iterator end = newScaledParent.end(); for (std::vector<double>::iterator it = newScaledParent.begin(); it != end; ++it) { (*it) *= alpha->getValue(); } RateMatrix_GTR rm = RateMatrix_GTR( newParent.size() ); RbPhylogenetics::Gtr::computeRateMatrix( exchangeabilityRates->getValue(), newParent, &rm ); uniqueBaseFrequencies.push_back( newParent ); uniqueMatrices.push_back( rm ); matrixIndex[nodeIndex] = uniqueMatrices.size()-1; values->insert(nodeIndex, rm); size_t numChildren = node.getNumberOfChildren(); if ( numChildren > 0 ) { for (size_t i = 0; i < numChildren; ++i) { const TopologyNode& child = node.getChild(i); recursiveSimulate(child,values,newScaledParent); } } } else { // no change size_t parentIndex = node.getParent().getIndex(); values->insert(nodeIndex, uniqueMatrices[ matrixIndex[ parentIndex ] ]); size_t numChildren = node.getNumberOfChildren(); if ( numChildren > 0 ) { for (size_t i = 0; i < numChildren; ++i) { const TopologyNode& child = node.getChild(i); recursiveSimulate(child,values,scaledParent); } } } }
void MultivariateBrownianPhyloProcess::recursiveCorruptAll(const TopologyNode& from) { dirtyNodes[from.getIndex()] = true; for (size_t i = 0; i < from.getNumberOfChildren(); ++i) { recursiveCorruptAll(from.getChild(i)); } }
void BrownianPhyloProcess::recursiveSimulate(const TopologyNode& from) { size_t index = from.getIndex(); if (! from.isRoot()) { // x ~ normal(x_up, sigma^2 * branchLength) size_t upindex = from.getParent().getIndex(); double standDev = sigma->getValue() * sqrt(from.getBranchLength()); double mean = (*value)[upindex] + drift->getValue() * from.getBranchLength(); // simulate the new Val RandomNumberGenerator* rng = GLOBAL_RNG; (*value)[index] = RbStatistics::Normal::rv( mean, standDev, *rng); } // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { recursiveSimulate(from.getChild(i)); } }
double AutocorrelatedLognormalRateBranchwiseVarDistribution::recursiveLnProb( const TopologyNode& n ) { // get the index size_t nodeIndex = n.getIndex(); double lnProb = 0.0; size_t numChildren = n.getNumberOfChildren(); double scale = scaleValue->getValue(); if ( numChildren > 0 ) { double parentRate = log( (*value)[nodeIndex] ); for (size_t i = 0; i < numChildren; ++i) { const TopologyNode& child = n.getChild(i); lnProb += recursiveLnProb(child); size_t childIndex = child.getIndex(); // compute the variance double variance = sigma->getValue()[childIndex] * child.getBranchLength() * scale; double childRate = (*value)[childIndex]; // the mean of the LN dist is parentRate = exp[mu + (variance / 2)], // where mu is the location param of the LN dist (see Kishino & Thorne 2001) double mu = parentRate - (variance * 0.5); double stDev = sqrt(variance); lnProb += RbStatistics::Lognormal::lnPdf(mu, stDev, childRate); } } return lnProb; }
void MultivariateBrownianPhyloProcess::recursiveSimulate(const TopologyNode& from) { size_t index = from.getIndex(); if (from.isRoot()) { std::vector<double>& val = (*value)[index]; for (size_t i=0; i<getDim(); i++) { val[i] = 0; } } else { // x ~ normal(x_up, sigma^2 * branchLength) std::vector<double>& val = (*value)[index]; sigma->getValue().drawNormalSampleCovariance((*value)[index]); size_t upindex = from.getParent().getIndex(); std::vector<double>& upval = (*value)[upindex]; for (size_t i=0; i<getDim(); i++) { val[i] += upval[i]; } } // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { recursiveSimulate(from.getChild(i)); } }
double BrownianPhyloProcess::recursiveLnProb( const TopologyNode& from ) { double lnProb = 0.0; size_t index = from.getIndex(); double val = (*value)[index]; if (! from.isRoot()) { // x ~ normal(x_up, sigma^2 * branchLength) size_t upindex = from.getParent().getIndex(); double standDev = sigma->getValue() * sqrt(from.getBranchLength()); double mean = (*value)[upindex] + drift->getValue() * from.getBranchLength(); lnProb += RbStatistics::Normal::lnPdf(val, standDev, mean); } // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { lnProb += recursiveLnProb(from.getChild(i)); } return lnProb; }
void PhyloWhiteNoiseProcess::recursiveSimulate(const TopologyNode& from) { if (! from.isRoot()) { // get the index size_t index = from.getIndex(); // compute the variance along the branch double mean = 1.0; double stdev = sigma->getValue() / sqrt(from.getBranchLength()); double alpha = mean * mean / (stdev * stdev); double beta = mean / (stdev * stdev); // simulate a new Val RandomNumberGenerator* rng = GLOBAL_RNG; double v = RbStatistics::Gamma::rv( alpha,beta, *rng); // we store this val here (*value)[index] = v; } // simulate the val for each child (if any) size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { const TopologyNode& child = from.getChild(i); recursiveSimulate(child); } }
void AutocorrelatedLognormalRateBranchwiseVarDistribution::recursiveSimulate(const TopologyNode& node, double parentRate) { // get the index size_t nodeIndex = node.getIndex(); // compute the variance along the branch double scale = scaleValue->getValue(); double variance = sigma->getValue()[nodeIndex] * node.getBranchLength() * scale; double mu = log(parentRate) - (variance * 0.5); double stDev = sqrt(variance); // simulate a new rate RandomNumberGenerator* rng = GLOBAL_RNG; double nodeRate = RbStatistics::Lognormal::rv( mu, stDev, *rng ); // we store this rate here (*value)[nodeIndex] = nodeRate; // simulate the rate for each child (if any) size_t numChildren = node.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { const TopologyNode& child = node.getChild(i); recursiveSimulate(child,nodeRate); } }
void RealNodeContainer::recursiveClampAt(const TopologyNode& from, const ContinuousCharacterData* data, size_t l) { if (from.isTip()) { // get taxon index size_t index = from.getIndex(); std::string taxon = tree->getTipNames()[index]; size_t dataindex = data->getIndexOfTaxon(taxon); if (data->getCharacter(dataindex,l) != -1000) { (*this)[index] = data->getCharacter(dataindex,l); clampVector[index] = true; //std::cerr << "taxon : " << index << '\t' << taxon << " trait value : " << (*this)[index] << '\n'; } else { std::cerr << "taxon : " << taxon << " is missing for trait " << l+1 << '\n'; } } // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { recursiveClampAt(from.getChild(i),data,l); } }
void RealNodeContainer::recursiveGetStats(const TopologyNode& from, double& e1, double& e2, int& n) const { double tmp = (*this)[from.getIndex()]; n++; e1 += tmp; e2 += tmp * tmp; // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { recursiveGetStats(from.getChild(i),e1,e2,n); } }
double AutocorrelatedBranchMatrixDistribution::recursiveLnProb( const TopologyNode& n ) { // get the index size_t nodeIndex = n.getIndex(); double lnProb = 0.0; size_t numChildren = n.getNumberOfChildren(); if ( numChildren > 0 ) { std::vector<double> parent = (*value)[nodeIndex].getStationaryFrequencies(); std::vector<double>::iterator end = parent.end(); for (std::vector<double>::iterator it = parent.begin(); it != end; ++it) { (*it) *= alpha->getValue(); } for (size_t i = 0; i < numChildren; ++i) { const TopologyNode& child = n.getChild(i); lnProb += recursiveLnProb(child); size_t childIndex = child.getIndex(); // RateMatrix& rm = (*value)[childIndex]; // compare if the child has a different matrix if ( matrixIndex[nodeIndex] == matrixIndex[childIndex] ) { // no change -> just the probability of no change lnProb += log( 1.0 - changeProbability->getValue() ); } else { // change: // probability of change lnProb += log( changeProbability->getValue() ); const std::vector<double>& descendant = (*value)[childIndex].getStationaryFrequencies(); // const std::vector<double>& descendant = uniqueMatrices[ matrixIndex[childIndex] ].getStationaryFrequencies(); // probability of new descendant values double p = RbStatistics::Dirichlet::lnPdf(parent, descendant); lnProb += p; } } } return lnProb; }
double MultivariateBrownianPhyloProcess::recursiveLnProb( const TopologyNode& from ) { double lnProb = 0.0; size_t index = from.getIndex(); std::vector<double> val = (*value)[index]; if (! from.isRoot()) { if (1) { // if (dirtyNodes[index]) { // x ~ normal(x_up, sigma^2 * branchLength) size_t upindex = from.getParent().getIndex(); std::vector<double> upval = (*value)[upindex]; const MatrixReal& om = sigma->getValue().getInverse(); double s2 = 0; for (size_t i = 0; i < getDim(); i++) { double tmp = 0; for (size_t j = 0; j < getDim(); j++) { tmp += om[i][j] * (val[j] - upval[j]); } s2 += (val[i] - upval[i]) * tmp; } double logprob = 0; logprob -= 0.5 * s2 / from.getBranchLength(); logprob -= 0.5 * (sigma->getValue().getLogDet() + sigma->getValue().getDim() * log(from.getBranchLength())); nodeLogProbs[index] = logprob; dirtyNodes[index] = false; } lnProb += nodeLogProbs[index]; } // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { lnProb += recursiveLnProb(from.getChild(i)); } return lnProb; }
void RealNodeContainer::recursiveGetTipValues(const TopologyNode& from, ContinuousCharacterData& nameToVal) const { if(from.isTip()) { double tmp = (*this)[from.getIndex()]; std::string name = tree->getTipNames()[from.getIndex()]; ContinuousTaxonData dataVec = ContinuousTaxonData(name); double contObs = tmp; dataVec.addCharacter( contObs ); nameToVal.addTaxonData( dataVec ); return; } // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { recursiveGetTipValues(from.getChild(i), nameToVal ); } }
std::string RealNodeContainer::recursiveGetNewick(const TopologyNode& from) const { std::ostringstream s; if (from.isTip()) { s << getTimeTree()->getTipNames()[from.getIndex()] << "_"; // std::cerr << from.getIndex() << '\t' << getTimeTree()->getTipNames()[from.getIndex()] << "_"; // std::cerr << (*this)[from.getIndex()] << '\n'; // exit(1); } else { s << "("; // propagate forward size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { s << recursiveGetNewick(from.getChild(i)); if (i < numChildren-1) { s << ","; } } s << ")"; } s << (*this)[from.getIndex()]; /* if (from.isTip() && (! isClamped(from.getIndex()))) { std::cerr << "leaf is not clamped\n"; // get taxon index size_t index = from.getIndex(); std::cerr << "index : " << index << '\n'; std::string taxon = tree->getTipNames()[index]; std::cerr << "taxon : " << index << '\t' << taxon << '\n'; std::cerr << " trait value : " << (*this)[index] << '\n'; exit(1); }*/ // if (!from.isRoot()) { s << ":"; s << getTimeTree()->getBranchLength(from.getIndex()); // } return s.str(); }
double PhyloWhiteNoiseProcess::recursiveLnProb(const TopologyNode &from) { double lnProb = 0.0; if (! from.isRoot()) { // compute the variance double mean = 1.0; double stdev = sigma->getValue() / sqrt(from.getBranchLength()); double alpha = mean * mean / (stdev * stdev); double beta = mean / (stdev * stdev); double v = (*value)[from.getIndex()]; lnProb += log( RbStatistics::Gamma::lnPdf(alpha,beta,v) ); } size_t numChildren = from.getNumberOfChildren(); for (size_t i = 0; i < numChildren; ++i) { const TopologyNode& child = from.getChild(i); lnProb += recursiveLnProb(child); } return lnProb; }
Tree* NewickConverter::convertFromNewick(std::string const &n) { // create and allocate the tree object Tree *tau = new Tree(); std::vector<TopologyNode*> nodes; std::vector<double> brlens; // create a string-stream and throw the string into it std::stringstream ss (std::stringstream::in | std::stringstream::out); ss << n; // ignore white spaces std::string trimmed = ""; char c; while ( ss.good() ) { c = ss.get(); if ( c != ' ') trimmed += c; } // construct the tree starting from the root TopologyNode *root = createNode( trimmed, nodes, brlens); // set up the tree tau->setRoot( root ); // set the branch lengths for (size_t i = 0; i < nodes.size(); ++i) { nodes[i]->setBranchLength(brlens[i]); } if(root->getNumberOfChildren() == 2) tau->setRooted(true); // return the tree, the caller is responsible for destruction return tau; }
/** Perform the move */ void RateAgeBetaShift::performMcmcMove( double lHeat, double pHeat ) { // Get random number generator RandomNumberGenerator* rng = GLOBAL_RNG; Tree& tau = tree->getValue(); RbOrderedSet<DagNode*> affected; tree->getAffectedNodes( affected ); double oldLnLike = 0.0; bool checkLikelihoodShortcuts = rng->uniform01() < 0.001; if ( checkLikelihoodShortcuts == true ) { for (RbOrderedSet<DagNode*>::iterator it = affected.begin(); it != affected.end(); ++it) { (*it)->touch(); oldLnLike += (*it)->getLnProbability(); } } // pick a random node which is not the root and neithor the direct descendant of the root TopologyNode* node; size_t nodeIdx = 0; do { double u = rng->uniform01(); nodeIdx = size_t( std::floor(tau.getNumberOfNodes() * u) ); node = &tau.getNode(nodeIdx); } while ( node->isRoot() || node->isTip() ); TopologyNode& parent = node->getParent(); // we need to work with the times double parent_age = parent.getAge(); double my_age = node->getAge(); double child_Age = node->getChild( 0 ).getAge(); if ( child_Age < node->getChild( 1 ).getAge()) { child_Age = node->getChild( 1 ).getAge(); } // now we store all necessary values storedNode = node; storedAge = my_age; storedRates[nodeIdx] = rates[nodeIdx]->getValue(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); storedRates[childIdx] = rates[childIdx]->getValue(); } // draw new ages and compute the hastings ratio at the same time double m = (my_age-child_Age) / (parent_age-child_Age); double a = delta * m + 1.0; double b = delta * (1.0-m) + 1.0; double new_m = RbStatistics::Beta::rv(a, b, *rng); double my_new_age = (parent_age-child_Age) * new_m + child_Age; // compute the Hastings ratio double forward = RbStatistics::Beta::lnPdf(a, b, new_m); double new_a = delta * new_m + 1.0; double new_b = delta * (1.0-new_m) + 1.0; double backward = RbStatistics::Beta::lnPdf(new_a, new_b, m); // set the age tau.getNode(nodeIdx).setAge( my_new_age ); // touch the tree so that the likelihoods are getting stored tree->touch(); // get the probability ratio of the tree double treeProbRatio = tree->getLnProbabilityRatio(); // set the rates double pa = node->getParent().getAge(); double my_new_rate =(pa - my_age) * storedRates[nodeIdx] / (pa - my_new_age); // now we set the new value // this will automcatically call a touch rates[nodeIdx]->setValue( new double( my_new_rate ) ); // get the probability ratio of the new rate double ratesProbRatio = rates[nodeIdx]->getLnProbabilityRatio(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); double a = node->getChild(i).getAge(); double child_new_rate = (my_age - a) * storedRates[childIdx] / (my_new_age - a); // now we set the new value // this will automcatically call a touch rates[childIdx]->setValue( new double( child_new_rate ) ); // get the probability ratio of the new rate ratesProbRatio += rates[childIdx]->getLnProbabilityRatio(); } if ( checkLikelihoodShortcuts == true ) { double lnProbRatio = 0; double newLnLike = 0; for (RbOrderedSet<DagNode*>::iterator it = affected.begin(); it != affected.end(); ++it) { double tmp = (*it)->getLnProbabilityRatio(); lnProbRatio += tmp; newLnLike += (*it)->getLnProbability(); } if ( fabs(lnProbRatio) > 1E-8 ) { double lnProbRatio2 = 0; double newLnLike2 = 0; for (RbOrderedSet<DagNode*>::iterator it = affected.begin(); it != affected.end(); ++it) { double tmp2 = (*it)->getLnProbabilityRatio(); lnProbRatio2 += tmp2; newLnLike2 += (*it)->getLnProbability(); } throw RbException("Likelihood shortcut computation failed in rate-age-proposal."); } } double hastingsRatio = backward - forward; double ln_acceptance_ratio = lHeat * pHeat * (treeProbRatio + ratesProbRatio) + hastingsRatio; if (ln_acceptance_ratio >= 0.0) { numAccepted++; tree->keep(); rates[nodeIdx]->keep(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->keep(); } } else if (ln_acceptance_ratio < -300.0) { reject(); tree->restore(); rates[nodeIdx]->restore(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->restore(); } } else { double r = exp(ln_acceptance_ratio); // Accept or reject the move double u = GLOBAL_RNG->uniform01(); if (u < r) { numAccepted++; //keep tree->keep(); rates[nodeIdx]->keep(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->keep(); } } else { reject(); tree->restore(); rates[nodeIdx]->restore(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->restore(); } } } }
/** Perform the move */ void RateAgeBetaShift::performMove( double heat, bool raiseLikelihoodOnly ) { // Get random number generator RandomNumberGenerator* rng = GLOBAL_RNG; TimeTree& tau = tree->getValue(); // pick a random node which is not the root and neithor the direct descendant of the root TopologyNode* node; size_t nodeIdx = 0; do { double u = rng->uniform01(); nodeIdx = size_t( std::floor(tau.getNumberOfNodes() * u) ); node = &tau.getNode(nodeIdx); } while ( node->isRoot() || node->isTip() ); TopologyNode& parent = node->getParent(); // we need to work with the times double parent_age = parent.getAge(); double my_age = node->getAge(); double child_Age = node->getChild( 0 ).getAge(); if ( child_Age < node->getChild( 1 ).getAge()) { child_Age = node->getChild( 1 ).getAge(); } // now we store all necessary values storedNode = node; storedAge = my_age; storedRates[nodeIdx] = rates[nodeIdx]->getValue(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); storedRates[childIdx] = rates[childIdx]->getValue(); } // draw new ages and compute the hastings ratio at the same time double m = (my_age-child_Age) / (parent_age-child_Age); double a = delta * m + 1.0; double b = delta * (1.0-m) + 1.0; double new_m = RbStatistics::Beta::rv(a, b, *rng); double my_new_age = (parent_age-child_Age) * new_m + child_Age; // compute the Hastings ratio double forward = RbStatistics::Beta::lnPdf(a, b, new_m); double new_a = delta * new_m + 1.0; double new_b = delta * (1.0-new_m) + 1.0; double backward = RbStatistics::Beta::lnPdf(new_a, new_b, m); // set the age tau.setAge( node->getIndex(), my_new_age ); tree->touch(); double treeProbRatio = tree->getLnProbabilityRatio(); // set the rates rates[nodeIdx]->setValue( new double((node->getParent().getAge() - my_age) * storedRates[nodeIdx] / (node->getParent().getAge() - my_new_age))); double ratesProbRatio = rates[nodeIdx]->getLnProbabilityRatio(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->setValue( new double((my_age - node->getChild(i).getAge()) * storedRates[childIdx] / (my_new_age - node->getChild(i).getAge()))); ratesProbRatio += rates[childIdx]->getLnProbabilityRatio(); } std::set<DagNode*> affected; tree->getAffectedNodes( affected ); double lnProbRatio = 0; for (std::set<DagNode*>::iterator it = affected.begin(); it != affected.end(); ++it) { (*it)->touch(); lnProbRatio += (*it)->getLnProbabilityRatio(); } if ( fabs(lnProbRatio) > 1E-6 ) { // throw RbException("Likelihood shortcut computation failed in rate-age-proposal."); std::cout << "Likelihood shortcut computation failed in rate-age-proposal." << std::endl; } double hastingsRatio = backward - forward; double lnAcceptanceRatio = treeProbRatio + ratesProbRatio + hastingsRatio; if (lnAcceptanceRatio >= 0.0) { numAccepted++; tree->keep(); rates[nodeIdx]->keep(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->keep(); } } else if (lnAcceptanceRatio < -300.0) { reject(); tree->restore(); rates[nodeIdx]->restore(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->restore(); } } else { double r = exp(lnAcceptanceRatio); // Accept or reject the move double u = GLOBAL_RNG->uniform01(); if (u < r) { numAccepted++; //keep tree->keep(); rates[nodeIdx]->keep(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->keep(); } } else { reject(); tree->restore(); rates[nodeIdx]->restore(); for (size_t i = 0; i < node->getNumberOfChildren(); i++) { size_t childIdx = node->getChild(i).getIndex(); rates[childIdx]->restore(); } } } }
void PhyloBrownianProcessREML::recursiveComputeLnProbability( const TopologyNode &node, size_t nodeIndex ) { // check for recomputation if ( node.isTip() == false && dirtyNodes[nodeIndex] ) { // mark as computed dirtyNodes[nodeIndex] = false; std::vector<double> &p_node = this->partialLikelihoods[this->activeLikelihood[nodeIndex]][nodeIndex]; std::vector<double> &mu_node = this->contrasts[this->activeLikelihood[nodeIndex]][nodeIndex]; // get the number of children size_t num_children = node.getNumberOfChildren(); for (size_t j = 1; j < num_children; ++j) { size_t leftIndex = nodeIndex; const TopologyNode *left = &node; if ( j == 1 ) { left = &node.getChild(0); leftIndex = left->getIndex(); recursiveComputeLnProbability( *left, leftIndex ); } const TopologyNode &right = node.getChild(j); size_t rightIndex = right.getIndex(); recursiveComputeLnProbability( right, rightIndex ); const std::vector<double> &p_left = this->partialLikelihoods[this->activeLikelihood[leftIndex]][leftIndex]; const std::vector<double> &p_right = this->partialLikelihoods[this->activeLikelihood[rightIndex]][rightIndex]; // get the per node and site contrasts const std::vector<double> &mu_left = this->contrasts[this->activeLikelihood[leftIndex]][leftIndex]; const std::vector<double> &mu_right = this->contrasts[this->activeLikelihood[rightIndex]][rightIndex]; // get the propagated uncertainties double delta_left = this->contrastUncertainty[this->activeLikelihood[leftIndex]][leftIndex]; double delta_right = this->contrastUncertainty[this->activeLikelihood[rightIndex]][rightIndex]; // get the scaled branch lengths double v_left = 0; if ( j == 1 ) { v_left = this->computeBranchTime(leftIndex, left->getBranchLength()); } double v_right = this->computeBranchTime(rightIndex, right.getBranchLength()); // add the propagated uncertainty to the branch lengths double t_left = v_left + delta_left; double t_right = v_right + delta_right; // set delta_node = (t_l*t_r)/(t_l+t_r); this->contrastUncertainty[this->activeLikelihood[nodeIndex]][nodeIndex] = (t_left*t_right) / (t_left+t_right); double stdev = sqrt(t_left+t_right); for (int i=0; i<this->numSites; i++) { mu_node[i] = (mu_left[i]*t_right + mu_right[i]*t_left) / (t_left+t_right); // get the site specific rate of evolution double standDev = this->computeSiteRate(i) * stdev; // compute the contrasts for this site and node double contrast = mu_left[i] - mu_right[i]; // compute the probability for the contrasts at this node double lnl_node = RbStatistics::Normal::lnPdf(0, standDev, contrast); // sum up the probabilities of the contrasts p_node[i] = lnl_node + p_left[i] + p_right[i]; } // end for-loop over all sites } // end for-loop over all children } // end if we need to compute something for this node. }