1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270
|
//
// partitionmodelplen.cpp
// iqtree
//
// Created by Olga on 04/05/17.
//
//
#include <stdio.h>
#include "partitionmodelplen.h"
#include "utils/timeutil.h"
/**********************************************************
* class PartitionModelPlen
**********************************************************/
//const double MIN_GENE_RATE = 0.001;
//const double MAX_GENE_RATE = 1000.0;
//const double TOL_GENE_RATE = 0.0001;
PartitionModelPlen::PartitionModelPlen()
: PartitionModel()
{
// optimizing_part = -1;
}
PartitionModelPlen::PartitionModelPlen(Params ¶ms, PhyloSuperTreePlen *tree, ModelsBlock *models_block)
: PartitionModel(params, tree, models_block)
{
// optimizing_part = -1;
}
PartitionModelPlen::~PartitionModelPlen()
{
}
void PartitionModelPlen::startCheckpoint() {
checkpoint->startStruct("PartitionModelPlen");
}
void PartitionModelPlen::saveCheckpoint() {
startCheckpoint();
PhyloSuperTreePlen *tree = (PhyloSuperTreePlen*)site_rate->getTree();
if (!tree->fixed_rates) {
int nrates = tree->part_info.size();
double *part_rates = new double[nrates];
for (int i = 0; i < nrates; i++)
part_rates[i] = tree->part_info[i].part_rate;
CKP_ARRAY_SAVE(nrates, part_rates);
delete [] part_rates;
}
endCheckpoint();
PartitionModel::saveCheckpoint();
}
void PartitionModelPlen::restoreCheckpoint() {
startCheckpoint();
PhyloSuperTreePlen *tree = (PhyloSuperTreePlen*)site_rate->getTree();
if (!tree->fixed_rates) {
int nrates = tree->part_info.size();
double *part_rates = new double[nrates];
if (CKP_ARRAY_RESTORE(nrates, part_rates)) {
for (int i = 0; i < nrates; i++)
tree->part_info[i].part_rate = part_rates[i];
tree->mapTrees();
}
delete [] part_rates;
}
endCheckpoint();
PartitionModel::restoreCheckpoint();
}
double PartitionModelPlen::optimizeParameters(int fixed_len, bool write_info, double logl_epsilon, double gradient_epsilon) {
PhyloSuperTreePlen *tree = (PhyloSuperTreePlen*)site_rate->getTree();
double tree_lh = 0.0, cur_lh = 0.0;
int ntrees = tree->size();
//tree->initPartitionInfo(); // FOR OLGA: needed here
for(int part = 0; part < ntrees; part++){
tree->part_info[part].cur_score = 0.0;
}
if (fixed_len == BRLEN_OPTIMIZE) {
tree_lh = tree->optimizeAllBranches(1);
} else {
tree_lh = tree->computeLikelihood();
}
cout<<"Initial log-likelihood: "<<tree_lh<<endl;
double begin_time = getRealTime();
int i;
for(i = 1; i < tree->params->num_param_iterations; i++){
cur_lh = 0.0;
if (tree->part_order.empty()) tree->computePartitionOrder();
#ifdef _OPENMP
#pragma omp parallel for reduction(+: cur_lh) schedule(dynamic) if(tree->num_threads > 1)
#endif
for (int partid = 0; partid < ntrees; partid++) {
int part = tree->part_order[partid];
// Subtree model parameters optimization
tree->part_info[part].cur_score = tree->at(part)->getModelFactory()->optimizeParametersOnly(i+1,
gradient_epsilon/min(min(i,ntrees),10), tree->part_info[part].cur_score);
if (tree->part_info[part].cur_score == 0.0)
tree->part_info[part].cur_score = tree->at(part)->computeLikelihood();
cur_lh += tree->part_info[part].cur_score;
// normalize rates s.t. branch lengths are #subst per site
double mean_rate = tree->at(part)->getRate()->rescaleRates();
if (fabs(mean_rate-1.0) > 1e-6) {
if (tree->fixed_rates) {
outError("Unsupported -spj. Please use proportion edge-linked partition model (-spp)");
}
tree->at(part)->scaleLength(mean_rate);
tree->part_info[part].part_rate *= mean_rate;
}
}
if (tree->params->link_alpha) {
cur_lh = optimizeLinkedAlpha(write_info, gradient_epsilon);
}
if (verbose_mode >= VB_MED)
cout << "LnL after optimizing individual models: " << cur_lh << endl;
if (cur_lh <= tree_lh - 1.0) {
// more info for ASSERTION
writeInfo(cout);
tree->printTree(cout, WT_BR_LEN+WT_NEWLINE);
}
ASSERT(cur_lh > tree_lh - 1.0 && "individual model opt reduces LnL");
tree->clearAllPartialLH();
// Optimizing gene rate
if(!tree->fixed_rates){
cur_lh = optimizeGeneRate(gradient_epsilon);
if (verbose_mode >= VB_MED) {
cout << "LnL after optimizing partition-specific rates: " << cur_lh << endl;
writeInfo(cout);
}
ASSERT(cur_lh > tree_lh - 1.0 && "partition rate opt reduces LnL");
}
// Optimizing branch lengths
int my_iter = min(5,i+1);
if (fixed_len == BRLEN_OPTIMIZE){
double new_lh = tree->optimizeAllBranches(my_iter, logl_epsilon);
ASSERT(new_lh > cur_lh - 1.0);
cur_lh = new_lh;
} else if (fixed_len == BRLEN_SCALE) {
double scaling = 1.0;
double new_lh = tree->optimizeTreeLengthScaling(MIN_BRLEN_SCALE, scaling, MAX_BRLEN_SCALE, gradient_epsilon);
ASSERT(new_lh > cur_lh - 1.0);
cur_lh = new_lh;
}
cout<<"Current log-likelihood at step "<<i<<": "<<cur_lh<<endl;
if(fabs(cur_lh-tree_lh) < logl_epsilon) {
tree_lh = cur_lh;
break;
}
// make sure that the new logl is not so bad compared with previous logl
ASSERT(cur_lh > tree_lh - 1.0 && "branch length opt reduces LnL");
tree_lh = cur_lh;
}
// cout <<"OPTIMIZE MODEL has finished"<< endl;
if (write_info)
writeInfo(cout);
cout << "Parameters optimization took " << i-1 << " rounds (" << getRealTime()-begin_time << " sec)" << endl << endl;
return tree_lh;
}
double PartitionModelPlen::optimizeParametersGammaInvar(int fixed_len, bool write_info, double logl_epsilon, double gradient_epsilon) {
outError("Thorough I+G parameter optimization does not work with edge-linked partition model yet");
return 0.0;
}
void PartitionModelPlen::writeInfo(ostream &out) {
PhyloSuperTreePlen *tree = (PhyloSuperTreePlen*)site_rate->getTree();
int ntrees = tree->size();
if (!tree->fixed_rates) {
out << "Partition-specific rates: ";
for(int part = 0; part < ntrees; part++){
out << " " << tree->part_info[part].part_rate;
}
out << endl;
}
}
double PartitionModelPlen::optimizeGeneRate(double gradient_epsilon)
{
PhyloSuperTreePlen *tree = (PhyloSuperTreePlen*)site_rate->getTree();
// BQM 22-05-2015: change to optimize individual rates
int i;
double score = 0.0;
double nsites = tree->getAlnNSite();
vector<DoubleVector> brlen;
brlen.resize(tree->branchNum);
tree->getBranchLengths(brlen);
double max_brlen = 0.0;
for (i = 0; i < brlen.size(); i++)
for (int j = 0; j < brlen[i].size(); j++)
if (brlen[i][j] > max_brlen)
max_brlen = brlen[i][j];
if (tree->part_order.empty()) tree->computePartitionOrder();
#ifdef _OPENMP
#pragma omp parallel for reduction(+: score) private(i) schedule(dynamic) if(tree->num_threads > 1)
#endif
for (int j = 0; j < tree->size(); j++) {
int i = tree->part_order[j];
double min_scaling = 1.0/tree->at(i)->getAlnNSite();
double max_scaling = nsites / tree->at(i)->getAlnNSite();
if (max_scaling < tree->part_info[i].part_rate)
max_scaling = tree->part_info[i].part_rate;
if (min_scaling > tree->part_info[i].part_rate)
min_scaling = tree->part_info[i].part_rate;
tree->part_info[i].cur_score = tree->at(i)->optimizeTreeLengthScaling(min_scaling, tree->part_info[i].part_rate, max_scaling, gradient_epsilon);
score += tree->part_info[i].cur_score;
}
// now normalize the rates
double sum = 0.0;
size_t nsite = 0;
for (i = 0; i < tree->size(); i++) {
sum += tree->part_info[i].part_rate * tree->at(i)->aln->getNSite();
if (tree->at(i)->aln->seq_type == SEQ_CODON && tree->rescale_codon_brlen)
nsite += 3*tree->at(i)->aln->getNSite();
else
nsite += tree->at(i)->aln->getNSite();
}
sum /= nsite;
if (sum > tree->params->max_branch_length / max_brlen) {
outWarning("Too high (saturated) partition rates for proportional partition model!");
// outWarning("Please switch to the edge-equal partition model via -q option instead of -spp");
// exit(EXIT_FAILURE);
}
tree->scaleLength(sum);
sum = 1.0/sum;
for (i = 0; i < tree->size(); i++)
tree->part_info[i].part_rate *= sum;
return score;
}
int PartitionModelPlen::getNParameters(int brlen_type) {
PhyloSuperTreePlen *tree = (PhyloSuperTreePlen*)site_rate->getTree();
int df = 0;
for (PhyloSuperTreePlen::iterator it = tree->begin(); it != tree->end(); it++) {
df += (*it)->getModelFactory()->model->getNDim() +
(*it)->getModelFactory()->model->getNDimFreq() +
(*it)->getModelFactory()->site_rate->getNDim();
}
df += tree->branchNum;
if(!tree->fixed_rates)
df += tree->size()-1;
if (linked_alpha > 0.0)
df ++;
return df;
}
int PartitionModelPlen::getNDim(){
PhyloSuperTreePlen *tree = (PhyloSuperTreePlen*)site_rate->getTree();
int ndim = tree->size() -1;
return ndim;
}
|