PageRenderTime 85ms CodeModel.GetById 15ms app.highlight 61ms RepoModel.GetById 1ms app.codeStats 1ms

/archive/source.archive/bm.only/trait.rjmcmc.11172010.R

http://github.com/eastman/auteur
R | 742 lines | 623 code | 83 blank | 36 comment | 113 complexity | 8559d1d1884098c4610b656a2f5f2cf4 MD5 | raw file
  1## Written by LJ HARMON (2008), AL HIPP (2008), and JM EASTMAN (2010)
  2
  3rjMCMC.trait<-function (phy, data, ngen=1000, sampleFreq=100, probMergeSplit=0.05, probRoot=0.01, lambdaK=log(2), constrainK=FALSE, jumpsize=2, fileBase="result", simplestart=FALSE) 
  4{
  5	model="BM"
  6	heat=1
  7#	require(ouch)
  8	require(geiger)
  9		
 10### prepare data for rjMCMC
 11	cur.model		<- model								
 12
 13	dataList		<- prepare.data(phy, data)			
 14	ape.tre			<- cur.tre <- dataList$ape.tre			
 15	orig.dat		<- dataList$orig.dat						
 16	nn				<- length(ape.tre$edge.length)	
 17	node.des		<- sapply(unique(c(ape.tre$edge[1,1],ape.tre$edge[,2])), function(x) get.descendants.of.node(x, ape.tre))
 18	names(node.des) <- c(ape.tre$edge[1,1], unique(ape.tre$edge[,2]))
 19	subtrees		<- subtrees(phy)
 20	subtrees.marker <- sapply(subtrees, function(x) return(min(x$node.label)))
 21
 22# initialize parameters
 23	if(is.numeric(constrainK) & (constrainK > length(ape.tre$edge) | constrainK < 1)) stop("Constraint on rate shifts is nonsensical. Ensure that constrainK is at least 1 and less than the number of available nodes in the tree.")
 24
 25	if(simplestart | is.numeric(constrainK) ) {
 26		if(is.numeric(constrainK)) {
 27			init.rate	<- generate.starting.point(data, ape.tre, node.des, theta=FALSE, K=constrainK, jumpsize=jumpsize)
 28		} else {
 29			init.rate	<- list(values=rep(fitContinuous(phy,data)$Trait1$beta,length(phy$edge.length)),delta=rep(0,length(phy$edge.length)))
 30		}
 31	} else {
 32		init.rate	<- generate.starting.point(data, ape.tre, node.des, theta=FALSE, K=constrainK, jumpsize=jumpsize )
 33	}
 34	
 35    cur.rates		<- init.rate$values			
 36	cur.delta.rates	<- init.rate$delta
 37	cur.root		<- adjust.value(mean(orig.dat), jumpsize)
 38	cur.tre			<- apply.BMM.to.apetree(ape.tre, cur.rates)
 39	
 40	mod.cur = new.bm.lik.fn(cur.rates, cur.root, cur.tre, orig.dat)
 41	
 42# proposal counts
 43	nRateProp =		0	
 44	nRateSwapProp = 0									
 45	nRcatProp =		0										
 46	nRootProp =		0
 47	nRateOK =		0											
 48	nRateSwapOK =	0
 49	nRcatOK =		0											
 50	nRootOK =		0
 51
 52	cOK=c(												
 53		nRateOK,
 54		nRateSwapOK,
 55		nRcatOK, 
 56		nRootOK 
 57	)
 58		
 59	tickerFreq=ceiling(ngen/30)
 60	
 61# file handling
 62	parmBase=paste(model, fileBase, "parameters/",sep=".")
 63	if(!file.exists(parmBase)) dir.create(parmBase)
 64	parlogger(model=model, init=TRUE, node.des, parmBase=parmBase)
 65	errorLog=file(paste(parmBase,paste(cur.model, fileBase, "rjTraitErrors.log",sep="."),sep="/"),open='w+')
 66	generate.error.message(i=NULL, mod.cur=NULL, mod.new=NULL, lnR=NULL, errorLog=errorLog, init=TRUE)
 67	runLog=file(paste(parmBase,paste(cur.model, fileBase, "rjTrait.log",sep="."),sep="/"),open='w+')
 68	generate.log(bundled.parms=NULL, cur.model, file=runLog, init=TRUE)
 69
 70
 71### Begin rjMCMC
 72    for (i in 1:ngen) {
 73        lnLikelihoodRatio <- lnHastingsRatio <- lnPriorRatio <- 0
 74		startparms = c(nRateProp, nRateSwapProp, nRcatProp, nRootProp)
 75
 76## BM IMPLEMENTATION ##
 77		while(1) {
 78			if (runif(1) < (2 * probMergeSplit) & !constrainK) {			# adjust rate categories
 79				nr=split.or.merge(cur.delta.rates, cur.rates, ape.tre, node.des, lambdaK)
 80				new.rates=nr$new.values
 81				new.delta.rates=nr$new.delta
 82				new.root=cur.root
 83				nRcatProp=nRcatProp+1
 84				lnHastingsRatio=nr$lnHastingsRatio
 85				lnPriorRatio=nr$lnPriorRatio
 86				break()
 87			} else { 
 88				if(runif(1)<probRoot) {										# adjust root
 89					new.root=adjust.value(cur.root, jumpsize)
 90					new.rates=cur.rates
 91					new.delta.rates=cur.delta.rates
 92					nRootProp=nRootProp+1
 93					break()
 94				} else {													
 95					if(runif(1)>0.1 & length(unique(cur.rates))>1) {		# swap rates classes
 96						nr=switcheroo(cur.rates)
 97						new.rates=nr$new.values ##
 98						new.delta.rates=cur.delta.rates 
 99						new.root=cur.root
100						nRateSwapProp=nRateSwapProp+1
101						break()
102					} else {												# adjust rates
103						new.rates=adjust.rate(cur.rates, jumpsize)
104						new.delta.rates=cur.delta.rates
105						new.root=cur.root
106						nRateProp=nRateProp+1
107						break()
108					} 
109				}
110			}
111		}
112		
113		new.tre=apply.BMM.to.apetree(ape.tre, new.rates)
114		mod.new=try(new.bm.lik.fn(new.rates, new.root, new.tre, orig.dat),silent=TRUE)
115		
116		if(inherits(mod.new, "try-error")) {mod.new=as.list(mod.new); mod.new$lnL=Inf}
117		if(!is.infinite(mod.new$lnL)) {
118			lnLikelihoodRatio = mod.new$lnL - mod.cur$lnL
119		} else {
120			lnLikelihoodRatio = -Inf
121			failures=paste("./failed", model, "parameters/",sep=".")
122			failed.trees=paste(failures, "trees", sep="/")
123			failed.plots=paste(failures, "pdf", sep="/")
124			lapply(list(failures, failed.trees, failed.plots), function(x) if(!file.exists(x)) dir.create(x))
125				
126			pdf(file=paste(failed.plots, paste(i,fileBase,"failed.model.pdf",sep="."),sep="/"))
127				plot(new.tre, tip=FALSE)
128			dev.off()
129			
130			write.tree(new.tre, file=paste(failed.trees, paste(i,fileBase,"tre",sep="."),sep="/"))
131			write.table(matrix(c(i, new.rates),nrow=1),file=paste(failures,"failed.rates.txt",sep="/"),append=TRUE,col=FALSE,row=FALSE,quote=FALSE)
132			write.table(matrix(c(i, new.root),nrow=1),file=paste(failures,"failed.root.txt",sep="/"),append=TRUE,col=FALSE,row=FALSE,quote=FALSE)
133		}
134		
135	# compare likelihoods
136		endparms = c(nRateProp=nRateProp, nRateSwapProp=nRateSwapProp, nRcatProp=nRcatProp, nRootProp=nRootProp)
137
138		r=assess.lnR((heat * lnLikelihoodRatio + heat * lnPriorRatio + lnHastingsRatio)->lnR)
139
140		if (runif(1) <= r$r) {			# adopt proposal
141			decision="adopt"
142			cur.root <- new.root
143			cur.rates <- new.rates
144			cur.delta.rates <- new.delta.rates
145			curr.lnL <- mod.new$lnL
146			cur.tre <- new.tre
147			mod.cur <- mod.new
148			cOK <- determine.accepted.proposal(startparms, endparms, cOK)
149		} else {						# deny proposal	
150			decision="reject"
151			curr.lnL <- mod.cur$lnL 
152		}
153		
154		if(is.infinite(curr.lnL)) stop("starting point has exceptionally poor likelihood")
155
156		if(r$error) generate.error.message(i, mod.cur, mod.new, lnR, errorLog)
157		
158		if(i%%tickerFreq==0) {
159			if(i==tickerFreq) cat("|",rep(" ",9),toupper("generations complete"),rep(" ",9),"|","\n")
160			cat(". ")
161		}
162		
163		if(i%%sampleFreq==0) {
164			bundled.parms=list(gen=i, mrate=exp(mean(log(cur.rates))), cats=sum(cur.delta.rates)+1, root=cur.root, lnL=curr.lnL)
165			generate.log(bundled.parms, cur.model, file=runLog)			
166			parlogger(model=model, init=FALSE, node.des=node.des, i=i, curr.lnL=curr.lnL, cur.root=cur.root, cur.rates=cur.rates, cur.delta.rates=cur.delta.rates, parmBase=parmBase)
167		}
168	}
169	
170# End rjMCMC
171	close(errorLog)
172	close(runLog)
173	summarize.run(cOK, endparms, cur.model)	
174}
175
176
177## AUXILIARY FUNCTIONS
178
179new.bm.lik.fn <- function(rates, root, ape.tre, orig.dat) { # from LJ HARMON	
180# mod 10.18.2010 JM Eastman: using determinant()$modulus rather than det() to stay in log-space
181	tree=ape.tre
182	y=orig.dat
183	
184	b <- vcv.phylo(ape.tre)
185	w <- rep(root, nrow(b))
186	num <- -t(y-w)%*%solve(b)%*%(y-w)/2
187	den <- 0.5*(length(y)*log(2*pi) + as.numeric(determinant(b)$modulus))
188	list(
189		 root=root,
190		 lnL = (num-den)[1,1]
191		 )	
192}
193
194split.or.merge <- function(cur.delta, cur.values, phy, node.des, lambda=lambda) { 
195# mod 11.17.2010 JM Eastman
196#	cur.delta=cur.delta.rates; cur.values=cur.rates; theta=FALSE; lambda=log(2); jumpsize=2
197	bb=cur.delta
198	vv=cur.values
199	names(vv)<-names(bb)<-phy$edge[,2]
200	new.bb=choose.one(bb)
201	new.vv=vv
202	
203	s=names(new.bb[bb!=new.bb])
204	all.shifts=as.numeric(names(bb[bb>0]))
205	all.D=node.des[[which(names(node.des)==s)]]
206	if(length(all.D)!=0) {
207		untouchable=unlist(lapply(all.shifts[all.shifts>s], function(x)node.des[[which(names(node.des)==x)]]))
208		remain.unchanged=union(all.shifts, untouchable)
209	} else {
210		untouchable=NULL
211		remain.unchanged=list()
212	}
213	
214	marker=match(s, names(new.vv))
215	nn=length(vv)
216	K=sum(bb)+1
217	N=Ntip(phy)
218	
219	ncat=sum(bb)
220	cur.vv=as.numeric(vv[marker])
221	ca.vv=length(which(vv==cur.vv))
222	
223	if(sum(new.bb)>sum(bb)) {			# add transition: SPLIT
224#		print("s")
225
226		n.desc=sum(!all.D%in%remain.unchanged)+1
227		n.split=sum(vv==cur.vv)-n.desc
228		if(runif(1)<0.5) {
229			u=runif(1,min=-0.5*n.desc*cur.vv, max=0.5*n.split*cur.vv)
230			nr.desc=cur.vv + u/n.desc
231			nr.split=cur.vv - u/n.split
232#			(n.desc*nr.desc+n.split*nr.split)/(n.desc+n.split); nr.desc; nr.split
233		} else {
234			u=runif(1,min=-0.5*n.split*cur.vv, max=0.5*n.desc*cur.vv)
235			nr.desc=cur.vv - u/n.desc
236			nr.split=cur.vv + u/n.split
237#			(n.desc*nr.desc+n.split*nr.split)/(n.desc+n.split); nr.desc; nr.split
238		}
239		new.vv[vv==cur.vv]=nr.split
240		if(length(remain.unchanged)==0) {	# assign new value to all descendants
241			new.vv[match(all.D,names(new.vv))] = nr.desc
242		} else {							# assign new value to all 'open' descendants 
243			new.vv[match(all.D[!(all.D%in%remain.unchanged)],names(new.vv))] = nr.desc
244		}
245		new.vv[match(s, names(new.vv))]=nr.desc
246		lnHastingsRatio = log((K+1)/(2*N-2-K)) ### from Drummond and Suchard 2010: where N is tips, K is number of local parms in tree
247		lnPriorRatio = log(ptpois(K+1,lambda,nn)/ptpois(K,lambda,nn))
248		
249	} else {							# drop transition: MERGE
250#		print("m")
251		
252		anc = get.ancestor.of.node(s, phy)
253		if(!is.root(anc, phy)) {			# base new rate on ancestral rate of selected branch
254			anc.vv=as.numeric(vv[match(anc,names(vv))])
255			na.vv=length(which(vv==anc.vv))
256			nr=(anc.vv*na.vv+cur.vv*ca.vv)/(ca.vv+na.vv)
257			new.vv[vv==cur.vv | vv==anc.vv]=nr
258		} else {							# if ancestor of selected node is root, base new rate on sister node
259			sister.tmp=get.desc.of.node(anc,phy)
260			sister=sister.tmp[sister.tmp!=s]
261			sis.vv=as.numeric(vv[match(sister,names(vv))])
262			ns.vv=length(which(vv==sis.vv))
263			nr=(sis.vv*ns.vv+cur.vv*ca.vv)/(ca.vv+ns.vv)
264			new.vv[vv==cur.vv | vv==sis.vv]=nr			
265		}
266		lnHastingsRatio = log((2*N-2-K+1)/K) ### from Drummond and Suchard 2010: where N is tips, K is number of local parms in tree
267		lnPriorRatio = log(ptpois(K-1,lambda,nn)/ptpois(K,lambda,nn))
268		
269	}
270	
271	return(list(new.delta=new.bb, new.values=new.vv, lnHastingsRatio=lnHastingsRatio, lnPriorRatio=lnPriorRatio))
272}
273
274switcheroo <- function(cur.values) { 
275# mod 11.17.2010 JM Eastman
276#	cur.delta=cur.delta.rates; cur.values=cur.rates; theta=FALSE; lambda=log(2); jumpsize=2
277
278	vv=cur.values
279	rates.sample=unique(vv)[sample(1:length(unique(vv)), 2)]
280	new.vv=vv
281	new.vv[which(vv==rates.sample[1])]=rates.sample[2]
282	new.vv[which(vv==rates.sample[2])]=rates.sample[1]
283	
284	return(list(new.values=new.vv))
285}
286
287generate.starting.point <- function(data, phy, node.des, theta=FALSE, K=FALSE, jumpsize) { 
288# updated JME 11.14.2010
289	nn=length(phy$edge.length)
290	ntip=Ntip(phy)
291	if(!K) nshifts=rtpois(1,log(ntip),nn) else nshifts=K-1
292	if(nshifts!=0) bb=sample(c(rep(1,nshifts),rep(0,nn-nshifts)),replace=FALSE) else bb=rep(0,nn)
293	names(bb)=phy$edge[,2]
294	shifts=as.numeric(names(bb)[bb==1])
295	if(theta) {
296		min.max=c(adjust.value(min(data), jumpsize), adjust.value(max(data), jumpsize))
297		values=runif(sum(bb)+1, min=min.max[min.max==min(min.max)], max=min.max[min.max==max(min.max)])
298	} else {
299		init.rate=fitContinuous(phy,data)[[1]]$beta
300		min.max=c(adjust.rate(init.rate, jumpsize), adjust.rate(init.rate, jumpsize))
301		values=runif(sum(bb)+1, min=min.max[min.max==min(min.max)], max=min.max[min.max==max(min.max)])
302	}
303	internal.shifts<-tip.shifts<-numeric(0)
304	internal.shifts=sort(shifts[shifts>ntip])
305	tip.shifts=shifts[shifts<=ntip]
306	
307	if(length(internal.shifts)==0 & length(tip.shifts)==0) {
308		vv=rep(values, nn)
309	} else {
310		vv=bb
311		vv[]=values[length(values)]
312		i=0
313		if(length(internal.shifts)!=0) {
314			for(i in 1:length(internal.shifts)) {
315				d=node.des[which(names(node.des)==internal.shifts[i])]
316				vv[match(c(internal.shifts[i], unlist(d)), names(vv))]=values[i]
317			}
318		}
319		if(length(tip.shifts)!=0) {
320			for(j in 1:length(tip.shifts)) {
321				vv[match(tip.shifts[j], names(vv))]=values[j+i]
322			}
323		}
324	}
325	return(list(delta=unname(bb), values=unname(vv)))
326}
327
328assign.root.theta <- function(cur.theta, phy) {
329	root=phy$edge[1,1]
330	desc=phy$edge[which(phy$edge[,1]==root),2]
331	root.theta=cur.theta[sample(which(phy$edge[,2]==desc),1)]
332	root.theta
333}
334
335get.descendants.of.node <- function(node, phy) {
336	storage=list()
337	if(is.null(node)) node=phy$edge[1,1]
338	if(node%in%phy$edge[,1]) {
339		desc=c(sapply(node, function(x) {return(phy$edge[which(phy$edge[,1]==x),2])}))
340		while(any(desc%in%phy$edge[,1])) {
341			desc.true=desc[which(desc%in%phy$edge[,1])]
342			desc.desc=lapply(desc.true, function(x) return(get.desc.of.node(x, phy)))
343			
344			storage=list(storage,union(desc,desc.desc))
345			
346			desc=unlist(desc.desc)
347		}
348		storage=sort(unique(unlist(union(desc,storage))))
349	}
350	return(storage)
351}
352
353get.desc.of.node <- function(node, phy) {
354	return(phy$edge[which(phy$edge[,1]==node),2])
355}
356
357get.ancestor.of.node<-function(node, phy) {
358	anc=phy$edge[which(phy$edge[,2]==node),1]
359	anc
360}
361
362choose.one <- function(cur.delta)
363{
364	bb=cur.delta
365	s=sample(1:length(bb),1)
366	bb[s]=1-bb[s]
367	bb
368}
369
370is.root<-function(node,phy) {
371	if(node==phy$edge[1,1]) return(TRUE) else return(FALSE)
372}
373
374assess.lnR <- function(lnR) {
375	if(is.na(lnR) || abs(lnR)==Inf) {
376		error=TRUE
377		r=0
378	} else {
379		if(lnR < -20) {
380			r=0 
381		} else if(lnR > 0) {
382			r=1 
383		} else {
384			r=exp(lnR)
385		}
386		error=FALSE
387	}
388	return(list(r=r, error=error))
389}
390
391apply.BMM.to.apetree<-function(apetree, rates) {
392	apetree$edge.length=apetree$edge.length*rates
393	return(apetree)
394}
395
396adjust.value <- function(value, jumpsize) {
397# mod 10.20.2010 JM Eastman
398	vv=value
399	rch <- runif(1, min = -jumpsize, max = jumpsize)
400	return(vv[]+rch)
401}
402
403adjust.rate <- function(rate, jumpsize) {
404# mod 10.20.2010 JM Eastman
405	vv=rate
406	v=log(vv)
407	rch <- runif(1, min = -jumpsize, max = jumpsize)
408	v=exp(v[]+rch)
409	v
410}
411
412prepare.data <- function(phy, data) {
413	td <- treedata(phy, data, sort = T)	
414	return(list(ape.tre=td$phy, orig.dat=td$data[,1]))
415}
416
417summarize.run<-function(cOK, endparms, cur.model) {
418#	end = c(nRateProp, nRateSwapProp, nRcatProp, nRootProp, nAlphaProp, nThetaProp, nThetaSwapProp, nTcatProp)
419
420	if(cur.model=="BM") {
421		names(endparms)<-names(cOK)<-c("rates", "rate.swap", "rate.catg", "root")
422	} else {
423		names(endparms)<-names(cOK)<-c("rates", "rate.swap", "rate.catg", "root", "alpha", "theta", "theta.swap", "theta.catg")
424	}
425	names(endparms)=paste("pr.",names(endparms),sep="")
426	names(cOK)=paste("ok.",names(cOK),sep="")
427	cat("\n")
428	if(cur.model=="BM") {
429		print(endparms[paste("pr.", c("rates", "rate.swap", "rate.catg", "root"),sep="")])
430		print(cOK[paste("ok.", c("rates", "rate.swap", "rate.catg", "root"),sep="")])
431	} else {
432		print(endparms[paste("pr.", c("rates", "rate.swap", "rate.catg", "alpha", "theta", "theta.swap", "theta.catg"),sep="")])
433		print(cOK[paste("ok.", c("rates", "rate.swap", "rate.catg", "alpha", "theta", "theta.swap", "theta.catg"),sep="")])
434	}
435}
436
437generate.log<-function(bundled.parms, cur.model, file, init=FALSE) {
438#	bundled.parms=list(gen=i, mrate=mean(cur.rates), cats=sum(cur.delta.rates)+1, root=cur.root, alpha=cur.alpha, reg=sum(cur.delta.theta)+1, theta=mean(cur.theta), lnL=curr.lnL)
439	res=bundled.parms
440
441	if(init) {
442		if(cur.model=="BM") msg=paste("gen", "model", "mean.rate", "rates", "root", "lnL", sep="\t") else msg=paste("gen", "model", "mean.rate", "rates", "alpha", "regimes", "mean.theta", "lnL", sep="\t")
443	} else {
444		if(cur.model=="BM") {
445			msg<-paste(res$gen, sQuote(cur.model), sprintf("%.3f", res$mrate), res$cats, sprintf("%.3f", res$root), sprintf("%.3f", res$lnL),sep="\t")
446		} else {
447			msg<-paste(res$gen, sQuote(cur.model), sprintf("%.3f", res$mrate), res$cats, sprintf("%.4f", res$alpha), res$reg, sprintf("%.3f", res$theta), sprintf("%.3f", res$lnL),sep="\t")
448		}
449	}
450	write(msg, file=file, append=TRUE)
451}
452
453parlogger<-function(model, init=FALSE, node.des, i, curr.lnL, cur.root=NULL, cur.rates=NULL, cur.delta.rates=NULL, cur.alpha=NULL, cur.theta=NULL, cur.delta.theta=NULL, parmBase) {	
454	if(init) {
455		msg.long=names(node.des[-1])
456		if(model=="BM") {
457			msg.short=paste("gen", "model", "lnL", sep="\t") 
458			parlogs=paste(parmBase, paste(c("summary", "root", "rates", "rate.shifts"), "txt", sep="."), sep="")
459			sapply(parlogs[1], function(x) write.table(msg.short, x, quote=FALSE, col.names=FALSE, row.names=FALSE))
460			sapply(parlogs[2], function(x) write.table(NULL, x, quote=FALSE, col.names=FALSE, row.names=FALSE))
461			sapply(parlogs[3:length(parlogs)], function(x) write.table(paste(msg.long,collapse=" "), x, quote=FALSE, col.names=FALSE, row.names=FALSE))
462		} else {
463			msg.short=paste("gen", "model", "lnL", sep="\t")
464			parlogs=paste(parmBase, paste(c("summary", "root", "alpha", "rates", "rate.shifts", "optima", "optima.shifts"), "txt", sep="."), sep="")
465			sapply(parlogs[1], function(x) write.table(msg.short, x, quote=FALSE, col.names=FALSE, row.names=FALSE))
466			sapply(parlogs[2:3], function(x) write.table(NULL, x, quote=FALSE, col.names=FALSE, row.names=FALSE))
467			sapply(parlogs[4:length(parlogs)], function(x) write.table(paste(msg.long,collapse=" "), x, quote=FALSE, col.names=FALSE, row.names=FALSE))
468		}
469	} else {
470		if(model=="BM") {
471			parlogs=paste(parmBase, paste(c("summary", "root", "rates", "rate.shifts"), "txt", sep="."), sep="")
472			outlist=list(paste(i, model, curr.lnL, sep="\t"), cur.root, cur.rates, cur.delta.rates) 
473			sapply(1:length(outlist), function(x) write.table(paste(outlist[[x]],collapse=" "), parlogs[[x]], quote=FALSE, col.names=FALSE, row.names=FALSE, append=TRUE)) 
474		} else {
475			parlogs=paste(parmBase, paste(c("summary", "root", "alpha", "rates", "rate.shifts", "optima", "optima.shifts"), "txt", sep="."), sep="")
476			outlist=list(paste(i, model, curr.lnL, sep="\t"), cur.root, cur.alpha, cur.rates, cur.delta.rates, cur.theta, cur.delta.theta) 
477			sapply(1:length(outlist), function(x) write.table(paste(outlist[[x]],collapse=" "), parlogs[[x]], quote=FALSE, col.names=FALSE, row.names=FALSE, append=TRUE)) 
478		}
479	}
480}
481
482generate.error.message<-function(i, mod.cur, mod.new, lnR, errorLog, init=FALSE) {
483	if(init) {
484		initmsg = write(paste("gen", "curr.lnL", "new.lnL", "lnR", sep="\t"), file=errorLog)
485	} else {
486		write(paste(i, sprintf("%.3f", mod.cur$lnL), sprintf("%.3f", mod.new$lnL), sprintf("%.3f", lnR), sep="\t"), file=errorLog)
487	}
488}
489
490determine.accepted.proposal<-function(startparms, endparms, cOK) {
491	which(startparms!=endparms)->marker
492	cOK[marker]=cOK[marker]+1
493	cOK
494}
495
496rtpois<-function(N, lambda, k) {
497	p=ppois(k, lambda, lower.tail=FALSE)
498	out=qpois(runif(N, min=0, max=1-p), lambda)
499	out
500}
501				 
502ptpois<-function(x, lambda, k) {
503	p.k=ppois(k, lambda, lower.tail=FALSE)
504	p.x=ppois(x, lambda, lower.tail=FALSE)
505	ptp=p.x/(1-p.k)
506	return(ptp)
507}
508				 
509## SUMMARIZING MCMC -- PLOTTING FUNCTIONS ##
510shifts.plot=function(phy, base.dir, burnin=0, level=0.03, internal.only=FALSE, paint.branches=TRUE, pdf.base=NULL, verbosity.base=NULL, ...) {
511	color.length=21
512	require(ape)
513	oldwd=getwd()
514	setwd(base.dir)
515	files=dir()
516#	if(optima) {
517#		shifts.file=files[grep("optima.shifts.txt", files)]
518#		estimates.file=files[grep("optima.txt", files)]
519#	} else {
520		shifts.file=files[grep("rate.shifts.txt", files)]
521		estimates.file=files[grep("rates.txt", files)]
522#	}
523	if(!is.null(pdf.base)) lab=paste(pdf.base, paste("bi", burnin, sep=""), paste("fq", level, sep=""), sep=".")
524
525	cat("READING shifts...\n")
526	results=read.table(shifts.file)
527	names(results)=as.numeric(results[1,])
528	results=results[-c(1:(burnin+1)),]
529	if(!all(names(results)%in%phy$edge[,2])) stop("Cannot process results: check that tree matches posterior summaries")
530	hits.tmp=apply(results,1,sum)
531	hits=length(hits.tmp[hits.tmp>0])
532	aa.tmp=apply(results, 2, function(x) sum(x))
533	aa.tmp=aa.tmp/hits
534	aa=aa.tmp[aa.tmp>=level]
535	aa=aa[order(aa, decreasing=TRUE)]
536	aa.nodes=aa[which(as.numeric(names(aa))>Ntip(phy))]
537	aa.tips=aa[which(as.numeric(names(aa))<=Ntip(phy))]
538	aa=aa.nodes
539	print(aa.tips)
540	nodes=as.numeric(names(aa))
541	nodes=nodes[nodes>Ntip(phy)]
542	desc=lapply(nodes, function(x) {
543				foo=get.descendants.of.node(x,phy)
544				if(length(foo)) return(phy$tip.label[foo[foo<=Ntip(phy)]]) else return(NULL)
545				}
546				)
547	shifts=apply(results, 1, sum)
548	cat("READING estimates...\n")
549	ests=read.table(estimates.file)
550	names(ests)=ests[1,]
551	ests=ests[-c(1:(burnin+1)),]
552	average.ests.tmp=apply(ests,2,mean)
553	average.ests=average.ests.tmp-median(average.ests.tmp)
554	if(paint.branches) {
555		require(colorspace)
556#		cce=diverge_hcl(color.length, c = c(100, 0), l = c(50, 90), power = 1.1)
557		cce=diverge_hcl(2*color.length, power = 0.5)
558		e.seq=seq(min(average.ests),max(average.ests),length=color.length)
559		color.gamut=seq(-max(abs(e.seq)), max(abs(e.seq)), length=length(cce))
560		colors.branches=sapply(average.ests, function(x) cce[which(min(abs(color.gamut-x))==abs(color.gamut-x))])
561		colors.branches=colors.branches[match(as.numeric(names(colors.branches)), phy$edge[,2])]
562	} else {
563		colors.branches=1
564	}
565	if(length(nodes)) {
566		require(colorspace)
567		cols=sapply(nodes, function(x) {
568						 a=get.ancestor.of.node(x, phy)
569						 comp=ests[,which(names(ests)==a)]
570						 this=ests[,which(names(ests)==x)]
571						 if(!length(comp)) { # dealing with first descendant of root
572						 d=get.desc.of.node(a, phy)
573						 d=d[which(d!=x)]
574						 
575						# find if sister descendant also experienced rate shift 
576						 d.shifts=results[, which(names(results)==d)]
577						 comp=ests[,which(names(ests)==d)]
578						 x.res=sapply(1:length(d.shifts), function(y) {
579									  if(d.shifts[y]==0) {
580									  zz=this[y]-comp[y]
581									  if(zz>0) {
582									  return(1) 
583									  } else {
584									  if(zz<0) return(-1) else return(0)
585									  }
586									  } else {
587									  return(0)
588									  }
589									  }
590									  )
591						 x.res=mean(x.res[x.res!=0])
592						 } else {
593						 yy=this-comp
594						 zz=yy
595						 zz[yy>0]=1
596						 zz[yy<0]=-1
597						 zz[yy==0]=0
598						 x.res=mean(zz[zz!=0])
599						 }
600						 return(x.res)
601						 }
602						 )
603		ccx=diverge_hcl(21, c = c(100, 0), l = c(50, 90), power = 1.1)
604		c.seq=round(seq(-1,1,by=0.1),digits=1)
605		colors.nodes=ccx[match(round(cols,1),c.seq)]
606		names(colors.nodes)=nodes
607		colors.nodes=colors.nodes[which(as.numeric(names(colors.nodes))>Ntip(phy))]
608	} else {
609		colors.nodes=NULL
610		cols=NULL
611	}
612		
613# send to pdf 
614	if(!is.null(pdf.base)) pdf(file=paste(oldwd, paste(lab,"pdf",sep="."),sep="/"))
615	plot(phy, cex=0.1, lwd=0.4, edge.width=0.5, edge.color=colors.branches, ...)
616	nn=(Ntip(phy)+1):max(phy$edge[,2])
617	ll<-cc<-rr<-rep(0,length(nn))
618	ll[nn%in%nodes]=1
619	if(length(nodes)) cc[nn%in%nodes]=aa[match(nn[nn%in%nodes],nodes)]/max(aa) 
620	if(is.null(colors.nodes)) {
621		nodelabels(pch=ifelse(ll==1, 21, NA), bg=ifelse(ll==1, "black", "white"), cex=2*cc)
622	} else {
623		rr[nn%in%nodes]=colors.nodes[order(names(colors.nodes))]
624		nodelabels(pch=ifelse(ll==1, 21, NA), cex=2*cc, col=rr, bg=rr, lwd=0.5)
625		if(length(aa.tips) & !internal.only) {
626			tt=rep(0,Ntip(phy))
627			tt[(1:Ntip(phy))[as.numeric(names(aa.tips))]]=1
628			tt.cex=ifelse(tt==1, aa.tips/max(aa.tmp), 0)
629			tt.col=sapply(names(aa.tips), function(x) {
630						a=get.ancestor.of.node(x, phy)
631						comp=ests[,which(names(ests)==a)]
632						this=ests[,which(names(ests)==x)]
633						if(!length(comp)) { # dealing with first descendant of root
634						d=get.desc.of.node(a, phy)
635						d=d[which(d!=x)]
636						
637						# find if sister descendant also experienced rate shift 
638						d.shifts=results[, which(names(results)==d)]
639						comp=ests[,which(names(ests)==d)]
640						x.res=sapply(1:length(d.shifts), function(y) {
641									 if(d.shifts[y]==0) {
642									 zz=this[y]-comp[y]
643									 if(zz>0) {
644									 return(1) 
645									 } else {
646									 if(zz<0) return(-1) else return(0)
647									 }
648									 } else {
649									 return(0)
650									 }
651									 }
652									 )
653						x.res=mean(x.res[x.res!=0])
654						} else {
655						yy=this-comp
656						zz=yy
657						zz[yy>0]=1
658						zz[yy<0]=-1
659						zz[yy==0]=0
660						x.res=mean(zz[zz!=0])
661						}
662						return(x.res)
663						}
664						)
665			tt.ccx=diverge_hcl(21, c = c(100, 0), l = c(50, 90), power = 1.1)
666			tt.c.seq=round(seq(-1,1,by=0.1),digits=1)
667			tt.colors.nodes=tt.ccx[match(round(tt.col,1),tt.c.seq)]
668			names(tt.colors.nodes)=as.numeric(names(aa.tips))
669			colors.tips.tmp=tt.colors.nodes[order(as.numeric(names(tt.colors.nodes)))]
670			colors.tips=rep(NA, Ntip(phy))
671			colors.tips[(1:Ntip(phy))[as.numeric(names(aa.tips))]]=colors.tips.tmp
672
673			tiplabels(pch=ifelse(tt==1, 21, NA), cex=2*tt.cex, col=colors.tips, bg=colors.tips, lwd=0.5)
674		}
675		legend("bottomleft", title=toupper("direction"), cex=0.5, pt.cex=1, text.col="darkgray", legend = sprintf("%6.1f", rev(c.seq[seq(1,length(c.seq),by=2)])), pch=21, ncol=1, col = "darkgray", pt.bg = rev(ccx[seq(1,length(c.seq),by=2)]), box.lty="blank", border="white")
676		if(length(nodes)) legend("bottomright", title=toupper("frequency"), cex=0.5, pt.cex=2*(seq(min(cc), max(cc), length=6)), text.col="darkgray", legend = sprintf("%10.3f", seq(round(min(aa),2), round(max(aa),2), length=6)), pch=21, ncol=1, col = "darkgray", pt.bg = "white", box.lty="blank", border="white")
677	}
678	if(!is.null(pdf.base)) dev.off()
679		
680	if(length(nodes)) names(cols)=names(aa)
681	setwd(oldwd)
682	names(desc)=paste("node",nodes,sep=".")
683	if(!is.null(verbosity.base)) {
684		summary.table=c(hits/nrow(results), mean(shifts), median(shifts))
685		names(summary.table)=c("shift.prob", "mean.shift","median.shifts")
686		write.table(summary.table, file=paste(verbosity.base, "run.summary.txt", sep="."), quote=FALSE, col.names=FALSE)
687		if(length(nodes)) {
688			for(nn in 1:length(nodes)) {
689				write.table(c(nodes[nn], desc[[nn]]), file=paste(verbosity.base, "shift.descendants.txt", sep="."), append=TRUE, quote=FALSE, row.names=FALSE, col.names=FALSE)
690			}
691		}
692		allres=merge(as.data.frame(aa),as.data.frame(cols),by=0)
693		allres=allres[order(allres$aa, decreasing=TRUE),]
694		names(allres)=c("node", "conditional.shift.probability", "relative.shift.direction")
695		write.table(allres, file=paste(verbosity.base, "estimates.summary.txt", sep="."), quote=FALSE, row.names=FALSE)
696	}
697	return(list(desc=desc, shift.prob=hits/nrow(results), mean.shifts=mean(shifts), median.shifts=median(shifts), shift.prop=aa, shift.direction=cols, mean.rates=average.ests.tmp, mean.shifts=aa.tmp/nrow(results)))
698}
699
700
701plot.quant<-function(phy, data, data.names, shrinker=1, relative=TRUE, cex=0.5, adj=c(15,5), lwd=0.5, font=1, ...){
702	std=function(x,max,min){return((x-min)/(max-min))}
703	if(relative) {
704		data=sapply(data, function(x) std(x, max(data), min(data)))
705	}
706	
707	f=match(phy$tip.label, data.names)
708	if(is.na(sum(f))) {
709		stop("tips cannot be matched to tree")
710	}
711	else {
712		data=data[f]
713		phy$trait=data
714		plot.phylo(phy, show.tip.label=T, label.offset=adj[1], cex=cex, font=font, ...)
715		tiplabels(text=NULL, pch=21, adj=adj[2], cex=c(phy$trait/shrinker), bg="white",lwd=lwd)
716		} 
717}
718
719branchcol.plot=function(phy, cur.rates, color.length=4, ...) {
720	color.length=color.length
721	require(ape)
722	require(colorspace)
723	ests=cur.rates
724	names(ests)=phy$edge[,2]
725	cce=diverge_hcl(color.length, power = 0.5)
726	e.seq=seq(min(ests),max(ests),length=color.length)
727	colors.branches=sapply(ests, function(x) cce[which(min(abs(e.seq-x))==abs(e.seq-x))])
728	names(colors.branches)=names(ests)
729	colors.branches=colors.branches[match(as.numeric(names(colors.branches)), phy$edge[,2])]	
730	plot(phy, cex=0.1, lwd=0.4, edge.width=0.5, edge.color=colors.branches, ...)
731}
732
733
734prior.posterior.plot=function(data, rpois.lambda) {
735# data should be a vector of number of shifts, for instance, for all interations across the MCMC sampling (post-burnin)
736	require(ggplot2)
737	dataset <- data.frame(variable = gl(2, length(data), labels = c("posterior", "prior")), value = c(data, rpois(length(data),rpois.lambda))) 
738	ggplot(data = dataset, aes(x = value, fill = variable)) + geom_histogram(position = "dodge") + scale_fill_manual(values = c("gray", "black"))	
739}
740
741
742