Loading [MathJax]/extensions/tex2jax.js
LALInference 4.1.9.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
LALInferenceNestedSampler.c
Go to the documentation of this file.
1/* Implementation of Nested Sampling for LALInference.
2 * (C) John Veitch, 2010
3 */
4
5#include <config.h>
6
7#include <stdio.h>
8#include <stdlib.h>
9#include <float.h>
10#include <signal.h>
11#include <time.h>
12#include <sys/times.h>
13#ifdef HAVE_UNISTD_H
14#include <unistd.h>
15#endif
16
17#include <lal/TimeDelay.h>
18#include <lal/LALInferenceConfig.h>
19#include <lal/LALStdlib.h>
20#include <lal/LALInferenceNestedSampler.h>
21#include <lal/LALInferencePrior.h>
22#include <lal/LALInferenceLikelihood.h>
23#include <lal/LALInferenceProposal.h>
24#include <lal/LALInferenceReadData.h>
25#include <lal/LALInferenceHDF5.h>
26#include <lal/LALInferencePriorVolumes.h>
27
28#include "logaddexp.h"
29
30#define PROGRAM_NAME "LALInferenceNestedSampler.c"
31#define CVS_ID_STRING "$Id$"
32#define CVS_REVISION "$Revision$"
33#define CVS_SOURCE "$Source$"
34#define CVS_DATE "$Date$"
35#define CVS_NAME_STRING "$Name$"
36
37#define MAX_MCMC 5000 /* Maximum chain length, set to be higher than expected from a reasonable run */
38#define ACF_TOLERANCE 0.01 /* Desired maximum correlation of MCMC samples */
39
40#ifdef __GNUC__
41#define UNUSED __attribute__ ((unused))
42#else
43#define UNUSED
44#endif
45
47
48/**
49 * structure holding internal state of the NS integrator
50 */
51typedef struct tagNSintegralState
52{
61
62static struct itimerval checkpoint_timer;
63
64/** Utility functions for the resume functionality */
65/** Write the current state to a checkpoint file the given filename */
66/** Read the given filename to populate a given LALInferenceRunState and NSintegralState */
69{
70 XLALH5FileAddScalarAttribute(group, "iteration", &(s->iteration), LAL_U4_TYPE_CODE );
71 XLALH5FileWriteREAL8Vector(group, "logZarray", s->logZarray);
72 XLALH5FileWriteREAL8Vector(group, "oldZarray", s->oldZarray);
73 XLALH5FileWriteREAL8Vector(group, "Harray", s->Harray);
74 XLALH5FileWriteREAL8Vector(group, "logwarray", s->logwarray);
75 XLALH5FileWriteREAL8Vector(group, "logtarray", s->logtarray);
76 XLALH5FileWriteREAL8Vector(group, "logt2array", s->logt2array);
77 return(0);
78}
79
82{
83 XLALH5FileQueryScalarAttributeValue(&(s->iteration), group, "iteration");
84 s->logZarray = XLALH5FileReadREAL8Vector(group, "logZarray");
85 s->oldZarray = XLALH5FileReadREAL8Vector(group, "oldZarray");
86 s->Harray = XLALH5FileReadREAL8Vector(group, "Harray");
87 s->logwarray = XLALH5FileReadREAL8Vector(group, "logwarray");
88 s->logtarray = XLALH5FileReadREAL8Vector(group, "logtarray");
89 s->logt2array = XLALH5FileReadREAL8Vector(group, "logt2array");
90 s->size = s->logZarray->length;
91 return(0);
92}
93
94
97
99{
100 LALH5File *h5file = XLALH5FileOpen(filename,"w");
101 int retcode;
102 if(!h5file)
103 {
104 fprintf(stderr,"Unable to save resume file %s!\n",filename);
105 return(1);
106 }
107 UINT4 Nlive=*(UINT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nlive");
109 XLAL_TRY(group = XLALH5GroupOpen(h5file,"lalinference"),retcode);
110 if(retcode!=XLAL_SUCCESS) return(retcode);
111 XLAL_TRY(group = XLALH5GroupOpen(group,"lalinferencenest_checkpoint"),retcode);
112 if(retcode!=XLAL_SUCCESS) return(retcode);
113 retcode = _saveNSintegralStateH5(group,s);
114 if(retcode) XLAL_ERROR(XLAL_EFAILED,"Unable to save integral state\n");
115 LALInferenceH5VariablesArrayToDataset(group, runState->livePoints, Nlive, "live_points");
116 INT4 N_output_array=0;
117 if(LALInferenceCheckVariable(runState->algorithmParams,"N_outputarray")) N_output_array=LALInferenceGetINT4Variable(runState->algorithmParams,"N_outputarray");
118 if(N_output_array>0)
119 {
120 LALInferenceVariables **output_array=NULL;
121 output_array=*(LALInferenceVariables ***)LALInferenceGetVariable(runState->algorithmParams,"outputarray");
122 LALInferenceH5VariablesArrayToDataset(group, output_array, N_output_array, "past_chain");
123 XLALH5FileAddScalarAttribute(group, "N_outputarray", &N_output_array, LAL_I4_TYPE_CODE);
124 }
125 struct tms tms_buffer;
126 if(times(&tms_buffer))
127 {
128 REAL8 execution_time = (REAL8)tms_buffer.tms_utime / (REAL8) sysconf (_SC_CLK_TCK);
129 REAL8 *prevtime=NULL;
130 if( (prevtime=LALInferenceGetVariable(runState->algorithmParams,"cpu_time"))) execution_time+=*prevtime;
131 XLALH5FileAddScalarAttribute(group, "cpu_time", &execution_time, LAL_D_TYPE_CODE);
132 }
134 XLALH5FileClose(h5file);
136 return(retcode);
137}
138
139static int CheckOutputFileContents(char *filename);
140/* Check what's in the output file (filename). Returns
141 * 0 : File does not exist or is not valid HDF5
142 * 1 : File contains a resume state
143 * 2 : File is valid HDF5 but not a resume file, probably complete
144 */
146{
147 int retcode = 0;
148 int result = 0;
149 LALH5File *h5file = NULL;
150 if (access(filename, F_OK)==-1) return(0);
152 XLAL_TRY(h5file = XLALH5FileOpen(filename,"r"),retcode);
153 if(retcode != XLAL_SUCCESS) return(0);
155 XLAL_TRY(group = XLALH5GroupOpen(h5file,"lalinference"),retcode);
156 if(retcode == XLAL_SUCCESS) result=2;
157 /* Look for checkpoint */
158 XLAL_TRY(XLALH5GroupOpen(group,"lalinferencenest_checkpoint"),retcode);
159 if(retcode == XLAL_SUCCESS) result=1;
160 XLALH5FileClose(h5file);
161 return(result);
162}
163
165{
166 int retcode;
167 LALH5File *h5file;
168 UINT4 Nlive;
169 if( access( filename, F_OK ) == -1 ) return(1);
171 XLAL_TRY(h5file = XLALH5FileOpen(filename,"r"),retcode);
172 if(retcode!=XLAL_SUCCESS) return(retcode);
174 XLAL_TRY(group = XLALH5GroupOpen(h5file,"lalinference"),retcode);
175 if(retcode!=XLAL_SUCCESS) return(retcode);
176 XLAL_TRY(group = XLALH5GroupOpen(group,"lalinferencenest_checkpoint"),retcode);
177 if(retcode!=XLAL_SUCCESS) return(retcode);
178 UINT4 N_outputarray;
179 LALInferenceVariables **outputarray;
180 XLALH5FileQueryScalarAttributeValue(&N_outputarray, group, "N_outputarray");
181 if(!h5file)
182 {
183 fprintf(stderr,"Unable to load resume file %s!\n",filename);
184 return(1);
185 }
186 printf("restoring nested sampling integral state\n");
188 if(retcode){
189 fprintf(stderr,"Unable to read nested sampling state - unable to resume!\n");
190 XLALH5FileClose(h5file);
191 return 1;
192 }
193 LALH5Dataset *liveGroup = XLALH5DatasetRead(group,"live_points");
194 retcode = LALInferenceH5DatasetToVariablesArray(liveGroup , &(runState->livePoints), &Nlive );
195 printf("restored %i live points\n",Nlive);
196 XLALH5DatasetFree(liveGroup);
197 if(N_outputarray>0)
198 {
199 printf("restoring %i past iterations\n",N_outputarray);
200 LALH5Dataset *outputGroup = XLALH5DatasetRead(group, "past_chain");
201 retcode |= LALInferenceH5DatasetToVariablesArray(outputGroup, &outputarray, &N_outputarray);
204 XLALH5DatasetFree(outputGroup);
205 }
206 double execution_time = 0.0;
207 if ( ( execution_time = XLALH5FileQueryScalarAttributeValue(&execution_time, group, "cpu_time") ) )
208 {
210 }
211
213 XLALH5FileClose(h5file);
214 printf("done restoring\n");
215 return(retcode);
216
217}
218
219/** Sync the live points to the differential evolution buffer */
221
222/* This is checked by the main loop to determine when to checkpoint */
223static volatile sig_atomic_t __ns_saveStateFlag = 0;
224/* This indicates the main loop should terminate */
225static volatile sig_atomic_t __ns_exitFlag = 0;
226
227/* Signal handler for SIGINT, which is produced by condor
228 * when evicting a job, or when the user presses Ctrl-C */
229static void catch_interrupt(UNUSED int sig, UNUSED siginfo_t *siginfo,UNUSED void *context);
230static void catch_interrupt(UNUSED int sig, UNUSED siginfo_t *siginfo,UNUSED void *context)
231{
232
235}
236
237/* Signal handler for SIGALRM, for periodic checkpointing */
238static void catch_alarm(UNUSED int sig, UNUSED siginfo_t *siginfo,UNUSED void *context);
239static void catch_alarm(UNUSED int sig, UNUSED siginfo_t *siginfo,UNUSED void *context)
240{
242}
243
244/**
245 * Install the signal handlers for checkpointing.
246 * If checkpoint_exit!=0, then install the catch_alarm_condor_exit_code
247 * handler to exit after checkpointing
248 */
249static void install_resume_handler(int checkpoint_exit);
250static void install_resume_handler(int checkpoint_exit)
251{
252 /* Install a periodic alarm that will trigger a checkpoint */
253 int sigretcode=0;
254 struct sigaction sa;
255 int checkpoint_time = 3600; /* Default 1 hour */
256 if (checkpoint_exit)
257 {
258 sa.sa_sigaction=catch_interrupt;
259 /* Increase time between checkpoints because there is overhead in restarting */
260 checkpoint_time = 3*3600; /* 3 hours */
261 }
262 else
263 {
264 sa.sa_sigaction=catch_alarm;
265 }
266 sa.sa_flags=SA_SIGINFO;
267 sigretcode=sigaction(SIGVTALRM,&sa,NULL);
268 if(sigretcode!=0) fprintf(stderr,"WARNING: Cannot establish checkpoint timer!\n");
269 /* Condor sends SIGUSR2 to checkpoint and continue */
270 sigretcode=sigaction(SIGUSR2,&sa,NULL);
271 if(sigretcode!=0) fprintf(stderr,"WARNING: Cannot establish checkpoint on SIGUSR2.\n");
272 checkpoint_timer.it_interval.tv_sec=checkpoint_time;
273 checkpoint_timer.it_interval.tv_usec=0;
274 checkpoint_timer.it_value=checkpoint_timer.it_interval;
275 setitimer(ITIMER_VIRTUAL,&checkpoint_timer,NULL);
276 /* Install the handler for the condor interrupt signal */
277 sa.sa_sigaction=catch_interrupt;
278 sigretcode=sigaction(SIGINT,&sa,NULL);
279 if(sigretcode!=0) fprintf(stderr,"WARNING: Cannot establish checkpoint on SIGINT.\n");
280 /* Condor sends SIGTERM to vanilla universe jobs to evict them */
281 sigretcode=sigaction(SIGTERM,&sa,NULL);
282 if(sigretcode!=0) fprintf(stderr,"WARNING: Cannot establish checkpoint on SIGTERM.\n");
283 /* Condor sends SIGTSTP to standard universe jobs to evict them.
284 *I think condor handles this, so didn't add a handler CHECK */
285}
286
287static UINT4 UpdateNMCMC(LALInferenceRunState *runState);
288/* Prototypes for private "helper" functions. */
289
290static REAL8 LALInferenceNSSample_logt(int Nlive,gsl_rng *RNG);
291
292//static void SamplePriorDiscardAcceptance(LALInferenceRunState *runState);
293static REAL8 mean(REAL8 *array,int N);
294
295
296/** Calculate covariance matrix from a collection of live points */
297static void LALInferenceNScalcCVM(gsl_matrix **cvm, LALInferenceVariables **Live, UINT4 Nlive);
298
299/* log( exp(a) - exp(b) ) */
300static double logsubexp(double a, double b);
301static double logsubexp(double a, double b)
302{
303 if(b>a)
304 {
305 fprintf(stderr,"Cannot take log of negative number %lf - %lf = %lf !\n",exp(a),exp(b),exp(a)-exp(b));
306 return -INFINITY;
307 }
308 else
309 {
310 return a + log1p(-exp(b-a));
311 }
312}
313
314static void SetupEigenProposals(LALInferenceRunState *runState);
315
316/**
317 * Update the internal state of the integrator after receiving the lowest logL
318 * value logL
319 */
320static REAL8 incrementEvidenceSamples(gsl_rng *GSLrandom, UINT4 Nlive, REAL8 logL, NSintegralState *s);
321static REAL8 incrementEvidenceSamples(gsl_rng *GSLrandom, UINT4 Nlive, REAL8 logL, NSintegralState *s)
322{
323 REAL8 Wtarray[s->size];
324 /* Update evidence array */
325 for(UINT4 j=0;j<s->size;j++){
326 s->oldZarray->data[j]=s->logZarray->data[j];
327 Wtarray[j]=s->logwarray->data[j]+logL+logsubexp(0,s->logt2array->data[j] + s->logtarray->data[j])-log(2.0);
328 s->logZarray->data[j]=logaddexp(s->logZarray->data[j],Wtarray[j]);
329 if(isfinite(s->oldZarray->data[j]) && isfinite(s->logZarray->data[j]))
330 s->Harray->data[j]= exp(Wtarray[j]-s->logZarray->data[j])*logL
331 + exp(s->oldZarray->data[j]-s->logZarray->data[j])*(s->Harray->data[j]+s->oldZarray->data[j])-s->logZarray->data[j];
332 REAL8 tmp=s->logtarray->data[j];
333 s->logtarray->data[j]=s->logt2array->data[j];
334 s->logt2array->data[j]=tmp;
335 s->logtarray->data[j]=LALInferenceNSSample_logt(Nlive,GSLrandom);
336 s->logwarray->data[j]+=s->logtarray->data[j];
337 }
338 s->iteration++;
339 return(mean(s->logZarray->data,s->logZarray->length));
340}
341
342static void printAdaptiveJumpSizes(FILE *file, LALInferenceThreadState *threadState);
344{
345 LALInferenceVariableItem *this=threadState->currentParams->head;
346 REAL8 *val=NULL;
347 char tmpname[1000]="";
348 int first=1;
349 while(this)
350 {
351 sprintf(tmpname,"%s_%s",this->name,ADAPTSUFFIX);
352 if(LALInferenceCheckVariable(threadState->proposalArgs,tmpname))
353 {
354 if(first) {fprintf(file,"Adaptive proposal step sizes:\n"); first=0;}
355 val=(REAL8 *)LALInferenceGetVariable(threadState->proposalArgs,tmpname);
356 fprintf(file,"%s: %lf\n",this->name,*val);
357 }
358 this=this->next;
359 }
360}
361
362
363/* Create Internal arrays for sampling the integral */
364static NSintegralState *initNSintegralState(UINT4 Nruns, UINT4 Nlive);
366{
368 s->iteration=0;
369 s->size=Nruns;
370 s->logZarray = XLALCreateREAL8Vector(Nruns);
371 s->oldZarray = XLALCreateREAL8Vector(Nruns);
372 s->Harray = XLALCreateREAL8Vector(Nruns);
373 s->logwarray = XLALCreateREAL8Vector(Nruns);
374 s->logtarray=XLALCreateREAL8Vector(Nruns);
375 s->logt2array=XLALCreateREAL8Vector(Nruns);
376 REAL8 logw=0;
377
378 if(s->logZarray==NULL || s->Harray==NULL || s->oldZarray==NULL || s->logwarray==NULL)
379 {fprintf(stderr,"Unable to allocate RAM\n"); exit(-1);}
380 for(UINT4 i=0;i<Nruns;i++) {s->logwarray->data[i]=logw; s->logZarray->data[i]=-INFINITY;
381 s->oldZarray->data[i]=-INFINITY; s->Harray->data[i]=0.0;s->logtarray->data[i]=-1.0/Nlive;
382 s->logt2array->data[i]=-1.0/Nlive;
383 }
384 return s;
385}
386
387static REAL8 mean(REAL8 *array,int N){
388 REAL8 sum=0.0;
389 int i;
390 for(i=0;i<N;i++) sum+=array[i];
391 return sum/((REAL8) N);
392}
393
394
395static REAL8 LALInferenceNSSample_logt(int Nlive,gsl_rng *RNG){
396 REAL8 t=0.0;
397 REAL8 a=0.0;
398 while((Nlive--)>1) {a=gsl_rng_uniform(RNG); t = t>a ? t : a;}
399 return(log(t));
400}
401
403 INT4 max = 0;
404 INT4 maxMCMC = MAX_MCMC;
405 /* Measure Autocorrelations if the Nmcmc is not over-ridden */
406 if(!LALInferenceGetProcParamVal(runState->commandLine,"--Nmcmc") && !LALInferenceGetProcParamVal(runState->commandLine,"--nmcmc")){
407 if(LALInferenceCheckVariable(runState->algorithmParams,"maxmcmc"))
408 maxMCMC = *(INT4 *)LALInferenceGetVariable(runState->algorithmParams,"maxmcmc");
409 if(LALInferenceCheckVariable(runState->algorithmParams,"Nmcmc")){ /* if already estimated the length */
410 if ( LALInferenceGetINT4Variable(runState->algorithmParams,"Nmcmc") != 0 ){
411 max=4 * *(INT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nmcmc"); /* We will use this to go out 4x last ACL */
412 }
413 else { max=4*maxMCMC; } /* otherwise use the MAX_MCMC */
414 }
415 else max=4*maxMCMC; /* otherwise use the MAX_MCMC */
416 if(max>4*maxMCMC) max=4*maxMCMC;
417 UINT4 Nlive = *(INT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nlive");
418 UINT4 rnd=gsl_rng_uniform_int(runState->GSLrandom,Nlive);
419 /* Single threaded here */
420 LALInferenceCopyVariables(runState->livePoints[rnd],runState->threads[0].currentParams);
422 max=10;
423 for(LALInferenceVariableItem *this=acls->head;this;this=this->next) {
424 if(LALInferenceCheckVariable(runState->algorithmParams,"verbose"))
425 fprintf(stdout,"Autocorrelation length of %s: %i\n",this->name,(INT4) *(REAL8 *)this->value);
426 if(*(REAL8 *)this->value>max) {
427 max=(INT4) *(REAL8 *)this->value;
428 }
429 }
431 XLALFree(acls);
432 if(max>maxMCMC){
433 fprintf(stderr,"Warning: Estimated chain length %i exceeds maximum %i!\n",max,maxMCMC);
434 max=maxMCMC;
435 }
436 LALInferenceSetVariable(runState->algorithmParams,"Nmcmc",&max);
437 }
438 /* Single threaded here */
439 if (LALInferenceGetProcParamVal(runState->commandLine,"--proposal-kde"))
441 return(max);
442}
443
444/* estimateCovarianceMatrix reads the list of live points,
445 and works out the covariance matrix of the varying parameters
446 - CIRCULAR parameters are wrapped around before the calculation and must be
447 scaled to the range 0 -> 2pi. Only works with REAL8 variables */
448static void LALInferenceNScalcCVM(gsl_matrix **cvm, LALInferenceVariables **Live, UINT4 Nlive)
449{
450 UINT4 i,j,k;
451 UINT4 ND=0;
452 LALInferenceVariableItem *item,*k_item,*j_item;
453 REAL8 *means, *ms, *mc,jval=0.,kval=0.;
454
455 /* Find the number of dimensions which vary in the covariance matrix */
456 for(item=Live[0]->head;item!=NULL;item=item->next)
458
459 /* Set up matrix if necessary */
460 if(*cvm==NULL)
461 {if(NULL==(*cvm=gsl_matrix_alloc(ND,ND))) {fprintf(stderr,"Unable to allocate matrix memory\n"); exit(1);}}
462 else {
463 if((*cvm)->size1!=(*cvm)->size2 || (*cvm)->size1!=ND)
464 { fprintf(stderr,"ERROR: Matrix wrong size. Something has gone wrong in LALInferenceNScalcCVM\n");
465 exit(1);
466 }
467 }
468 /* clear the matrix */
469 for(i=0;i<(*cvm)->size1;i++) for(j=0;j<(*cvm)->size2;j++) gsl_matrix_set(*cvm,i,j,0.0);
470
471 /* Find the means */
472 if(NULL==(means = XLALMalloc((size_t)ND*sizeof(REAL8)))){fprintf(stderr,"Can't allocate RAM"); exit(-1);}
473 if(NULL==(ms = XLALMalloc((size_t)ND*sizeof(REAL8)))){fprintf(stderr,"Can't allocate RAM"); exit(-1);}
474 if(NULL==(mc = XLALMalloc((size_t)ND*sizeof(REAL8)))){fprintf(stderr,"Can't allocate RAM"); exit(-1);}
475 for(i=0;i<ND;i++){
476 means[i]=0.0;
477 ms[i] = 0.;
478 mc[i] = 0.;
479 }
480 for(item=Live[0]->head,j=0;item;item=item->next){
482 for(i=0;i<Nlive;i++){
483 void *ptr = LALInferenceGetVariable(Live[i],item->name);
484 switch(item->vary)
485 {
487 {
488 means[j]+=*(REAL8 *)ptr;
489 break;
490 }
492 {
493 ms[j] += sin(*(REAL8 *)ptr);
494 mc[j] += cos(*(REAL8 *)ptr);
495 break;
496 }
497 default:
498 {
499 break;
500 }
501 }
502 }
503 j++;
504 }
505
506 for(item=Live[0]->head,j=0;item;item=item->next){
508 means[j]/=(REAL8)Nlive;
509 j++;
510 }
512 ms[j]/=(REAL8)Nlive;
513 mc[j]/=(REAL8)Nlive;
514
515 means[j] = atan2(ms[j], mc[j]);
516 means[j] = means[j]<0? 2.0*LAL_PI + means[j] : means[j];
517
518 j++;
519 }
520 }
521
522 XLALFree(ms);
523 XLALFree(mc);
524
525 /* Iterate over the matrix elements one at a time */
526 j_item=k_item=Live[0]->head;
527 for(j=0,j_item=Live[0]->head;j_item;j_item=j_item->next)
528 {
529 /* Skip non-varying parameters */
530 if((j_item->vary!=LALINFERENCE_PARAM_LINEAR && j_item->vary!=LALINFERENCE_PARAM_CIRCULAR )||j_item->type!=LALINFERENCE_REAL8_t ) {continue;}
531
532 for(k_item=Live[0]->head,k=0;k_item &&k<=j;k_item=k_item->next)
533 {
534 /* Skip non-varying parameters */
535 if((k_item->vary!=LALINFERENCE_PARAM_LINEAR && k_item->vary!=LALINFERENCE_PARAM_CIRCULAR)||k_item->type!=LALINFERENCE_REAL8_t) {continue;}
536
537 /* Loop over live points updating covariance elements */
538 for(i=0;i<Nlive;i++)
539 {
540 void *jptr=LALInferenceGetVariable(Live[i],j_item->name);
541 void *kptr=LALInferenceGetVariable(Live[i],k_item->name);
542
543 /* Check for linear or circular */
544 if( j_item->vary==LALINFERENCE_PARAM_CIRCULAR )
545 jval = LALInferenceAngularDistance(*(REAL8 *)jptr, means[j]);
546 else if( j_item->vary==LALINFERENCE_PARAM_LINEAR )
547 jval = *(REAL8 *)jptr - means[j];
548
549 if( k_item->vary==LALINFERENCE_PARAM_CIRCULAR )
550 kval = LALInferenceAngularDistance(*(REAL8 *)kptr, means[k]);
551 else if( k_item->vary==LALINFERENCE_PARAM_LINEAR )
552 kval = *(REAL8 *)kptr - means[k];
553
554 gsl_matrix_set(*cvm,j,k, gsl_matrix_get(*cvm,j,k) + kval*jval);
555
556 }
557
558 k++;
559 } /* end loop over k */
560 j++;
561 }/* end loop over j */
562
563 /* Normalise */
564 for(i=0;i<ND;i++) for(j=0;j<ND;j++) gsl_matrix_set(*cvm,i,j,gsl_matrix_get(*cvm,i,j)/((REAL8) Nlive));
565 XLALFree(means);
566
567 /* the other half */
568 for(i=0;i<ND;i++)
569 for(j=0;j<i;j++)
570 gsl_matrix_set(*cvm,j,i,gsl_matrix_get(*cvm,i,j));
571
572 return;
573}
574
576{
577 char help[]="\
578 ----------------------------------------------\n\
579 --- Nested Sampling Algorithm Parameters -----\n\
580 ----------------------------------------------\n\
581 --Nlive N Number of live points to use\n\
582 (--Nmcmc M) Over-ride auto chain length determination and use <M> MCMC samples\n\
583 (--maxmcmc M) Use at most M MCMC points when autodetermining the chain (5000)\n\
584 (--Nmcmcinitial M) Use M MCMC points when initially resampling from the prior\n\
585 (otherwise default is to use maxmcmc)\n\
586 (--sloppyratio S) Number of sub-samples of the prior for every sample from the\n\
587 limited prior\n\
588 (--Nruns R) Number of parallel samples from logt to use(1)\n\
589 (--tolerance dZ) Tolerance of nested sampling algorithm (0.1)\n\
590 (--randomseed seed) Random seed of sampling distribution\n\
591 (--prior ) Set the prior to use (InspiralNormalised,SkyLoc,malmquist)\n\
592 (default: InspiralNormalised)\n\
593 (--sampleprior N) For Testing: Draw N samples from the prior, will not perform the\n\
594 nested sampling integral\n\
595 (--progress) Output some progress information at each iteration\n\
596 (--verbose) Output more info. N=1: errors, N=2 (default): warnings, N=3: info\n\
597 (--resume) Allow non-condor checkpointing every 4 hours. If given will check \n\
598 or OUTFILE_resume and continue if possible\n\
599 (--checkpoint-exit-code N) Exit with code N when checkpoint is complete.\n\
600 For use with condor's +SuccessCheckpointExitCode option\n\
601 \n";
602
604 /* Print command line arguments if help requested */
605 if(runState == NULL || LALInferenceGetProcParamVal(runState->commandLine,"--help"))
606 {
607 fprintf(stdout,"%s",help);
608 return;
609 }
610 ProcessParamsTable *commandLine=runState->commandLine;
611
612 INT4 verbose=0;
613 INT4 x=0;
614 ppt=LALInferenceGetProcParamVal(commandLine,"--verbose");
615 if(ppt) {
616 if(ppt->value[0]){
617 x=atoi(ppt->value);
618 switch(x){
619 case 0:
620 verbose=LALNDEBUG; /* Nothing */
621 break;
622 case 1:
623 verbose=LALMSGLVL1; /* Only errors */
624 break;
625 case 2:
626 verbose=LALMSGLVL2; /* Errors and warnings */
627 break;
628 case 3:
629 verbose=LALMSGLVL3; /* Errors, warnings and info */
630 break;
631 default:
633 break;
634 }
635 }
636 else verbose=LALMSGLVL2; /* Errors and warnings */
639 }
640 INT4 tmpi=0;
641 REAL8 tmp=0;
642
643 /* Single thread only */
644 LALInferenceThreadState *threadState = &runState->threads[0];
645
646 /* Set up the appropriate functions for the nested sampling algorithm */
649
650 /* use the ptmcmc proposal to sample prior */
652 REAL8 temp=1.0;
654
656
657 /* Number of live points */
658 ppt=LALInferenceGetProcParamVal(commandLine,"--Nlive");
659 if (!ppt) ppt=LALInferenceGetProcParamVal(commandLine,"--nlive");
660 if(ppt)
661 tmpi=atoi(ppt->value);
662 else {
663 fprintf(stderr,"Error, must specify number of live points\n");
664 exit(1);
665 }
667
668 /* Number of points in MCMC chain */
669 ppt=LALInferenceGetProcParamVal(commandLine,"--Nmcmc");
670 if(!ppt) ppt=LALInferenceGetProcParamVal(commandLine,"--nmcmc");
671 if(ppt){
672 tmpi=atoi(ppt->value);
673 LALInferenceAddVariable(runState->algorithmParams,"Nmcmc",&tmpi,
675 printf("set number of MCMC points, over-riding auto-determination!\n");
676 }
677
678 /* Maximum number of points in MCMC chain */
679 ppt=LALInferenceGetProcParamVal(commandLine,"--maxmcmc");
680 if(ppt){
681 tmpi=atoi(ppt->value);
682 LALInferenceAddVariable(runState->algorithmParams,"maxmcmc",&tmpi,
684 }
685
686 /* Set fraction for sloppy sampling */
687 if((ppt=LALInferenceGetProcParamVal(commandLine,"--sloppyfraction")))
688 tmp=atof(ppt->value);
689 else tmp=0.0;
690 LALInferenceAddVariable(runState->algorithmParams,"sloppyfraction",&tmp,
692
693 /* Optionally specify number of parallel runs */
694 ppt=LALInferenceGetProcParamVal(commandLine,"--Nruns");
695 if(ppt) {
696 tmpi=atoi(ppt->value);
698 }
699
700 printf("set tolerance.\n");
701 /* Tolerance of the Nested sampling integrator */
702 ppt=LALInferenceGetProcParamVal(commandLine,"--tolerance");
703 if(ppt){
704 tmp=strtod(ppt->value,(char **)NULL);
707 }
708 double zero=0.0;
710
711 return;
712
713}
714
715
716/* NestedSamplingAlgorithm implements the nested sampling algorithm,
717 see e.g. Sivia & Skilling "Data Analysis: A Bayesian Tutorial, 2nd edition.
718 REQUIREMENTS:
719 Calling routine must have set up runState->livePoints already to
720 contain samples from the prior distribution.
721 runState->algorithmParams must contain a variable "logLikelihoods"
722 which contains a REAL8 array of likelihood values for the live
723 points.
724 */
725
727{
728 UINT4 iter=0,i,j,minpos;
729 /* Single thread here */
730 LALInferenceThreadState *threadState = &runState->threads[0];
731 UINT4 HDFOUTPUT=1;
732 UINT4 Nlive=*(UINT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nlive");
733 UINT4 Nruns=100;
734 REAL8 *logZarray,*Harray,*logwarray,*logtarray;
735 REAL8 TOLERANCE=0.1;
736 REAL8 logZ,logZnew,logLmin,logLmax=-INFINITY,logLtmp,logw,H,logZnoise,dZ=0;
738 FILE *fpout=NULL;
739 REAL8 neginfty=-INFINITY;
740 REAL8 zero=0.0;
741 REAL8 *logLikelihoods=NULL;
742 UINT4 verbose=0;
743 REAL8 sloppyfrac;
744 UINT4 displayprogress=0;
746 UINT4 samplePrior=0; //If this flag is set to a positive integer, code will just draw this many samples from the prior
748 int CondorExitCode=0;
749
750 if(!runState->logsample) runState->logsample=LALInferenceLogSampleToArray;
751
752 if((ppt=LALInferenceGetProcParamVal(runState->commandLine,"--checkpoint-exit-code")))
753 CondorExitCode=atoi(ppt->value);
754
755 if ( !LALInferenceCheckVariable(runState->algorithmParams, "logZnoise" ) ){
756 logZnoise=LALInferenceNullLogLikelihood(runState->data);
757
759 }
760 else{
761 logZnoise =
762 *(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"logZnoise");
763 }
764
765 logLikelihoods=(REAL8 *)(*(REAL8Vector **)LALInferenceGetVariable(runState->algorithmParams,"logLikelihoods"))->data;
766
767 if((ppt=LALInferenceGetProcParamVal(runState->commandLine,"--sampleprior")))
768 {
769 samplePrior=atoi(ppt->value);
770 fprintf(stdout,"Generating %i samples from the prior\n",samplePrior);
771 }
772
774 displayprogress=verbose;
775
776 /* Operate on parallel runs if requested */
777 if(LALInferenceCheckVariable(runState->algorithmParams,"Nruns"))
778 Nruns = *(UINT4 *) LALInferenceGetVariable(runState->algorithmParams,"Nruns");
779
780 /* Create workspace for arrays */
781 NSintegralState *s=NULL;
782
783 if(LALInferenceCheckVariable(runState->algorithmParams,"tolerance"))
784 TOLERANCE = *(REAL8 *) LALInferenceGetVariable(runState->algorithmParams,"tolerance");
785
786 /* Check that necessary parameters are created */
787 if(!LALInferenceCheckVariable(runState->algorithmParams,"logLmin"))
789
790 if(!LALInferenceCheckVariable(runState->algorithmParams,"accept_rate"))
792 if(!LALInferenceCheckVariable(runState->algorithmParams,"sub_accept_rate"))
794 if(!LALInferenceCheckVariable(runState->algorithmParams,"sloppyfraction"))
796
797 /* Open output file */
798 ppt = NULL;
799 ppt=LALInferenceGetProcParamVal(runState->commandLine,"--outfile");
800 if(!ppt){
801 fprintf(stderr,"Must specify --outfile <filename.hdf5>\n");
802 exit(1);
803 }
804 char *outfile=ppt->value;
805 /* Check if the output file has hdf5 extension */
806 if(strstr(outfile,".h5") || strstr(outfile,".hdf")) HDFOUTPUT=1;
807 else HDFOUTPUT=0;
808
809 double logvolume=0.0;
810 if ( LALInferenceCheckVariable( runState->livePoints[0], "chirpmass" ) ){
811 /* If a cbc run, calculate the mass-distance volume and store it to file*/
812 /* Do it before algorithm starts so that we can kill the run and still get this */
813 int errnum=0;
814 XLAL_TRY(logvolume=log(LALInferenceMassDistancePriorVolume(runState)), errnum);
815 if(XLAL_IS_REAL8_FAIL_NAN(logvolume))
816 {
817 XLALPrintWarning("Could not compute log volume: %s\n",XLALErrorString(errnum));
818 logvolume=0.0;
819 }
820 }
821
822 if(LALInferenceGetProcParamVal(runState->commandLine,"--progress"))
823 displayprogress=1;
824
825 minpos=0;
826
827 logw=log(1.0-exp(-1.0/Nlive));
828 i=0;
829 /* sort points for consistent order before creating the matrix */
830 for(i=0;i<Nlive;i++)
832
833 /* Set up eigenvector proposals */
834 SetupEigenProposals(runState);
835
836 /* Use the live points as differential evolution points */
837 /* Single thread here */
838 syncLivePointsDifferentialPoints(runState,threadState);
839 threadState->differentialPointsSkip=1;
840
841 if(!LALInferenceCheckVariable(runState->algorithmParams,"Nmcmc")){
842 INT4 tmp=MAX_MCMC;
843 if(LALInferenceGetProcParamVal(runState->commandLine,"--Nmcmcinitial")){
844 tmp=atoi(LALInferenceGetProcParamVal(runState->commandLine,"--Nmcmcinitial")->value);
845 }
847 }
848 s=initNSintegralState(Nruns,Nlive);
849
850 /* Check if output/resume file exists as a valid HDF5 file */
851 int filetest = CheckOutputFileContents(outfile);
852 int retcode=1;
853 if (filetest==2) /* Run is complete, do not overwrite */
854 {
855 printf("Output file %s contains complete run, not over-writing\n",outfile);
856 exit(1);
857 }
858 else if(filetest==1) /* Run contains a resume file */
859 {
860 /* Check for an interrupted run */
861 if(LALInferenceGetProcParamVal(runState->commandLine,"--resume")){
862 fprintf(stderr,"Resuming from %s\n",outfile);
863 retcode=ReadNSCheckPointH5(outfile,runState,s);
864 if(retcode==0){
865 for(i=0;i<Nlive;i++) logLikelihoods[i]=*(REAL8 *)LALInferenceGetVariable(runState->livePoints[i],"logL");
866 iter=s->iteration;
867 }
868 }
869 }
870 /* If we did not resume, start a new run */
871 if(retcode!=0)
872 {
873
874 if(LALInferenceGetProcParamVal(runState->commandLine,"--resume"))
875 fprintf(stdout,"Starting a new run.\n");
876 /* Sprinkle points */
877 LALInferenceSetVariable(runState->algorithmParams,"logLmin",&neginfty);
878 int sprinklewarning=0;
879 for(i=0;i<Nlive;i++) {
880 threadState->currentParams=runState->livePoints[i];
881 j=0;
882 do
883 {
884 runState->evolve(runState);
885 logLikelihoods[i]=runState->likelihood(runState->livePoints[i],runState->data,threadState->model);
886 if(!sprinklewarning && j++==100) { fprintf(stderr,"Warning: having difficulty sampling prior, check your prior bounds\n"); sprinklewarning=1;}
887 }while(isnan(logLikelihoods[i]) || isinf(logLikelihoods[i]));
888 if(XLALPrintProgressBar((double)i/(double)Nlive)) fprintf(stderr,"\n");
889 }
890 }
891
892 logZarray=s->logZarray->data;
893 logtarray=s->logtarray->data;
894 logwarray=s->logwarray->data;
895 Harray=s->Harray->data;
896
897 /* Find maximum likelihood and sanity check */
898 for(i=0;i<Nlive;i++)
899 {
901 logLtmp=logLikelihoods[i];
902 logLmax=logLtmp>logLmax? logLtmp : logLmax;
903 if(isnan(logLikelihoods[i]) || isinf(logLikelihoods[i])) {
904 fprintf(stderr,"Detected logL[%i]=%lf! Sanity checking...\n",i,logLikelihoods[i]);
905 if(LALInferenceSanityCheck(runState))
906 exit(1);
907 }
908 }
909 /* sort points for consistent order after attaching delta parameters */
910 for(i=0;i<Nlive;i++)
912
913 /* Update the covariance matrix for proposal distribution */
914 SetupEigenProposals(runState);
915
916 /* Reset proposal stats before starting */
918
919 /* Set the number of MCMC points */
920 UpdateNMCMC(runState);
921 /* Output some information */
922 if(verbose){
923 LALInferencePrintProposalStatsHeader(stdout,threadState->cycle);
924 LALInferencePrintProposalStats(stdout,threadState->cycle);
926 printAdaptiveJumpSizes(stdout, threadState);
927 }
928 /* Write out names of parameters */
929 if(!HDFOUTPUT)
930 {
931 FILE *lout=NULL;
932 char param_list[FILENAME_MAX+16];
933 snprintf(param_list,sizeof(param_list),"%s_params.txt",outfile);
934 lout=fopen(param_list,"w");
936 fclose(lout);
937 }
938 minpos=0;
939 threadState->currentParams=currentVars;
940 fprintf(stdout,"Starting nested sampling loop!\n");
941 /* Install interrupt handler for resuming */
942 if(LALInferenceGetProcParamVal(runState->commandLine,"--resume"))
943 {
944 install_resume_handler(CondorExitCode);
945 }
946 /* Iterate until termination condition is met */
947 do {
948 /* Find minimum likelihood sample to replace */
949 minpos=0;
950 for(i=1;i<Nlive;i++){
951 if(logLikelihoods[i]<logLikelihoods[minpos])
952 minpos=i;
953 }
954 logLmin=logLikelihoods[minpos];
955 if(samplePrior) logLmin=-INFINITY;
956
957 logZnew=incrementEvidenceSamples(runState->GSLrandom, Nlive, logLikelihoods[minpos], s);
958 //deltaZ=logZnew-logZ; - set but not used
959 H=mean(Harray,Nruns);
960 logZ=logZnew;
961 if(runState->logsample) runState->logsample(runState->algorithmParams,runState->livePoints[minpos]);
962 UINT4 itercounter=0;
963
964 /* Generate a new live point */
965 do{ /* This loop is here in case it is necessary to find a different sample */
966 /* Clone an old live point and evolve it */
967 while((j=gsl_rng_uniform_int(runState->GSLrandom,Nlive))==minpos){};
968 LALInferenceCopyVariables(runState->livePoints[j],threadState->currentParams);
969 threadState->currentLikelihood = logLikelihoods[j];
970 LALInferenceSetVariable(runState->algorithmParams,"logLmin",(void *)&logLmin);
971 runState->evolve(runState);
972 itercounter++;
973 }while( threadState->currentLikelihood<=logLmin || *(REAL8*)LALInferenceGetVariable(runState->algorithmParams,"accept_rate")==0.0);
974
975 LALInferenceCopyVariables(threadState->currentParams,runState->livePoints[minpos]);
976 logLikelihoods[minpos]=threadState->currentLikelihood;
977
978 if (threadState->currentLikelihood>logLmax)
979 logLmax=threadState->currentLikelihood;
980
981 logw=mean(logwarray,Nruns);
983 dZ=logaddexp(logZ,logLmax-((double) iter)/((double)Nlive))-logZ;
984 sloppyfrac=*(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"sloppyfraction");
985 if(displayprogress) fprintf(stderr,"%i: accpt: %1.3f Nmcmc: %i sub_accpt: %1.3f slpy: %2.1f%% H: %3.2lf nats logL:%.3lf ->%.3lf logZ: %.3lf deltalogLmax: %.2lf dZ: %.3lf Zratio: %.3lf \n",\
986 iter,\
987 *(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"accept_rate")/(REAL8)itercounter,\
988 *(INT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nmcmc"),\
989 *(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"sub_accept_rate"),\
990 100.0*sloppyfrac,\
991 H,\
992 logLmin,\
993 threadState->currentLikelihood,\
994 logZ,\
995 (logLmax - LALInferenceGetREAL8Variable(runState->algorithmParams,"logZnoise")), \
996 dZ,\
997 ( logZ - LALInferenceGetREAL8Variable(runState->algorithmParams,"logZnoise"))\
998 );
999 iter++;
1000
1001 /* Save progress */
1002 if(__ns_saveStateFlag!=0)
1003 {
1004 if(__ns_exitFlag) fprintf(stdout,"Saving state to %s.\n",outfile);
1005 WriteNSCheckPointH5(outfile,runState,s);
1006 fflush(fpout);
1008 }
1009 /* Have we been told to quit? */
1010 if(__ns_exitFlag) {
1011 exit(CondorExitCode);
1012 }
1013
1014 /* Update the proposal */
1015 if(!(iter%(Nlive/10))) {
1016 /* Update the covariance matrix */
1017 if ( LALInferenceCheckVariable( threadState->proposalArgs,"covarianceMatrix" ) ){
1018 SetupEigenProposals(runState);
1019 }
1020
1021 /* Update NMCMC from ACF */
1022 UpdateNMCMC(runState);
1023
1024 /* Sync the live points to differential points */
1025 syncLivePointsDifferentialPoints(runState,threadState);
1026
1027 /* Output some information */
1028 if(verbose){
1029 LALInferencePrintProposalStatsHeader(stdout,threadState->cycle);
1030 LALInferencePrintProposalStats(stdout,threadState->cycle);
1032 printAdaptiveJumpSizes(stdout, threadState);
1033 }
1034 }
1035
1036 }
1037 while(samplePrior?((Nlive+iter)<samplePrior):( iter <= Nlive || dZ> TOLERANCE)); /* End of NS loop! */
1038
1039 /* Sort the remaining points (not essential, just nice)*/
1040 for(i=0;i<Nlive-1;i++){
1041 minpos=i;
1042 logLmin=logLikelihoods[i];
1043 for(j=i+1;j<Nlive;j++){
1044 if(logLikelihoods[j]<logLmin)
1045 {
1046 minpos=j;
1047 logLmin=logLikelihoods[j];
1048 }
1049 }
1050 temp=runState->livePoints[minpos]; /* Put the minimum remaining point in the current position */
1051 runState->livePoints[minpos]=runState->livePoints[i];
1052 runState->livePoints[i]=temp;
1053 logLikelihoods[minpos]=logLikelihoods[i];
1054 logLikelihoods[i]=logLmin;
1055 }
1056 /* final corrections */
1057 for(i=0;i<Nlive;i++){
1058 logZ=incrementEvidenceSamples(runState->GSLrandom, Nlive-i, logLikelihoods[i], s);
1059 if(runState->logsample) runState->logsample(runState->algorithmParams,runState->livePoints[i]);
1060 }
1061
1062 LALInferenceVariables **output_array=NULL;
1063 UINT4 N_output_array=0;
1064 if(LALInferenceCheckVariable(runState->algorithmParams,"outputarray")
1065 &&LALInferenceCheckVariable(runState->algorithmParams,"N_outputarray") )
1066 {
1067 output_array=*(LALInferenceVariables ***)LALInferenceGetVariable(runState->algorithmParams,"outputarray");
1068 N_output_array=*(UINT4 *)LALInferenceGetVariable(runState->algorithmParams,"N_outputarray");
1069 }
1070 double logB=logZ-logZnoise;
1071 /* Pass output back through algorithmparams */
1075
1076 /* Write out the evidence */
1077 if(!HDFOUTPUT)
1078 {
1079 fpout = fopen(outfile,"w");
1080 for(i=0;i<N_output_array;i++)
1081 {
1082 LALInferencePrintSample(fpout,output_array[i]);
1083 fprintf(fpout,"\n");
1084 }
1085 fclose(fpout);
1086
1087 char bayesfile[FILENAME_MAX+10];
1088 sprintf(bayesfile,"%s_B.txt",outfile);
1089 fpout=fopen(bayesfile,"w");
1090 fprintf(fpout,"%lf %lf %lf %lf\n",logZ-logZnoise,logZ,logZnoise,logLmax);
1091 fclose(fpout);
1092
1093
1094
1095 }
1096 /* Write HDF5 file */
1097 if(HDFOUTPUT)
1098 {
1099
1100 LALH5File *h5file=XLALH5FileOpen(outfile, "w");
1101 // Create group heirarchy
1102 char runID[2048];
1103 if((ppt=LALInferenceGetProcParamVal(runState->commandLine,"--runid")))
1104 snprintf(runID,sizeof(runID),"%s_%s","lalinference_nest",ppt->value);
1105 else
1106 snprintf(runID,sizeof(runID),"lalinference_nest");
1107
1108 LALH5File *groupPtr = LALInferenceH5CreateGroupStructure(h5file, "lalinference", runID);
1109 /* Create run identifier group */
1111 /* TODO: Write metadata */
1112 XLALH5FileAddScalarAttribute(groupPtr, "log_evidence", &logZ, LAL_D_TYPE_CODE);
1113 XLALH5FileAddScalarAttribute(groupPtr, "log_bayes_factor", &logB, LAL_D_TYPE_CODE);
1114 XLALH5FileAddScalarAttribute(groupPtr, "information_nats", &H, LAL_D_TYPE_CODE);
1115 XLALH5FileAddScalarAttribute(groupPtr, "log_noise_evidence", &logZnoise, LAL_D_TYPE_CODE );
1116 XLALH5FileAddScalarAttribute(groupPtr, "log_max_likelihood", &logLmax , LAL_D_TYPE_CODE);
1117 XLALH5FileAddScalarAttribute(groupPtr, "number_live_points", &Nlive, LAL_U4_TYPE_CODE);
1118 XLALH5FileAddScalarAttribute(groupPtr, "log_prior_volume", &logvolume, LAL_D_TYPE_CODE);
1119 /*Add the whole command line to the output headers */
1120 char *cl=NULL;
1122 XLALH5FileAddStringAttribute(groupPtr,"CommandLine",cl);
1123
1124 /* Get the cpu usage */
1125 struct tms tms_buffer;
1126 if(times(&tms_buffer))
1127 {
1128 REAL8 execution_time = (REAL8)tms_buffer.tms_utime / (REAL8) sysconf (_SC_CLK_TCK) + LALInferenceGetREAL8Variable(runState->algorithmParams,"cpu_time");
1129 XLALH5FileAddScalarAttribute(groupPtr, "cpu_time", &execution_time, LAL_D_TYPE_CODE);
1130 }
1131
1132 LALInferenceVariables *injParams = NULL;
1133 if ( (injParams=LALInferencePrintInjectionSample(runState)) )
1134 {
1135 LALInferenceH5VariablesArrayToDataset(groupPtr, &injParams, 1, "injection_params");
1136 LALInferenceClearVariables(injParams);
1137 XLALFree(injParams);
1138 }
1139 XLALH5FileClose(h5file);
1141 }
1142
1143 if(output_array) {
1144 for(i=0;i<N_output_array;i++){
1145 LALInferenceClearVariables(output_array[i]);
1146 XLALFree(output_array[i]);
1147 }
1148 XLALFree(output_array);
1149 }
1150
1151 /* Free memory */
1152 XLALFree(logtarray); XLALFree(logwarray); XLALFree(logZarray);
1153}
1154
1155/* Calculate the autocorrelation function of the sampler (runState->evolve) for each parameter
1156 * Evolves the sample starting with the value passed in temp, with a maximum of max_iterations steps.
1157 Return the ACL for each parameter as a LALInferenceVariables */
1159{
1160 /* Single threaded here */
1161 LALInferenceThreadState *threadState = &runState->threads[0];
1162 ProcessParamsTable *ppt=NULL;
1163 char chainfilename[2048]="";
1164 char acf_file_name[2048]="";
1165 FILE *chainfile=NULL;
1166 FILE *acffile=NULL;
1167 UINT4 i,j;
1168 UINT4 nPar=0; // = LALInferenceGetVariableDimensionNonFixed(runState->currentParams);
1169 REAL8 **data_array=NULL;
1170 REAL8 **acf_array=NULL;
1172 INT4 thinning=1;
1173 max_iterations/=thinning;
1174 /* Find the number and names of variables */
1175 for(this=threadState->currentParams->head;this;this=this->next) if(this->vary!=LALINFERENCE_PARAM_FIXED && this->vary!=LALINFERENCE_PARAM_OUTPUT && this->type==LALINFERENCE_REAL8_t) nPar++;
1176 char **param_names=XLALCalloc(nPar,sizeof(char *));
1177 for(i=0,this=threadState->currentParams->head;this;this=this->next) if(this->vary!=LALINFERENCE_PARAM_FIXED && this->vary!=LALINFERENCE_PARAM_OUTPUT && this->type==LALINFERENCE_REAL8_t) param_names[i++]=this->name;
1178
1179 REAL8 ACF,ACL,max=0;
1181
1182 /* Back up the algorithm state and replace with a clean version for logSampletoarray */
1183 LALInferenceVariables myAlgParams,*oldAlgParams=runState->algorithmParams;
1184 LALInferenceVariables myCurrentParams,*oldCurrentParams=threadState->currentParams;
1185 memset(&myAlgParams,0,sizeof(LALInferenceVariables));
1186 memset(&myCurrentParams,0,sizeof(LALInferenceVariables));
1187 LALInferenceCopyVariables(oldAlgParams,&myAlgParams);
1188 if(LALInferenceCheckVariable(&myAlgParams,"outputarray")) LALInferenceRemoveVariable(&myAlgParams,"outputarray");
1189 if(LALInferenceCheckVariable(&myAlgParams,"N_outputarray")) LALInferenceRemoveVariable(&myAlgParams,"N_outputarray");
1190 if(LALInferenceCheckVariable(&myAlgParams,"outfile")) LALInferenceRemoveVariable(&myAlgParams,"outfile");
1191 LALInferenceRemoveVariable(&myAlgParams,"Nmcmc");
1193
1194 LALInferenceSortVariablesByName(&myCurrentParams);
1195 runState->algorithmParams=&myAlgParams;
1196 threadState->currentParams=&myCurrentParams;
1197 LALInferenceVariables **livePoints=runState->livePoints;
1198 UINT4 Nlive = *(INT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nlive");
1199 UINT4 BAILOUT=100; /* this should be the same as the bailout in the sampler */
1200 REAL8 accept=0.0;
1201 /* We can record write the MCMC chain to a file too */
1202 ppt=LALInferenceGetProcParamVal(runState->commandLine,"--acf-chainfile");
1203 if(ppt){
1204 sprintf(chainfilename,"%s.%i",ppt->value,__chainfile_iter);
1205 chainfile=fopen(chainfilename,"w");
1206 LALInferenceCopyVariables(livePoints[0],&myCurrentParams);
1207 LALInferenceSortVariablesByName(&myCurrentParams);
1208 for(this=myCurrentParams.head;this;this=this->next) fprintf(chainfile,"%s ",this->name);
1209 fprintf(chainfile,"\n");
1211 }
1212 ppt=LALInferenceGetProcParamVal(runState->commandLine,"--acf-file");
1213 if(ppt){
1214 sprintf(acf_file_name,"%s.%i",ppt->value,__chainfile_iter);
1215 acffile=fopen(acf_file_name,"w");
1216 }
1218 do{ /* Pick a random sample that isn't trapped in some corner*/
1219 UINT4 idx=gsl_rng_uniform_int(runState->GSLrandom,Nlive);
1220 /* Copy the variable to avoid over-writing one of the live points */
1221 LALInferenceCopyVariables(livePoints[idx],&myCurrentParams);
1222 threadState->currentParams=&myCurrentParams;
1223 i=0;
1224 do {
1225 evolve(runState);
1226 i++;
1227 accept=*(REAL8*)LALInferenceGetVariable(runState->algorithmParams,"accept_rate");
1228 }
1229 while(accept==0.0 && i<BAILOUT);
1230 }
1231 while(0.==accept);
1232 /* log the first sample*/
1234 /* Evolve the initial sample (i starts at 1)*/
1235 for(i=1;i<max_iterations;i++)
1236 {
1237 evolve(runState);
1239 }
1240
1241 /* Get the location of the sample array */
1242 LALInferenceVariables **variables_array=*(LALInferenceVariables ***)LALInferenceGetVariable(runState->algorithmParams,"outputarray");
1243
1244 /* Convert to a 2D array for ACF calculation */
1245 data_array=XLALCalloc(nPar,sizeof(REAL8 *));
1246 acf_array=XLALCalloc(nPar,sizeof(REAL8 *));
1247 for (i=0;i<(UINT4)nPar;i++){
1248 data_array[i]=XLALCalloc(max_iterations,sizeof(REAL8));
1249 acf_array[i]=XLALCalloc(max_iterations/2,sizeof(REAL8));
1250 }
1251 /* Measure autocorrelation in each dimension */
1252 /* Not ideal, should be measuring something like the det(autocorrelation-crosscorrelation matrix) */
1253 for (i=0;i<max_iterations;i++){
1254 for(j=0;j<nPar;j++) data_array[j][i]=*(REAL8 *)LALInferenceGetVariable(variables_array[i],param_names[j]);
1255 }
1256 this=myCurrentParams.head;
1257 for(i=0;i<(UINT4)nPar;i++){
1258 /* Subtract the mean */
1259 REAL8 this_mean = gsl_stats_mean(&data_array[i][0], 1, max_iterations);
1260 for(j=0;j<max_iterations;j++) data_array[i][j]-=this_mean;
1261 ACL=1;
1262 int startflag=1;
1263 ACF=1.;
1264 /* Use GSL to compute the ACF */
1265 for(UINT4 lag=0;ACF>=ACF_TOLERANCE&&lag<max_iterations/2;lag++){
1266 ACF=(REAL8) gsl_stats_correlation(&data_array[i][0], 1, &data_array[i][lag], 1, max_iterations-lag);
1267 if(isnan(ACF)) ACF=1.;
1268 acf_array[i][lag]=ACF;
1269 ACL+=2.0*ACF;
1270 if((ACF<ACF_TOLERANCE && startflag) || lag==max_iterations/2-1){
1271 startflag=0;
1272 ACL*=(REAL8)thinning;
1273 if(ACL>max) max=ACL;
1275 break;
1276 }
1277 }
1278 if(LALInferenceCheckVariable(runState->algorithmParams,"verbose")) fprintf(stdout,"%s: mean= %lf, ACL=%lf\n",param_names[i],this_mean,ACL);
1279 do{this=this->next;}while(this && (this->vary==LALINFERENCE_PARAM_FIXED || this->vary==LALINFERENCE_PARAM_OUTPUT || this->type!=LALINFERENCE_REAL8_t));
1280 }
1281 if(acffile){
1282 /* Write out the ACF */
1283 for(j=0;j<(UINT4)nPar;j++) fprintf(acffile,"%s ",param_names[j]);
1284 fprintf(acffile,"\n");
1285 for(i=0;i<max_iterations/2;i++){
1286 for(j=0;j<(UINT4)nPar;j++) fprintf(acffile,"%f ",acf_array[j][i]);
1287 fprintf(acffile,"\n");
1288 }
1289 }
1290
1291 /* Cache the samples */
1292 if(LALInferenceCheckVariable(oldAlgParams,"proposalcache"))
1293 {
1294 INT4 *Ncache=(INT4 *)LALInferenceGetVariable(oldAlgParams,"proposalcachesize");
1295 LALInferenceVariables **cache_ptr=(LALInferenceVariables **)LALInferenceGetVariable(oldAlgParams,"proposalcache");
1296 LALInferenceVariables *cache=*cache_ptr;
1297 UINT4 Nnew=max_iterations/(UINT4)(max/thinning);
1298 INT4 stride=max/thinning;
1299 if(LALInferenceCheckVariable(runState->algorithmParams,"verbose")) fprintf(stderr,"Caching %i samples\n",Nnew);
1300
1301 /* Copy independent samples */
1302 REAL8 oldLogL=-INFINITY;
1303 for(i=stride,j=*Ncache;j<Nnew+*Ncache&&i<max_iterations;i+=stride,j++)
1304 {
1305 REAL8 newlogL=*(REAL8 *)LALInferenceGetVariable(variables_array[i],"logL");
1306 if(newlogL==oldLogL) {j--; continue;}
1308 if(!cache) fprintf(stderr,"ERROR!!! Could not resize cache to %i!\n",j+1);
1309 memset(&(cache[j]),0,sizeof(LALInferenceVariables));
1310 LALInferenceCopyVariables(variables_array[i],&(cache[j]));
1311 oldLogL=newlogL;
1312 }
1313
1314 /* Update the state variables */
1315 *Ncache=j;
1316 *cache_ptr=cache;
1317 }
1318
1319 /* Clean up */
1320 for(i=0;i<(UINT4)nPar;i++) {XLALFree(data_array[i]); XLALFree(acf_array[i]);}
1321 for (i=0;i<max_iterations;i++){
1322 LALInferenceClearVariables(variables_array[i]);
1323 XLALFree(variables_array[i]);
1324 }
1325 XLALFree(variables_array);
1326 XLALFree(data_array); XLALFree(acf_array);
1327 LALInferenceClearVariables(&myAlgParams);
1328 LALInferenceClearVariables(&myCurrentParams);
1329 threadState->currentParams=oldCurrentParams;
1330 runState->algorithmParams=oldAlgParams;
1331 if(chainfile) fclose(chainfile);
1332 if(acffile) fclose(acffile);
1333 XLALFree(param_names);
1334 return(acls);
1335}
1336
1337/* Perform one MCMC iteration on runState->currentParams. Return 1 if accepted or 0 if not */
1339{
1340 /* Single threaded here */
1341 LALInferenceThreadState * threadState=&runState->threads[0];
1342 UINT4 outOfBounds=0;
1343 UINT4 adaptProp=0;
1344 //LALInferenceVariables tempParams;
1345 REAL8 logProposalRatio=0.0;
1346 //LALInferenceVariables *oldParams=&tempParams;
1347 LALInferenceVariables proposedParams;
1348 memset(&proposedParams,0,sizeof(proposedParams));
1349 REAL8 logLmin=*(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"logLmin");
1350 REAL8 thislogL=-INFINITY;
1351 UINT4 accepted=0;
1352
1353 REAL8 logPriorOld=*(REAL8 *)LALInferenceGetVariable(threadState->currentParams,"logPrior");
1354 if(adaptProp)
1355 {
1356 thislogL=runState->likelihood(threadState->currentParams,runState->data,threadState->model);
1357 if (logLmin<thislogL) outOfBounds=0;
1358 }
1359 LALInferenceCopyVariables(threadState->currentParams,&proposedParams);
1360
1361 logProposalRatio = threadState->proposal(threadState,threadState->currentParams,&proposedParams);
1362 REAL8 logPriorNew=runState->prior(runState, &proposedParams, threadState->model);
1363 if(isinf(logPriorNew) || isnan(logPriorNew) || log(gsl_rng_uniform(runState->GSLrandom)) > (logPriorNew-logPriorOld) + logProposalRatio)
1364 {
1365 /* Reject - don't need to copy new params back to currentParams */
1366 /*LALInferenceCopyVariables(oldParams,runState->currentParams); */
1367 }
1368 else {
1369 accepted=1;
1370 //printf("Accepted line %i\n",__LINE__);
1371 LALInferenceCopyVariables(&proposedParams,threadState->currentParams);
1372 LALInferenceSetVariable(threadState->currentParams,"logPrior",&logPriorNew);
1373 }
1374 LALInferenceClearVariables(&proposedParams);
1375
1376 if((!outOfBounds) && adaptProp)
1377 {
1378 thislogL=runState->likelihood(threadState->currentParams,runState->data,threadState->model);
1379 if(logLmin<thislogL) threadState->accepted = accepted;
1380 else threadState->accepted=accepted=0;
1381 LALInferenceUpdateAdaptiveJumps(threadState, 0.35);
1382 }
1383 //printf("logLnew = %lf, logPriorNew = %lf, logProposalRatio = %lf\n",thislogL,logPriorNew,logProposalRatio);
1384 threadState->accepted=accepted;
1386 //printf("Accepted = %i\n",accepted);
1387 return(accepted);
1388}
1389
1390/* Sample the prior N times, returns number of acceptances */
1392{
1393 UINT4 i=0;
1394 UINT4 Naccepted=0;
1395 for(i=0;i<N;i++) Naccepted+=LALInferenceMCMCSamplePrior(runState);
1396 return(Naccepted);
1397}
1398
1400{
1401 LALInferenceVariableItem *proposeIterator = params->head;
1402 UINT4 j=0,i=0;
1403 UINT4 N=eigenvectors->size1;
1404
1405 if(!*projection) *projection=XLALCreateREAL8Vector(N);
1406 if((*projection)->length==0) *projection=XLALCreateREAL8Vector(N);
1407
1408
1409 if (proposeIterator == NULL) {
1410 fprintf(stderr, "Bad proposed params in %s, line %d\n",
1411 __FILE__, __LINE__);
1412 exit(1);
1413 }
1414 for(i=0;i<N;i++){
1415 j=0;
1416 proposeIterator = params->head;
1417 (*projection)->data[i]=0.0;
1418 do {
1419 if (proposeIterator->vary != LALINFERENCE_PARAM_FIXED && proposeIterator->vary != LALINFERENCE_PARAM_OUTPUT) {
1420 (*projection)->data[i]+= *((REAL8 *)proposeIterator->value) * gsl_matrix_get(eigenvectors, j, i);
1421 j++;
1422 }
1423 }while ((proposeIterator = proposeIterator->next) != NULL && j < N);
1424 }
1425
1426}
1427
1428/* Cache wrapper around another sampler */
1430{
1431 INT4 Naccept=0;
1432 /* Single thread here */
1433 LALInferenceThreadState *threadState = &runState->threads[0];
1434 if(!LALInferenceCheckVariable(runState->algorithmParams,"proposalcache") || !LALInferenceCheckVariable(runState->algorithmParams,"proposalcachesize"))
1435 {
1436 fprintf(stderr,"Adding cache variables in the sampler\n");
1437 /* Add space for the proposal cache */
1439 INT4 newNcache=0;
1442 }
1443
1444 INT4 *Ncache=(INT4 *)LALInferenceGetVariable(runState->algorithmParams,"proposalcachesize");
1446 LALInferenceVariables *cache=*cache_ptr;
1447
1448 if(*Ncache==0 || LALInferenceGetProcParamVal(runState->commandLine,"--no-cache")){
1449 Naccept = LALInferenceNestedSamplingSloppySample(runState);
1450 return Naccept;
1451 }
1452 REAL8 logL=-INFINITY;
1453 REAL8 logLmin=*(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"logLmin");
1454
1455 /* Draw the last sample from the cache and reduce the size of the cache by one
1456 until we find one that has a high enough likelihood */
1457 do {
1458 LALInferenceVariables *new=&(cache[*Ncache-1]);
1459 logL=*(REAL8 *)LALInferenceGetVariable(new,"logL");
1460 if(logL>logLmin){
1461 threadState->currentLikelihood=logL;
1462 LALInferenceCopyVariables(new,threadState->currentParams);
1463 }
1465 new=NULL;
1466 cache=XLALRealloc(cache,sizeof(LALInferenceVariables) * (*Ncache-1));
1467 (*Ncache)--;
1468 } while (logL<=logLmin&&*Ncache>0);
1469
1470 *cache_ptr=cache;
1471 /* If we didn't get any acceptable samples, call the main sampler */
1472 if(*Ncache==0 && logL<=logLmin)
1473 {
1474 Naccept = LALInferenceNestedSamplingSloppySample(runState);
1475 }
1476
1477 return Naccept;
1478}
1479
1480/* Sample the limited prior distribution using the MCMC method as usual, but
1481 only check the likelihood bound x fraction of the time. Always returns a fulled checked sample.
1482 x=LALInferenceGetVariable(runState->algorithmParams,"sloppyfraction")
1483 */
1484
1486{
1487 LALInferenceVariables oldParams;
1488 /* Single thread here */
1489 LALInferenceThreadState *threadState = &runState->threads[0];
1490 LALInferenceIFOData *data=runState->data;
1491 REAL8 tmp;
1492 REAL8 Target=0.3;
1493 char tmpName[320];
1494 REAL8 logLold=*(REAL8 *)LALInferenceGetVariable(threadState->currentParams,"logL");
1495 memset(&oldParams,0,sizeof(oldParams));
1496 LALInferenceCopyVariables(threadState->currentParams,&oldParams);
1497 REAL8 logLmin=*(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"logLmin");
1498 UINT4 Nmcmc=*(UINT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nmcmc");
1499 REAL8 maxsloppyfraction=((REAL8)Nmcmc-1)/(REAL8)Nmcmc ;
1500 REAL8 sloppyfraction=maxsloppyfraction/2.0;
1501 REAL8 minsloppyfraction=0.;
1502 if(Nmcmc==1) maxsloppyfraction=minsloppyfraction=0.0;
1503 if (LALInferenceCheckVariable(runState->algorithmParams,"sloppyfraction"))
1504 sloppyfraction=*(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"sloppyfraction");
1505 UINT4 mcmc_iter=0,Naccepted=0,sub_accepted=0;
1506 UINT4 sloppynumber=(UINT4) (sloppyfraction*(REAL8)Nmcmc);
1507 UINT4 testnumber=Nmcmc-sloppynumber;
1508 /* +1 for the last iteration which we do check */
1509 UINT4 subchain_length=(sloppynumber/testnumber) +1;
1510 REAL8 logLnew=0.0;
1511 UINT4 sub_iter=0;
1512 UINT4 ifo=0;
1513 REAL8 counter=1.;
1514 UINT4 BAILOUT=100*testnumber; /* If no acceptance after 100 tries, will exit and the sampler will try a different starting point */
1515 const char *extra_names[]={"logL","optimal_snr","matched_filter_snr","deltalogL"}; /* Names for parameters to be stripped when sampling prior */
1516 UINT4 Nnames = 4;
1517 if ( Nmcmc ){
1518 do{
1519 counter=counter-1.;
1520 subchain_length=0;
1521 for(UINT4 i=0;i<Nnames;i++)
1522 {
1523 if(LALInferenceCheckVariable(threadState->currentParams,extra_names[i]))
1524 LALInferenceRemoveVariable(threadState->currentParams,extra_names[i]);
1525 }
1526 /* Draw an independent sample from the prior */
1527 do{
1528
1529 sub_accepted+=LALInferenceMCMCSamplePrior(runState);
1530 subchain_length++;
1531 counter+=(1.-sloppyfraction);
1532 }while(counter<1);
1533 /* Check that there was at least one accepted point */
1534 if(sub_accepted==0) {
1535 sub_iter+=subchain_length;
1536 mcmc_iter++;
1537 LALInferenceCopyVariables(&oldParams,threadState->currentParams);
1538 threadState->currentLikelihood=logLold;
1539 continue;
1540 }
1541 mcmc_iter++;
1542 sub_iter+=subchain_length;
1543 if(isfinite(logLmin)) logLnew=runState->likelihood(threadState->currentParams,runState->data,threadState->model);
1544 if(logLnew>logLmin || isinf(logLmin)) /* Accept */
1545 {
1546 Naccepted++;
1547 /* Update information to pass back out */
1549 if(LALInferenceCheckVariable(runState->algorithmParams,"logZnoise")){
1550 tmp=logLnew-*(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"logZnoise");
1552 }
1553 ifo=0;
1554 data=runState->data;
1555 while(data)
1556 {
1557 if(!threadState->model->ifo_loglikelihoods) break;
1558 tmp=threadState->model->ifo_loglikelihoods[ifo] - data->nullloglikelihood;
1559 sprintf(tmpName,"deltalogl%s",data->name);
1561 ifo++;
1562 data=data->next;
1563 }
1564 LALInferenceCopyVariables(threadState->currentParams,&oldParams);
1565 logLold=logLnew;
1566 threadState->currentLikelihood=logLnew;
1567 }
1568 else /* reject */
1569 {
1570 LALInferenceCopyVariables(&oldParams,threadState->currentParams);
1571 threadState->currentLikelihood=logLold;
1572 }
1573 }while((mcmc_iter<testnumber||threadState->currentLikelihood<=logLmin||Naccepted==0)&&(mcmc_iter<BAILOUT));
1574 }
1575 /* Make sure likelihood is filled in if it wasn't done during sampling */
1576 if(logLnew==0.0){
1577 logLnew=runState->likelihood(threadState->currentParams,runState->data,threadState->model);
1578 threadState->currentLikelihood=logLnew;
1580 if(LALInferenceCheckVariable(runState->algorithmParams,"logZnoise")){
1581 tmp=logLnew-*(REAL8 *)LALInferenceGetVariable(runState->algorithmParams,"logZnoise");
1583 }
1584 ifo=0;
1585 data=runState->data;
1586 while(data && threadState->model->ifo_loglikelihoods)
1587 {
1588 tmp=threadState->model->ifo_loglikelihoods[ifo] - data->nullloglikelihood;
1589 sprintf(tmpName,"deltalogl%s",data->name);
1591 ifo++;
1592 data=data->next;
1593 }
1594 }
1595
1596 /* Compute some statistics for information */
1597 REAL8 sub_accept_rate=(REAL8)sub_accepted/(REAL8)sub_iter;
1598 REAL8 accept_rate=(REAL8)Naccepted/(REAL8)testnumber;
1599 LALInferenceSetVariable(runState->algorithmParams,"accept_rate",&accept_rate);
1600 LALInferenceSetVariable(runState->algorithmParams,"sub_accept_rate",&sub_accept_rate);
1601 /* Adapt the sloppy fraction toward target acceptance of outer chain */
1602 if(isfinite(logLmin)){
1603 if((REAL8)accept_rate>Target) { sloppyfraction+=5.0/(REAL8)Nmcmc;}
1604 else { sloppyfraction-=5.0/(REAL8)Nmcmc;}
1605 if(sloppyfraction>maxsloppyfraction) sloppyfraction=maxsloppyfraction;
1606 if(sloppyfraction<minsloppyfraction) sloppyfraction=minsloppyfraction;
1607
1608 LALInferenceSetVariable(runState->algorithmParams,"sloppyfraction",&sloppyfraction);
1609 }
1610 /* Cleanup */
1611 LALInferenceClearVariables(&oldParams);
1612
1613 return Naccepted;
1614}
1615
1616
1617/* Evolve nested sampling algorithm by one step, i.e.
1618 evolve runState->currentParams to a new point with higher
1619 likelihood than currentLikelihood. Uses the MCMC method with sloppy sampling.
1620 */
1622{
1624 return Naccept;
1625}
1626
1628 /* Set up initial basket of live points, drawn from prior,
1629 by copying runState->currentParams to all entries in the array*/
1630
1631 /* Single thread here */
1632 LALInferenceThreadState *threadState = &runState->threads[0];
1633
1634 UINT4 Nlive=(UINT4)*(INT4 *)LALInferenceGetVariable(runState->algorithmParams,"Nlive");
1635 UINT4 i;
1636 REAL8Vector *logLs;
1637 REAL8 logPrior=0.0;
1638
1639 /* Allocate the array */
1640 /* runState->livePoints=XLALCalloc(Nlive,sizeof(LALVariables *)); */
1641 runState->livePoints=XLALCalloc(Nlive,sizeof(LALInferenceVariables *));
1642 if(runState->livePoints==NULL)
1643 {
1644 fprintf(stderr,"Unable to allocate memory for %i live points\n",Nlive);
1645 exit(1);
1646 }
1647 threadState->differentialPoints=runState->livePoints;
1648 threadState->differentialPointsLength=(size_t) Nlive;
1649 logLs=XLALCreateREAL8Vector(Nlive);
1650
1652 fprintf(stdout,"Sprinkling %i live points, may take some time\n",Nlive);
1653 LALInferenceVariables *curParsBackup=threadState->currentParams;
1654 for(i=0;i<Nlive;i++)
1655 {
1656 /* Clone the currentParams into LivePoints[i] */
1657 runState->livePoints[i]=XLALCalloc(1,sizeof(LALInferenceVariables));
1658 /* Copy the param structure */
1659 LALInferenceCopyVariables(threadState->currentParams,runState->livePoints[i]);
1660
1661 /* Sprinkle the varying points among prior */
1662 do{
1663 LALInferenceDrawFromPrior( runState->livePoints[i], runState->priorArgs, runState->GSLrandom );
1664 logPrior=runState->prior(runState,runState->livePoints[i],threadState->model);
1665 }while(isinf(logPrior) || isnan(logPrior));
1666 /* Populate log likelihood */
1667 logLs->data[i]=runState->likelihood(runState->livePoints[i],runState->data,threadState->model);
1670 }
1671 threadState->currentParams=curParsBackup;
1672 if(!threadState->currentParams) threadState->currentParams=XLALCalloc(1,sizeof(LALInferenceVariables));
1673}
1674
1675
1677{
1678 /* Single thread here */
1679 LALInferenceThreadState *threadState = &runState->threads[0];
1680 gsl_matrix *eVectors=NULL;
1681 gsl_vector *eValues =NULL;
1682 REAL8Vector *eigenValues=NULL;
1683 /* Check for existing covariance matrix */
1684 gsl_matrix **cvm=NULL;
1685 if(LALInferenceCheckVariable(threadState->proposalArgs,"covarianceMatrix"))
1686 LALInferenceRemoveVariable(threadState->proposalArgs,"covarianceMatrix");
1687// cvm=(gsl_matrix **)LALInferenceGetVariable(threadState->proposalArgs,"covarianceMatrix");
1688 cvm=XLALCalloc(1,sizeof(gsl_matrix *));
1689
1690 /* Add the covariance matrix for proposal distribution */
1691 UINT4 *Nlive=LALInferenceGetVariable(runState->algorithmParams,"Nlive");
1692 /* Sort the variables to ensure consistent order */
1693 for(UINT4 i=0;i<*Nlive;i++) LALInferenceSortVariablesByName(runState->livePoints[i]);
1694 LALInferenceNScalcCVM(cvm,runState->livePoints,*Nlive);
1695 UINT4 N=(*cvm)->size1;
1696
1697
1698 /* Check for the eigenvectors and values */
1699 if(LALInferenceCheckVariable(threadState->proposalArgs,"covarianceEigenvectors"))
1700 eVectors=*(gsl_matrix **)LALInferenceGetVariable(threadState->proposalArgs,"covarianceEigenvectors");
1701 else
1702 eVectors=gsl_matrix_alloc(N,N);
1703
1704 if(LALInferenceCheckVariable(threadState->proposalArgs,"covarianceEigenvalues"))
1705 eigenValues=*(REAL8Vector **)LALInferenceGetVariable(threadState->proposalArgs,"covarianceEigenvalues");
1706 else
1707 eigenValues=XLALCreateREAL8Vector(N);
1708
1709 /* Set up eigenvectors and eigenvalues. */
1710 gsl_matrix *covCopy = gsl_matrix_alloc(N,N);
1711 eValues = gsl_vector_alloc(N);
1712 gsl_eigen_symmv_workspace *ws = gsl_eigen_symmv_alloc(N);
1713 int gsl_status;
1714 gsl_matrix_memcpy(covCopy, *cvm);
1715
1716 if ((gsl_status = gsl_eigen_symmv(covCopy, eValues, eVectors, ws)) != GSL_SUCCESS) {
1717 XLALPrintError("Error in gsl_eigen_symmv (in %s, line %d): %d: %s\n", __FILE__, __LINE__, gsl_status, gsl_strerror(gsl_status));
1719 }
1720
1721 for (UINT4 i = 0; i < N; i++) {
1722 eigenValues->data[i] = gsl_vector_get(eValues,i);
1723 }
1724
1725 if(!LALInferenceCheckVariable(threadState->proposalArgs,"covarianceEigenvectors"))
1726 LALInferenceAddVariable(threadState->proposalArgs, "covarianceEigenvectors", &eVectors, LALINFERENCE_gslMatrix_t, LALINFERENCE_PARAM_FIXED);
1727 if(!LALInferenceCheckVariable(threadState->proposalArgs,"covarianceEigenvalues"))
1728 LALInferenceAddVariable(threadState->proposalArgs, "covarianceEigenvalues", &eigenValues, LALINFERENCE_REAL8Vector_t, LALINFERENCE_PARAM_FIXED);
1730
1731 gsl_matrix_free(covCopy);
1732 gsl_vector_free(eValues);
1733 gsl_eigen_symmv_free(ws);
1734 XLALFree(cvm);
1735}
1736
1737
1739{
1742
1743 for(INT4 i=0;i<N;i++)
1744 {
1745 if(!thread->differentialPoints[i]) thread->differentialPoints[i]=XLALCalloc(1,sizeof(LALInferenceVariables));
1747 }
1749 return(XLAL_SUCCESS);
1750}
struct tagLALH5File LALH5File
struct tagLALH5Dataset LALH5Dataset
#define TOLERANCE
int LALInferenceCheckNonEmptyFile(char *filename)
Returns 1 if a non-empty file exists, 0 otherwise.
int LALInferencePrintCheckpointFileInfo(char *filename)
Prints the size of the file.
int LALInferenceH5DatasetToVariablesArray(LALH5Dataset *dataset, LALInferenceVariables ***varsArray, UINT4 *N)
int LALInferenceH5VariablesArrayToDataset(LALH5File *h5file, LALInferenceVariables *const *const varsArray, UINT4 N, const char *TableName)
LALH5File * LALInferenceH5CreateGroupStructure(LALH5File *h5file, const char *codename, const char *runID)
Create a HDF5 heirarchy in the given LALH5File reference /codename/runID/ Returns a LALH5File pointer...
const char LALInferenceHDF5NestedSamplesDatasetName[]
static int WriteNSCheckPointH5(char *filename, LALInferenceRunState *runState, NSintegralState *s)
static int ReadNSCheckPointH5(char *filename, LALInferenceRunState *runState, NSintegralState *s)
static void catch_interrupt(UNUSED int sig, UNUSED siginfo_t *siginfo, UNUSED void *context)
static NSintegralState * initNSintegralState(UINT4 Nruns, UINT4 Nlive)
#define ACF_TOLERANCE
static int __chainfile_iter
static REAL8 mean(REAL8 *array, int N)
static void catch_alarm(UNUSED int sig, UNUSED siginfo_t *siginfo, UNUSED void *context)
static REAL8 incrementEvidenceSamples(gsl_rng *GSLrandom, UINT4 Nlive, REAL8 logL, NSintegralState *s)
Update the internal state of the integrator after receiving the lowest logL value logL.
static void printAdaptiveJumpSizes(FILE *file, LALInferenceThreadState *threadState)
static UINT4 UpdateNMCMC(LALInferenceRunState *runState)
static double logsubexp(double a, double b)
static int syncLivePointsDifferentialPoints(LALInferenceRunState *state, LALInferenceThreadState *thread)
Sync the live points to the differential evolution buffer.
static int _loadNSintegralStateH5(LALH5File *group, NSintegralState *s)
static volatile sig_atomic_t __ns_saveStateFlag
static void install_resume_handler(int checkpoint_exit)
Install the signal handlers for checkpointing.
static volatile sig_atomic_t __ns_exitFlag
#define MAX_MCMC
static struct itimerval checkpoint_timer
static REAL8 LALInferenceNSSample_logt(int Nlive, gsl_rng *RNG)
static void LALInferenceNScalcCVM(gsl_matrix **cvm, LALInferenceVariables **Live, UINT4 Nlive)
Calculate covariance matrix from a collection of live points.
static void SetupEigenProposals(LALInferenceRunState *runState)
static int _saveNSintegralStateH5(LALH5File *group, NSintegralState *s)
Utility functions for the resume functionality.
static int CheckOutputFileContents(char *filename)
double LALInferenceMassDistancePriorVolume(LALInferenceRunState *state)
#define max(a, b)
ProcessParamsTable * ppt
int j
ProcessParamsTable * ptr
int k
#define fprintf
const char *const name
int s
const double H
sigmaKerr data[0]
REAL8Vector * XLALH5FileReadREAL8Vector(LALH5File *file, const char *name)
int XLALH5FileWriteREAL8Vector(LALH5File *file, const char *name, REAL8Vector *vector)
void XLALH5FileClose(LALH5File UNUSED *file)
void XLALH5DatasetFree(LALH5Dataset UNUSED *dset)
LALH5File * XLALH5GroupOpen(LALH5File UNUSED *file, const char UNUSED *name)
int XLALH5FileQueryScalarAttributeValue(void UNUSED *value, LALH5File UNUSED *file, const char UNUSED *key)
int XLALH5FileAddScalarAttribute(LALH5File UNUSED *file, const char UNUSED *key, const void UNUSED *value, LALTYPECODE UNUSED dtype)
LALH5File * XLALH5FileOpen(const char UNUSED *path, const char UNUSED *mode)
int XLALH5FileAddStringAttribute(LALH5File UNUSED *file, const char UNUSED *key, const char UNUSED *value)
LALH5Dataset * XLALH5DatasetRead(LALH5File UNUSED *file, const char UNUSED *name)
#define LAL_PI
double REAL8
uint32_t UINT4
int32_t INT4
LAL_D_TYPE_CODE
LAL_I4_TYPE_CODE
LAL_U4_TYPE_CODE
LALNDEBUG
LALMSGLVL2
LALMSGLVL1
LALMSGLVL3
REAL8 LALInferenceAngularDistance(REAL8 a1, REAL8 a2)
Calculate shortest angular distance between a1 and a2 (modulo 2PI)
void LALInferencePrintSample(FILE *fp, LALInferenceVariables *sample)
Output the sample to file *fp, in ASCII format.
Definition: LALInference.c:923
int LALInferenceFprintParameterHeaders(FILE *out, LALInferenceVariables *params)
Print the parameter names to a file as a tab-separated ASCII line.
void LALInferenceAddVariable(LALInferenceVariables *vars, const char *name, const void *value, LALInferenceVariableType type, LALInferenceParamVaryType vary)
Add a variable named name to vars with initial value referenced by value.
Definition: LALInference.c:395
REAL8 LALInferenceGetREAL8Variable(LALInferenceVariables *vars, const char *name)
void LALInferenceRemoveVariable(LALInferenceVariables *vars, const char *name)
Remove name from vars Frees the memory for the name structure and its contents.
Definition: LALInference.c:450
char * LALInferencePrintCommandLine(ProcessParamsTable *procparams)
Output the command line based on the ProcessParamsTable procparams.
void LALInferenceCopyVariables(LALInferenceVariables *origin, LALInferenceVariables *target)
Deep copy the variables from one to another LALInferenceVariables structure.
Definition: LALInference.c:558
void * LALInferenceGetVariable(const LALInferenceVariables *vars, const char *name)
Return a pointer to the memory the variable vars is stored in specified by name User must cast this p...
Definition: LALInference.c:238
INT4 LALInferenceGetINT4Variable(LALInferenceVariables *vars, const char *name)
ProcessParamsTable * LALInferenceGetProcParamVal(ProcessParamsTable *procparams, const char *name)
Returns the element of the process params table with "name".
int LALInferencePrintProposalStatsHeader(FILE *fp, LALInferenceProposalCycle *cycle)
Output proposal statistics header to file *fp.
void LALInferencePrintProposalStats(FILE *fp, LALInferenceProposalCycle *cycle)
Output proposal statistics to file *fp.
void LALInferenceClearVariables(LALInferenceVariables *vars)
Delete the variables in this structure.
Definition: LALInference.c:532
void LALInferenceSortVariablesByName(LALInferenceVariables *vars)
Sorts the variable structure by name.
INT4 LALInferenceSanityCheck(LALInferenceRunState *state)
Sanity check the structures in the given state.
INT4(* LALInferenceEvolveOneStepFunction)(struct tagLALInferenceRunState *runState)
Perform one step of an algorithm, replaces runState ->currentParams.
Definition: LALInference.h:491
void LALInferenceLogSampleToArray(LALInferenceVariables *algorithmParams, LALInferenceVariables *vars)
Append the sample to an array which can be later processed by the user.
int LALInferenceCheckVariable(LALInferenceVariables *vars, const char *name)
Checks for name being present in vars returns 1(==true) or 0.
Definition: LALInference.c:525
void LALInferenceSetVariable(LALInferenceVariables *vars, const char *name, const void *value)
Set a variable named name in vars with a value.
Definition: LALInference.c:352
@ LALINFERENCE_PARAM_OUTPUT
A parameter that never changes, functions should respect this.
Definition: LALInference.h:131
@ LALINFERENCE_PARAM_FIXED
A parameter that is cyclic, such as an angle between 0 and 2pi.
Definition: LALInference.h:130
@ LALINFERENCE_PARAM_CIRCULAR
A parameter that simply has a maximum and a minimum.
Definition: LALInference.h:129
@ LALINFERENCE_PARAM_LINEAR
Definition: LALInference.h:128
@ LALINFERENCE_REAL8_t
Definition: LALInference.h:109
@ LALINFERENCE_INT4_t
Definition: LALInference.h:105
@ LALINFERENCE_REAL8Vector_t
Definition: LALInference.h:113
@ LALINFERENCE_void_ptr_t
Definition: LALInference.h:119
@ LALINFERENCE_gslMatrix_t
Definition: LALInference.h:112
REAL8 LALInferenceNullLogLikelihood(LALInferenceIFOData *data)
Identical to LALInferenceFreqDomainNullLogLikelihood, but returns the likelihood of a null template.
void LALInferenceNestedSamplingAlgorithm(LALInferenceRunState *runState)
NestedSamplingAlgorithm implements the nested sampling algorithm, see e.g.
INT4 LALInferenceNestedSamplingOneStep(LALInferenceRunState *runState)
A single iteration of the NS algorithm.
void LALInferenceSetupLivePointsArray(LALInferenceRunState *runState)
Setup the live points by calling runState->initVariables on each of them if it is specified.
LALInferenceVariables * LALInferenceComputeAutoCorrelation(LALInferenceRunState *runState, UINT4 max_iterations, LALInferenceEvolveOneStepFunction evolve)
Compute the autocorrelation length from the sampler at the current global iteration.
INT4 LALInferenceNestedSamplingSloppySample(LALInferenceRunState *runState)
Sample the limited prior distribution using the MCMC method as usual, but run a sub-chain of x iterat...
INT4 LALInferenceNestedSamplingCachedSampler(LALInferenceRunState *runState)
void LALInferenceNestedSamplingAlgorithmInit(LALInferenceRunState *runState)
Initialise the nested sampling algorithm by reading from the commandLine and setting up algorithmPara...
UINT4 LALInferenceMCMCSamplePriorNTimes(LALInferenceRunState *runState, UINT4 N)
Sample the prior N times, returns number of acceptances.
void LALInferenceProjectSampleOntoEigenvectors(LALInferenceVariables *params, gsl_matrix *eigenvectors, REAL8Vector **projection)
Project the sample in params onto the eigenvectors given in eigenvectors.
UINT4 LALInferenceMCMCSamplePrior(LALInferenceRunState *runState)
Perform one MCMC iteration on runState->currentParams.
void LALInferenceDrawFromPrior(LALInferenceVariables *output, LALInferenceVariables *priorArgs, gsl_rng *rdm)
Draw variables from the prior ranges.
#define ADAPTSUFFIX
void LALInferenceUpdateAdaptiveJumps(LALInferenceThreadState *thread, REAL8 targetAcceptance)
Update the adaptive proposal.
void LALInferenceTrackProposalAcceptance(LALInferenceThreadState *thread)
Update proposal statistics if tracking.
void LALInferenceSetupClusteredKDEProposalFromDEBuffer(LALInferenceThreadState *thread)
Setup a clustered-KDE proposal from the differential evolution buffer.
void LALInferenceZeroProposalStats(LALInferenceProposalCycle *cycle)
REAL8 LALInferenceCyclicProposal(LALInferenceThreadState *thread, LALInferenceVariables *currentParams, LALInferenceVariables *proposedParams)
Proposes a jump from the next proposal in the proposal cycle.
LALInferenceVariables * LALInferencePrintInjectionSample(LALInferenceRunState *runState)
Function to output a sample with logL values etc for the injection, if one is made.
void * XLALMalloc(size_t n)
void * XLALCalloc(size_t m, size_t n)
void * XLALRealloc(void *p, size_t n)
void XLALFree(void *p)
static const INT4 a
REAL8Vector * XLALCreateREAL8Vector(UINT4 length)
#define XLAL_ERROR_VOID(...)
const char * XLALErrorString(int errnum)
#define XLAL_ERROR(...)
#define XLAL_TRY(statement, errnum)
int XLALPrintError(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
int int XLALPrintWarning(const char *fmt,...) _LAL_GCC_PRINTF_FORMAT_(1
int XLALPrintProgressBar(double)
#define XLAL_IS_REAL8_FAIL_NAN(val)
XLAL_SUCCESS
XLAL_EFAILED
static double logaddexp(double x, double y)
Definition: logaddexp.h:37
help
type
LALCache * cache
Structure to contain IFO data.
Definition: LALInference.h:625
REAL8 * ifo_loglikelihoods
Network SNR at params
Definition: LALInference.h:445
Structure containing inference run state This includes pointers to the function types required to run...
Definition: LALInference.h:592
ProcessParamsTable * commandLine
Definition: LALInference.h:593
LALInferenceVariables * proposalArgs
The data from the interferometers.
Definition: LALInference.h:602
LALInferenceVariables * algorithmParams
Any special arguments for the prior function.
Definition: LALInference.h:604
struct tagLALInferenceIFOData * data
Log sample, i.e.
Definition: LALInference.h:601
LALInferenceVariables * priorArgs
Common arguments needed by proposals, to be copied to thread->cycle.
Definition: LALInference.h:603
LALInferenceThreadState * threads
Definition: LALInference.h:612
LALInferenceAlgorithm algorithm
A function that returns a new set of variables for the model.
Definition: LALInference.h:595
LALInferenceEvolveOneStepFunction evolve
The algorithm function.
Definition: LALInference.h:596
LALInferencePriorFunction prior
The algorithm's single iteration function.
Definition: LALInference.h:597
LALInferenceLogFunction logsample
The likelihood function.
Definition: LALInference.h:600
LALInferenceVariables ** livePoints
Parameters which control the running of the algorithm.
Definition: LALInference.h:605
LALInferenceLikelihoodFunction likelihood
MultiNest prior for the parameters.
Definition: LALInference.h:599
Structure containing chain-specific variables.
Definition: LALInference.h:541
LALInferenceVariables * currentParams
Prior boundaries, etc.
Definition: LALInference.h:556
LALInferenceModel * model
Cycle of proposals to call.
Definition: LALInference.h:548
INT4 accepted
This should be removed, can be given as an algorithmParams entry.
Definition: LALInference.h:574
size_t differentialPointsLength
Array of points for differential evolution.
Definition: LALInference.h:560
LALInferenceProposalFunction proposal
Step counter for this thread.
Definition: LALInference.h:546
LALInferenceProposalCycle * cycle
The proposal cycle.
Definition: LALInference.h:547
LALInferenceVariables * proposalArgs
Definition: LALInference.h:551
size_t differentialPointsSkip
Size of the differentialPoints memory block (must be >= length of differential points).
Definition: LALInference.h:566
LALInferenceVariables ** differentialPoints
Parameters proposed.
Definition: LALInference.h:559
The LALInferenceVariableItem list node structure This should only be accessed using the accessor func...
Definition: LALInference.h:154
char name[VARNAME_MAX]
Definition: LALInference.h:155
LALInferenceVariableType type
Definition: LALInference.h:157
struct tagVariableItem * next
Definition: LALInference.h:159
LALInferenceParamVaryType vary
Definition: LALInference.h:158
The LALInferenceVariables structure to contain a set of parameters Implemented as a linked list of LA...
Definition: LALInference.h:170
LALInferenceVariableItem * head
Definition: LALInference.h:171
structure holding internal state of the NS integrator
CHAR value[LIGOMETA_VALUE_MAX]
REAL8 * data
char * outfile