ADMB Documentation  -a65f1c97
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
hybmcmc.cpp
Go to the documentation of this file.
1 
8 #include <sstream>
9 using std::istringstream;
10 
11 # include <df1b2fun.h>
12 # include <adrndeff.h>
13 #include <admodel.h>
14 
15 #if defined(_MSC_VER)
16  #include <conio.h>
17 #endif
18 
19 #ifndef OPT_LIB
20  #include <cassert>
21  #include <climits>
22 #endif
23 
24 double better_rand(long int&);
25 void store_mcmc_values(const ofstream& ofs);
26 void set_labels_for_mcmc(void);
27 
28 void print_hist_data(const dmatrix& hist, const dmatrix& values,
29  const dvector& h, dvector& m, const dvector& s, const dvector& parsave,
30  long int iseed, double size_scale);
31 
32 int minnz(const dvector& x);
33 int maxnz(const dvector& xa);
34 
35 void read_hessian_matrix_and_scale1(int nvar, const dmatrix& _SS,double s,
36  int mcmc2_flag);
37 
38 int read_hist_data(const dmatrix& hist, const dvector& h,
39  dvector& m, const dvector& s, const dvector& parsave,long int& iseed,
40  const double& size_scale);
41 
42 void make_preliminary_hist(const dvector& s, const dvector& m,int nsim,
43  const dmatrix& values, dmatrix& hist, const dvector& h,int slots,
44  double total_spread,int probflag=0);
45 
46 void add_hist_values(const dvector& s, const dvector& m, const dmatrix& hist,
47  dvector& mcmc_values,double llc, const dvector& h,int nslots,
48  double total_spreadd,int probflag=0);
49 
50 void write_empirical_covariance_matrix(int ncor, const dvector& s_mean,
51  const dmatrix& s_covar, adstring& prog_name);
52 
53 void read_empirical_covariance_matrix(int nvar, const dmatrix& S,
54  const adstring& prog_name);
55 
56 void read_hessian_matrix_and_scale(int nvar, const dmatrix& S,
57  const dvector& pen_vector);
58 
60  const dvector& b1, dmatrix& ch, const double& wght,double pprobe,
62 
64  const dvector& b1, dmatrix& ch, const double& wght, const dvector& _y,
65  double pprobe, random_number_generator& rng);
66 
67 //void newton_raftery_bayes_estimate(double cbf,int ic, const dvector& lk,
68 //double d);
69 
71  (int feval,int iter,double fval,double gmax);
72 
73 void ad_update_function_minimizer_report(int feval,int iter,int phase,
74  double fval,double gmax,const char * cbuf);
75 void ad_update_mcmc_report(dmatrix& m,int i,int j,int ff=0);
76 void ad_update_mcmchist_report(dmatrix& mcmc_values,ivector& number_offsets,
77  dvector& mean_mcmc_values,dvector& h,int ff=0);
78 
83 void function_minimizer::hybrid_mcmc_routine(int nmcmc,int iseed0,double dscale,
84  int restart_flag)
85 {
87  uostream * pofs_psave=NULL;
88  dmatrix mcmc_display_matrix;
89  //int mcmc_save_index=1;
90  //int mcmc_wrap_flag=0;
91 
92  int on2=-1;
93  //int nvar1=0;
94  if ( (on2=option_match(ad_comm::argc,ad_comm::argv,"-nosdmcmc"))>-1)
95  {
96  //int no_sd_mcmc = 1;
97  }
98  if (mcmc2_flag==1)
99  {
101  //nvar1=initial_params::nvarcalc(); // get the number of active parameters
102  }
103 
105  {
106  cerr << " You must declare at least one object of type sdreport "
107  << endl << " to do the mcmc calculations" << endl;
108  //return;
109  }
110  {
111  //ofstream of_bf("testbf");
112 
113  //if (adjm_ptr) set_labels_for_mcmc();
114 
115  ivector number_offsets;
116  dvector lkvector;
117  //double current_bf=0;
118  //double lcurrent_bf=0;
119  //double size_scale=1.0;
120  //double total_spread=200;
121  //double total_spread=2500;
122  //uostream * pofs_sd = NULL;
123 
125  int nvar_x=initial_params::nvarcalc();
127  int nvar_re=initial_params::nvarcalc();
128 
129  int nvar=initial_params::nvarcalc(); // get the number of active parameters
130  dmatrix s_covar;
131  dvector s_mean;
132  //int ncsim=25000;
133  //int nslots=800;
134  //int nslots=3600;
135  //int initial_nsim=4800;
136  //int ntmp=0;
137  //int ncor=0;
138  //double bfsum=0;
139  //int ibfcount=0;
140  //double lbmax;
141 
142  s_covar.allocate(1,nvar,1,nvar);
143  s_mean.allocate(1,nvar);
144  s_mean.initialize();
145  s_covar.initialize();
146 
147  int ndvar=stddev_params::num_stddev_calc();
148  //int numdvar=stddev_params::num_stddev_number_calc();
149 
150  if (mcmc2_flag==0)
151  {
153  nvar=initial_params::nvarcalc(); // get the number of active parameters
154  }
155 
156  independent_variables parsave(1,nvar_re);
158 
159  dvector x(1,nvar);
160  //dvector scale(1,nvar);
161  dmatrix values;
162  //int have_hist_flag=0;
164  dvector pen_vector(1,nvar);
165  {
166  initial_params::reset(dvar_vector(x),pen_vector);
167  //cout << pen_vector << endl << endl;
168  }
169 
171  //initial_params::stddev_scale(scale,x);
173  dvector bmn(1,nvar);
174  dvector mean_mcmc_values(1,ndvar);
175  dvector s(1,ndvar);
176  dvector h(1,ndvar);
177  //dvector h;
178  dvector square_mcmc_values(1,ndvar);
179  square_mcmc_values.initialize();
180  mean_mcmc_values.initialize();
181  bmn.initialize();
182  int use_empirical_flag=0;
183  int diag_option=0;
184  //int topt=0;
185  int hybnstep=10;
186  double hybeps=0.1;
187  double _hybeps=-1.0;
188  int old_Hybrid_bounded_flag=-1;
189 
190  int on,nopt = 0;
191  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-hyeps",nopt))>-1)
192  {
193  if (!nopt)
194  {
195  cerr << "Usage -hyeps option needs number -- ignored" << endl;
196  }
197  else
198  {
199  istringstream ist(ad_comm::argv[on+1]);
200  ist >> _hybeps;
201 
202  if (_hybeps<=0)
203  {
204  cerr << "Usage -hyeps option needs positive number -- ignored\n";
205  _hybeps=-1.0;
206  }
207  }
208  }
209  if (_hybeps>0.0)
210  hybeps=_hybeps;
211 
212  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-hynstep",nopt))>-1)
213  {
214  if (nopt)
215  {
216  int iii=atoi(ad_comm::argv[on+1]);
217  if (iii < 1 )
218  {
219  cerr << " -hynstep argument must be integer between > 0 --"
220  " using default of 10" << endl;
221  }
222  else
223  {
224  hybnstep=iii;
225  }
226  }
227  }
228 
229  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcdiag"))>-1)
230  {
231  diag_option=1;
232  cout << " Setting covariance matrix to diagonal with entries " << dscale
233  << endl;
234  }
235  int mcrestart_flag=option_match(ad_comm::argc,ad_comm::argv,"-mcr");
236  dmatrix S(1,nvar,1,nvar);
237  dvector old_scale(1,nvar);
238  int old_nvar;
239  if (!diag_option)
240  {
241  int rescale_bounded_flag=0;
242  double rescale_bounded_power=0.5;
243  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcrb",nopt))>-1)
244  {
245  if (nopt)
246  {
247  int iii=atoi(ad_comm::argv[on+1]);
248  if (iii < 1 || iii > 9)
249  {
250  cerr << " -mcrb argument must be integer between 1 and 9 --"
251  " using default of 5" << endl;
252  rescale_bounded_power=0.5;
253  }
254  else
255  rescale_bounded_power=iii/10.0;
256  }
257  else
258  {
259  rescale_bounded_power=0.5;
260  }
261  rescale_bounded_flag=1;
262  }
263  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcec"))>-1)
264  {
265  use_empirical_flag=1;
266  }
267  if (use_empirical_flag)
268  {
270  }
271  else if (!rescale_bounded_flag)
272  {
273  if (mcmc2_flag==0)
274  {
275  read_covariance_matrix(S,nvar,old_Hybrid_bounded_flag,old_scale);
276  }
277  else
278  {
279  int tmp_nvar = 0;
280  adstring tmpstring = ad_comm::adprogram_name + ".bgs";
281  uistream uis((char*)(tmpstring));
282  if (!uis)
283  {
284  cerr << "error opening file " << tmpstring << endl;
285  ad_exit(1);
286  }
287  uis >> tmp_nvar;
288  if (!uis)
289  {
290  cerr << "error reading from file " << tmpstring << endl;
291  ad_exit(1);
292  }
293  if (tmp_nvar != nvar)
294  {
295  cerr << "size error reading from " << tmpstring << endl;
296  ad_exit(1);
297  }
298  uis >> S;
299  if (!uis)
300  {
301  cerr << "error reading from file " << tmpstring << endl;
302  ad_exit(1);
303  }
304  dvector tmp=read_old_scale(old_nvar);
305  old_scale=1.0;
306  old_scale(1,old_nvar)=tmp;
307  }
308  }
309  else
310  {
311  read_hessian_matrix_and_scale1(nvar,S,rescale_bounded_power,
312  mcmc2_flag);
313  //read_hessian_matrix_and_scale(nvar,S,pen_vector);
314  }
315 
316  // maybe we dont want to do this?
317  /*
318  { // scale covariance matrix for model space
319  dmatrix tmp(1,nvar,1,nvar);
320  for (int i=1;i<=nvar;i++)
321  {
322  tmp(i,i)=S(i,i)*(scale(i)*scale(i));
323  for (int j=1;j<i;j++)
324  {
325  tmp(i,j)=S(i,j)*(scale(i)*scale(j));
326  tmp(j,i)=tmp(i,j);
327  }
328  }
329  S=tmp;
330  }
331  */
332  }
333  else
334  {
335  S.initialize();
336  for (int i=1;i<=nvar;i++)
337  {
338  S(i,i)=dscale;
339  }
340  }
341 
342  // for hybrid mcmc option always save output
343  //if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcsave"))>-1)
344  if ( mcrestart_flag>-1)
345  {
346  // check that nvar is correct
347  uistream uis((char*)(ad_comm::adprogram_name + adstring(".psv")));
348  if (!uis)
349  {
350  cerr << "Error trying to open file" <<
351  ad_comm::adprogram_name + adstring(".psv") <<
352  " for mcrestart" << endl;
353  ad_exit(1);
354  } else {
355  int nv1 = 0;
356  uis >> nv1;
357  if (nv1 !=nvar) {
358  cerr << "wrong number of independent variables in" <<
360  << " I found " << nv1 << " instead of " << nvar << endl;
361  ad_exit(1);
362  }
363  // get last x vector from file
364  std::streamoff sz = parsave.size() * sizeof(double);
365  // backup from end of file
366  uis.seekg(-sz, ios::end);
367  uis >> parsave;
368  cout << " restored " << parsave(parsave.indexmin()) << " "
369  << parsave(parsave.indexmax()) << endl;
370  if (!uis)
371  {
372  cerr << "error resotring last mcmc poistion from file "
373  << ad_comm::adprogram_name + adstring(".psv") << endl;
374  ad_exit(1);
375  }
376  int ii=1;
377  dvector x0(1,nvar);
379  //x0.initialize();
380  //dvector pen_vector(1,nvar);
381  //initial_params::reset(dvar_vector(parsave),pen_vector);
382  //initial_params::xinit(x0);
383  //cout << " x0 " << x0(x0.indexmin()) << endl;
384  }
385  }
386 
387 
388  if ( mcrestart_flag>-1)
389  {
390  pofs_psave= new uostream( (char*)(ad_comm::adprogram_name
391  + adstring(".psv")),ios::app);
392  /*
393  pofs_psave->seekp(0,std::ios::end);
394  *pofs_psave << 123;
395  delete pofs_psave;
396  uistream uis((char*)(ad_comm::adprogram_name + adstring(".psv")));
397  if (!uis)
398  {
399  cerr << "Error trying to open file" <<
400  ad_comm::adprogram_name + adstring(".psv") <<
401  " for mcrestart" << endl;
402  ad_exit(1);
403  } else {
404  int nv1;
405  uis >> nv1;
406  if (nv1 !=nvar) {
407  cerr << "wrong number of independent variables in" <<
408  ad_comm::adprogram_name + adstring(".psv") << endl
409  << " I found " << nv1 << " instead of " << nvar << endl;
410  ad_exit(1);
411  }
412  }
413  */
414  }
415  else
416  {
417  pofs_psave=
418  new uostream((char*)(ad_comm::adprogram_name + adstring(".psv")));
419  }
420 
421  if (!pofs_psave|| !(*pofs_psave))
422  {
423  cerr << "Error trying to open file" <<
425  ad_exit(1);
426  }
427  if (mcrestart_flag == -1 )
428  {
429  (*pofs_psave) << nvar;
430  }
431 
432  // need to rescale the hessian
433  // get the current scale
434  dvector x0(1,nvar);
435  dvector current_scale(1,nvar);
437  // cout << "starting with " << x0(x0.indexmin()) << " "
438  // << x0(x0.indexmax()) << endl;
439  int mctmp=initial_params::mc_phase;
441  initial_params::stddev_scale(current_scale,x);
443  for (int i=1;i<=nvar;i++)
444  {
445  for (int j=1;j<=nvar;j++)
446  {
447  S(i,j)*=old_scale(i)*old_scale(j);
448  }
449  }
450  for (int i=1;i<=nvar;i++)
451  {
452  for (int j=1;j<=nvar;j++)
453  {
454  S(i,j)/=current_scale(i)*current_scale(j);
455  }
456  }
457 
458  //cout << sort(eigenvalues(S)) << endl;
459  dmatrix chd = choleski_decomp(S);
460  //dmatrix tchd = trans(chd);
461  //dmatrix chdinv=inv(chd);
462  //int sgn;
463  // ***************************************************************
464  // ***************************************************************
465  // NEW HYBRID MCMC
466  // ***************************************************************
467  // ***************************************************************
468  independent_variables z(1,nvar);
469  z=x0;
471  dvector g(1,nvar);
472  cout << "initial ll value " << get_hybrid_monte_carlo_value(nvar,z,g)
473  << endl;
474  dvector y(1,nvar);
475  y.initialize();
476 
477  if (mcmc2_flag==0)
478  {
480  }
481 
482  dvector p(1,nvar); // momentum
483  int iseed=2197;
484  if (iseed0) iseed=iseed0;
485  if ( mcrestart_flag>-1)
486  {
487  ifstream ifs("hybrid_seed");
488  ifs >> iseed;
489  if (!ifs)
490  {
491  cerr << "error reading random number seed" << endl;
492  }
493  }
494  random_number_generator rng(iseed);
496  ofstream ogs("sims");
498 
499  z=x0+chd*y;
500  /*double llc=*/get_hybrid_monte_carlo_value(nvar,z,g);
501  /*double llbest=*/get_hybrid_monte_carlo_value(nvar,z,g);
502  //lbmax=llbest;
503 
504 
505  int number_sims;
506  if (nmcmc<=0)
507  {
508  number_sims= 100000;
509  }
510  else
511  {
512  number_sims= nmcmc;
513  }
514  //double hybeps2=0.5*hybeps;
515  double beginprior=get_hybrid_monte_carlo_value(nvar,z,g);
516  dvector Fbegin=g*chd;
517  // use trand(chd) ?
518  //dvector Fbegin=tchd*g;
519  dvector F(1,nvar);
520  F=Fbegin;
521  p.fill_randn(rng);
522  if (robust_hybrid_flag)
523  {
524  double choice=randu(rng);
525  if (choice<0.05)
526  {
527  p*=3.0;
528  }
529  }
530  dmatrix xvalues(1,number_sims,1,nvar);
531  dvector yold(1,nvar);
532  yold=y;
533  double pprob;
534  if (robust_hybrid_flag==0)
535  {
536  pprob=0.5*norm2(p);
537  }
538  else
539  {
540  double r2=0.5*norm2(p);
541  pprob=-log(0.95*exp(-r2)+0.05/3.0*exp(-r2/9.0));
542  }
543  double Hbegin=beginprior+pprob;
544  double tmpprior = 0;
545  int ii=1;
547  // detmine whether to go forward or backward
548  cout << "Starting from y=" << y << endl;
549  cout << "Starting from z=" << z << endl;
550  cout << "Starting from gr=" << g << endl;
551  cout << "Starting from nll=" << beginprior << endl;
552  // cout << "Starting from chd=" << chd << endl;
553 
554  double iaccept=0.0;
555  for (int is=1;is<=number_sims;is++)
556  {
557  int forflag=1;
558  //double rnd=randu(rng);
559  //if (rnd<0.5) forflag=0;
560  double hstep,hstep2;
561  //if (forflag)
562  hstep=hybeps;
563  //else
564  // hstep=-hybeps;
565  hstep2=0.5*hstep;
566  // randomize the number of steps
567  double rnd2=randn(rng);
568 #ifdef OPT_LIB
569  int hnsteps = (int)(exp(0.2 * rnd2) * hybnstep);
570 #else
571  double _hnsteps=exp(0.2 * rnd2) * hybnstep;
572  assert(_hnsteps > 0 && _hnsteps <= (double)INT_MAX);
573  int hnsteps = (int)_hnsteps;
574 #endif
575  for (int i=1;i<=hnsteps;i++)
576  {
577  if (forflag==1)
578  {
579  dvector phalf=p-hstep2*F;
580  if (robust_hybrid_flag==0)
581  {
582  y+=hstep*phalf;
583  }
584  else
585  {
586  //pprob=-log(0.95*exp(-r2)+0.05/3.0*exp(-r2/9.0));
587  double r2=0.5*norm2(phalf);
588  double z=0.95*exp(-r2)+0.05/3.0*exp(-r2/9.0);
589  double xx=(0.95*exp(-r2)+0.05/27.0*exp(-r2/9.0))/z;
590  dvector zz=xx*phalf;
591  y+=hstep*zz;
592  }
593  z=x0+chd*y;
594  tmpprior=get_hybrid_monte_carlo_value(nvar,z,g);
595  F=g*chd;
596  //F=tchd*g;
597  p=phalf-hstep2*F;
598  }
599  else
600  {
601  dvector phalf=p+hstep2*F;
602  if (robust_hybrid_flag==0)
603  {
604  y-=hstep*phalf;
605  }
606  else
607  {
608  double r2=0.5*norm2(phalf);
609  double z=0.95*exp(-r2)+0.05/3.0*exp(-r2/9.0);
610  dvector zz=(0.95*exp(-r2)+0.05/27.0*exp(-r2/9.0))/z*phalf;
611  y-=hstep*phalf;
612  }
613  z=x0+chd*y;
614  tmpprior=get_hybrid_monte_carlo_value(nvar,z,g);
615  F=g*chd;
616  //F=tchd*g;
617  p=phalf+hstep2*F;
618  }
619  }
620  if (robust_hybrid_flag==0)
621  {
622  pprob=0.5*norm2(p);
623  }
624  else
625  {
626  double r2=0.5*norm2(p);
627  pprob=-log(0.95*exp(-r2)+0.05/3.0*exp(-r2/9.0));
628  }
629  double Ham=tmpprior+pprob;
630  double rr=randu(rng);
631  double pp=exp(Hbegin-Ham);
632  p.fill_randn(rng);
633  //p*=1.2;
634  if (robust_hybrid_flag)
635  {
636  double choice=randu(rng);
637  if (choice<0.05)
638  {
639  p*=3.0;
640  }
641  }
642  if (robust_hybrid_flag==0)
643  {
644  pprob=0.5*norm2(p);
645  }
646  else
647  {
648  double r2=0.5*norm2(p);
649  pprob=-log(0.95*exp(-r2)+0.05/3.0*exp(-r2/9.0));
650  }
651  if ((is%50)==1)
652  // cout << iaccept/is << " " << Hbegin-Ham << " " << Ham << endl;
653  cout << " hybrid sim " << is << " accept rate " << iaccept/is
654  << " Hbegin-Ham " << Hbegin-Ham << " Ham " << Ham << endl;
655  if (rr<pp)
656  {
657  iaccept++;
658  yold=y;
659  beginprior=tmpprior;
660  Hbegin=beginprior+pprob;
661  Fbegin=F;
662  ii=1;
664  }
665  else
666  {
667  y=yold;
668  z=x0+chd*y;
669  Hbegin=beginprior+pprob;
670  F=Fbegin;
671  }
672  (*pofs_psave) << parsave;
673  }
674  // cout << " saved " << parsave(parsave.indexmin()) << " "
675  // << parsave(parsave.indexmax()) << endl;
676  //double ll=get_hybrid_monte_carlo_value(nvar,parsave,g);
677  //cout << "ll " << ll << endl;
678  // ***************************************************************
679  // ***************************************************************
680  // ***************************************************************
681  ofstream ofs("hybrid_seed");
682  int seed=(int) (10000*randu(rng));
683  ofs << seed;
684  if (pofs_psave)
685  {
686  delete pofs_psave;
687  pofs_psave=NULL;
688  }
689  }
690 }
691 
void ad_update_mcmchist_report(dmatrix &mcmc_values, ivector &number_offsets, dvector &mean_mcmc_values, dvector &h, int ff=0)
static void restore_all_values(const dvector &x, const int &ii)
Definition: model3.cpp:105
void make_preliminary_hist(const dvector &s, const dvector &m, int nsim, const dmatrix &values, dmatrix &hist, const dvector &h, int slots, double total_spread, int probflag=0)
Definition: xxmcmc.cpp:1241
void read_hessian_matrix_and_scale(int nvar, const dmatrix &S, const dvector &pen_vector)
Definition: xxmcmc.cpp:1362
static void set_NO_DERIVATIVES(void)
Disable accumulation of derivative information.
Definition: gradstrc.cpp:641
Description not yet available.
static void set_active_random_effects(void)
Definition: model.cpp:267
void allocate(void)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: dmat0.cpp:8
dvector read_old_scale(int &old_nvar)
Definition: xxmcmc.cpp:1667
dvector new_probing_bounded_multivariate_normal(int nvar, const dvector &a1, const dvector &b1, dmatrix &ch, const double &wght, double pprobe, random_number_generator &rng)
Definition: prmonte.cpp:14
void ad_update_mcmc_stats_report(int feval, int iter, double fval, double gmax)
int read_hist_data(const dmatrix &hist, const dvector &h, dvector &m, const dvector &s, const dvector &parsave, long int &iseed, const double &size_scale)
#define x
Vector of double precision numbers.
Definition: dvector.h:50
static char ** argv
Definition: fvar.hpp:8866
void allocate(int ncl, int ncu)
Allocate memory for a dvector.
Definition: dvector.cpp:409
void print_hist_data(const dmatrix &hist, const dmatrix &values, const dvector &h, dvector &m, const dvector &s, const dvector &parsave, long int iseed, double size_scale)
void hybrid_mcmc_routine(int, int, double, int)
Description not yet available.
Definition: hybmcmc.cpp:83
exitptr ad_exit
Definition: gradstrc.cpp:53
void set_labels_for_mcmc(void)
static dvariable reset(const dvar_vector &x)
Definition: model.cpp:345
void ad_update_function_minimizer_report(int feval, int iter, int phase, double fval, double gmax, const char *cbuf)
ADMB variable vector.
Definition: fvar.hpp:2172
static void restore_start_phase(void)
Definition: model.cpp:275
static adstring adprogram_name
Definition: fvar.hpp:8860
static int num_stddev_calc(void)
Definition: model2.cpp:51
int atoi(adstring &s)
Returns a integer converted from input s.
Definition: atoi.cpp:20
df1_one_matrix choleski_decomp(const df1_one_matrix &MM)
Definition: df11fun.cpp:606
static int nvarcalc()
Definition: model.cpp:152
double randu(const random_number_generator &rng)
Uniform random number generator.
Definition: rngen.cpp:198
void fill_randn(long int &n)
Fill vector with random numbers.
Definition: ranfill.cpp:231
Description not yet available.
Definition: fvar.hpp:7951
void ad_update_mcmc_report(dmatrix &m, int i, int j, int ff=0)
void write_empirical_covariance_matrix(int ncor, const dvector &s_mean, const dmatrix &s_covar, adstring &prog_name)
Writes the covariance matrix out to a file, which is prog_name prepended to the extension &#39;...
Definition: xxmcmc.cpp:917
prnstream & endl(prnstream &)
Description not yet available.
Definition: fvar.hpp:1937
Array of integers(int) with indexes from index_min to indexmax.
Definition: ivector.h:50
double get_hybrid_monte_carlo_value(int nvar, const independent_variables &y, dvector &g)
Written by Dave, commented by Cole starting 8/31/2016 Description not yet available.
Description not yet available.
Definition: fvar.hpp:3398
static int argc
Definition: fvar.hpp:8863
Description not yet available.
static void xinit(const dvector &x)
Definition: model.cpp:226
d3_array exp(const d3_array &arr3)
Returns d3_array results with computed exp from elements in arr3.
Definition: d3arr2a.cpp:28
static void set_inactive_random_effects(void)
Definition: model.cpp:288
Description not yet available.
Definition: fvar.hpp:2819
void read_covariance_matrix(const dmatrix &S, int nvar, int &hbf, dvector &sscale)
Definition: xxmcmc.cpp:1320
void initialize(void)
Initialze all elements of dvector to zero.
Definition: dvect5.cpp:10
double norm2(const d3_array &a)
Return sum of squared elements in a.
Definition: d3arr2a.cpp:167
values
Definition: adjson.h:22
static void copy_all_values(const dvector &x, const int &ii)
Definition: model3.cpp:9
void read_empirical_covariance_matrix(int nvar, const dmatrix &S, const adstring &prog_name)
Reads the covariance matrix from a file, which is the program name prepended to the extension &#39;...
Definition: xxmcmc.cpp:961
unsigned int size() const
Get number of elements in array.
Definition: dvector.h:209
int option_match(int argc, char *argv[], const char *string)
Checks if the program has been invoked with a particular command line argument (&quot;string&quot;).
Definition: optmatch.cpp:25
void store_mcmc_values(const ofstream &ofs)
double randn(const random_number_generator &rng)
Normal number generator.
Definition: rngen.cpp:183
int maxnz(const dvector &xa)
Definition: xxmcmc.cpp:1161
Description not yet available.
void new_probing_bounded_multivariate_normal_mcmc(int nvar, const dvector &a1, const dvector &b1, dmatrix &ch, const double &wght, const dvector &_y, double pprobe, random_number_generator &rng)
Definition: prmonte.cpp:167
int minnz(const dvector &x)
Definition: xxmcmc.cpp:1144
static void set_YES_DERIVATIVES(void)
Enable accumulation of derivative information.
Definition: gradstrc.cpp:650
double better_rand(long int &idum)
Description not yet available.
Definition: bet_rand.cpp:18
Description not yet available.
Definition: fvar.hpp:3516
void read_hessian_matrix_and_scale1(int nvar, const dmatrix &_SS, double s, int mcmc2_flag)
Definition: xxmcmc.cpp:1453
static int stddev_scale(const dvector &d, const dvector &x)
Definition: model.cpp:202
static int mc_phase
Definition: admodel.h:846
void initialize(void)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: dmat7.cpp:12
void add_hist_values(const dvector &s, const dvector &m, const dmatrix &hist, dvector &mcmc_values, double llc, const dvector &h, int nslots, double total_spreadd, int probflag=0)
Definition: xxmcmc.cpp:1178
d3_array log(const d3_array &arr3)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: d3arr2a.cpp:13
static int num_stddev_params
Definition: admodel.h:2219