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 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421
|
/* RAxML-VI-HPC (version 2.2) a program for sequential and parallel estimation of phylogenetic trees
* Copyright August 2006 by Alexandros Stamatakis
*
* Partially derived from
* fastDNAml, a program for estimation of phylogenetic trees from sequences by Gary J. Olsen
*
* and
*
* Programs of the PHYLIP package by Joe Felsenstein.
*
* This program is free software; you may redistribute it and/or modify its
* under the terms of the GNU General Public License as published by the Free
* Software Foundation; either version 2 of the License, or (at your option)
* any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
* or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
* for more details.
*
*
* For any other enquiries send an Email to Alexandros Stamatakis
* Alexandros.Stamatakis@epfl.ch
*
* When publishing work that is based on the results from RAxML-VI-HPC please cite:
*
* Alexandros Stamatakis:"RAxML-VI-HPC: maximum likelihood-based phylogenetic analyses with thousands
* of taxa and mixed models".
* Bioinformatics 2006; doi: 10.1093/bioinformatics/btl446
*/
#ifndef WIN32
#include <unistd.h>
#endif
#include <math.h>
#include <time.h>
#include <stdlib.h>
#include <stdio.h>
#include <ctype.h>
#include <string.h>
#include "axml.h"
extern char run_id[128];
extern char workdir[1024];
extern double masterTime;
/*
function below not very interesting, standard RAxML stuff
*/
static double testInsertThorough(tree *tr, nodeptr r, nodeptr q, boolean useVector)
{
double
result,
qz[NUM_BRANCHES],
z[NUM_BRANCHES];
nodeptr
x = q->back,
s = r->back;
int
j;
for(j = 0; j < tr->numBranches; j++)
{
qz[j] = q->z[j];
z[j] = sqrt(qz[j]);
if(z[j] < zmin)
z[j] = zmin;
if(z[j] > zmax)
z[j] = zmax;
}
hookup(r->next, q, z, tr->numBranches);
hookup(r->next->next, x, z, tr->numBranches);
hookupDefault(r, s, tr->numBranches);
newviewGeneric(tr, r);
localSmooth(tr, r, smoothings);
if(useVector)
result = evaluateGenericVector(tr, r);
else
result = evaluateGeneric(tr, r);
hookup(q, x, qz, tr->numBranches);
r->next->next->back = r->next->back = (nodeptr) NULL;
return result;
}
/*
structure to store likelihood and insertion node
for each window position
*/
typedef struct
{
double lh;
nodeptr p;
} positionData;
/*
traverse the tree recursively and insert taxon into each position
*/
static void traverseBias(nodeptr p, nodeptr q, tree *tr, int *branchCounter, positionData *pd, int windowSize)
{
double
sum = 0.0,
result = 0.0;
int
i;
/*
compute the likelihood of inserting the tip attached to
p between q and q->back
Actually in testInsertThorough() we compute the per site
log likes which are then stored in an array tr->perSiteLL[i]
*/
result = testInsertThorough(tr, p, q, TRUE);
/*
stuff below can be removed at some point
it just makes me feel better
*/
for(i = 0; i < tr->cdta->endsite; i++)
sum += tr->perSiteLL[i];
assert(ABS(sum - result) < 0.001);
/*************************************/
/*
for each window position just compute the likelihood over
the window.
If its is better than the current best one, store the likelihood
and the insertion node in the data structure
*/
for(i = 0; i < tr->cdta->endsite - windowSize; i++)
{
int
j;
for(j = i, sum = 0.0; j < i + windowSize; j++)
sum += tr->perSiteLL[j];
if(sum > pd[i].lh)
{
pd[i].lh = sum;
pd[i].p = q;
}
}
*branchCounter = *branchCounter + 1;
/* and here comes the recursion */
if(!isTip(q->number, tr->rdta->numsp))
{
traverseBias(p, q->next->back, tr, branchCounter, pd, windowSize);
traverseBias(p, q->next->next->back, tr, branchCounter, pd, windowSize);
}
}
/*
functions to compute the node distances between inferred and true
placement positions. I think that they are correct.
*/
static int findRec(nodeptr ref, nodeptr best, int mxtips, int value)
{
if(isTip(ref->number, mxtips))
{
if(ref == best || ref->back == best)
return value;
else
return 0;
}
else
{
int
d1,
d2;
if(ref == best || ref->back == best)
return value;
d1 = findRec(ref->next->back, best, mxtips, value + 1);
d2 = findRec(ref->next->next->back, best, mxtips, value + 1);
assert((d1 > 0 && d2 == 0) || (d2 > 0 && d1 == 0) || (d1 == 0 && d2 == 0));
return (d1 + d2);
}
}
static double getNodeDistance(nodeptr ref, nodeptr best, int mxtips)
{
int
d1 = findRec(ref, best, mxtips, 0),
d2 = findRec(ref->back, best, mxtips, 0);
assert((d1 > 0 && d2 == 0) || (d2 > 0 && d1 == 0) || (d1 == 0 && d2 == 0));
return ((double)(d1 + d2));
}
void computePlacementBias(tree *tr, analdef *adef)
{
int
windowSize = adef->slidingWindowSize,
k,
i,
tips,
numTraversalBranches = (2 * (tr->mxtips - 1)) - 3; /* compute number of branches into which we need to insert once we have removed a taxon */
char
fileName[1024];
FILE
*outFile;
/* data for each sliding window starting position */
positionData
*pd = (positionData *)rax_malloc(sizeof(positionData) * (tr->cdta->endsite - windowSize));
double
*nodeDistances = (double*)rax_calloc(tr->cdta->endsite - windowSize, sizeof(double)), /* array to store node distnces ND for every sliding window position */
*distances = (double*)rax_calloc(tr->cdta->endsite, sizeof(double)); /* array to store avg distances for every site */
strcpy(fileName, workdir);
strcat(fileName, "RAxML_SiteSpecificPlacementBias.");
strcat(fileName, run_id);
outFile = myfopen(fileName, "w");
printBothOpen("Likelihood of comprehensive tree %f\n\n", tr->likelihood);
if(windowSize > tr->cdta->endsite)
{
printBothOpen("The size of your sliding window is %d while the number of sites in the alignment is %d\n\n", windowSize, tr->cdta->endsite);
exit(-1);
}
if(windowSize >= (int)(0.9 * tr->cdta->endsite))
printBothOpen("WARNING: your sliding window of size %d is only slightly smaller than you alignment that has %d sites\n\n", windowSize, tr->cdta->endsite);
printBothOpen("Sliding window size: %d\n\n", windowSize);
/* prune and re-insert on tip at a time into all branches of the remaining tree */
for(tips = 1; tips <= tr->mxtips; tips++)
{
nodeptr
myStart,
p = tr->nodep[tips]->back, /* this is the node at which we are prunung */
p1 = p->next->back,
p2 = p->next->next->back;
double
pz[NUM_BRANCHES],
p1z[NUM_BRANCHES],
p2z[NUM_BRANCHES];
int
branchCounter = 0;
/* reset array values for this tip */
for(i = 0; i < tr->cdta->endsite; i++)
{
pd[i].lh = unlikely;
pd[i].p = (nodeptr)NULL;
}
/* store the three branch lengths adjacent to the position at which we prune */
for(i = 0; i < tr->numBranches; i++)
{
p1z[i] = p1->z[i];
p2z[i] = p2->z[i];
pz[i] = p->z[i];
}
/* prune the taxon, optimizing the branch between p1 and p2 */
removeNodeBIG(tr, p, tr->numBranches);
printBothOpen("Pruning taxon Number %d [%s]\n", tips, tr->nameList[tips]);
/* find any tip to start traversing the tree */
myStart = findAnyTip(p1, tr->mxtips);
/* insert taxon, compute likelihood and remove taxon again from all branches */
traverseBias(p, myStart->back, tr, &branchCounter, pd, windowSize);
assert(branchCounter == numTraversalBranches);
/* for every sliding window position calc ND to the true/correct position at p */
for(i = 0; i < tr->cdta->endsite - windowSize; i++)
nodeDistances[i] = getNodeDistance(p1, pd[i].p, tr->mxtips);
/* now analyze */
for(i = 0; i < tr->cdta->endsite; i++)
{
double
d = 0.0;
int
s = 0;
/*
check site position, i.e., doe we have windowSize data points available
or fewer because we are at the start or the end of the alignment
*/
/*
for each site just accumulate the node distances we have for all
sliding windows that passed over this site
*/
if(i < windowSize)
{
for(k = 0; k <= i; k++, s++)
d += nodeDistances[k];
}
else
{
if(i < tr->cdta->endsite - windowSize)
{
for(k = i - windowSize + 1; k <= i; k++, s++)
d += nodeDistances[k];
}
else
{
for(k = i - windowSize; k < (tr->cdta->endsite - windowSize); k++, s++)
d += nodeDistances[k];
}
}
/*
now just divide the accumultaed ND distance by
the number of distances we have for this position and then add it to the acc
distances over all taxa.
I just realized that the version on which I did the tests
I sent to Simon I used
distances[i] = d / ((double)s);
instead of
distances[i] += d / ((double)s);
gamo tin poutana mou
*/
distances[i] += (d / ((double)s));
}
/*
re-connect taxon to its original position
*/
hookup(p->next, p1, p1z, tr->numBranches);
hookup(p->next->next, p2, p2z, tr->numBranches);
hookup(p, p->back, pz, tr->numBranches);
/*
fix likelihood vectors
*/
newviewGeneric(tr, p);
}
/*
now just compute the average ND over all taxa
*/
for(i = 0; i < tr->cdta->endsite; i++)
{
double
avg = distances[i] / ((double)tr->mxtips);
fprintf(outFile, "%d %f\n", i, avg);
}
printBothOpen("\nTime for EPA-based site-specific placement bias calculation: %f\n", gettime() - masterTime);
printBothOpen("Site-specific placement bias statistics written to file %s\n", fileName);
fclose(outFile);
exit(0);
}
|