ADMB Documentation  -a65f1c97
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
df1b2lap.cpp
Go to the documentation of this file.
1 /*
2  * $Id$
3  *
4  * Author: David Fournier
5  * Copyright (c) 2008-2012 Regents of the University of California
6  */
11 #include <sstream>
12 using std::istringstream;
13 #include <admodel.h>
14 #include <df1b2fun.h>
15 #include <adrndeff.h>
16 #ifndef OPT_LIB
17  #include <cassert>
18  #include <climits>
19 #endif
20 
21 static int no_stuff=0;
22 
23 #ifdef DIAG
24 int fcount =0;
25 static int write_sparse_flag=0;
26  static void trapper(void)
27  {
28  //int x=5;
29  }
30 #endif
32 unsigned int global_nvar=0;
33 
35 
36 double evaluate_function(const dvector& x,function_minimizer * pfmin);
37 void get_newton_raphson_info(int xs,int us,const init_df1b2vector _y,
38  dmatrix& Hess, dvector& grad,
40 
41 //dvariable AD_uf_inner(const dvector& x,const dvar_vector& u);
42 void get_second_ders(int xs,int us,const init_df1b2vector y,dmatrix& Hess,
45 
47 
53 
57 
63  (const dvector& x,function_minimizer * pfmin)
64 {
65  std::ostream& output_stream = get_output_stream();
66  //int on,nopt;
67  pfmin->inner_opt_flag=1;
68  double f=0.0;
69  double fb=1.e+100;
70  dvector g(1,usize);
71  dvector ub(1,usize);
72  independent_variables u(1,usize);
73  gradcalc(0,g);
74  fmc1.itn=0;
75  fmc1.ifn=0;
76  fmc1.ireturn=0;
77  initial_params::xinit(u); // get the initial values into the
78  initial_params::xinit(ubest); // get the initial values into the
79  fmc1.ialph=0;
80  fmc1.ihang=0;
81  fmc1.ihflag=0;
82  fmc1.use_control_c=0;
83 
84  if (init_switch)
85  {
86  u.initialize();
87  }
88  else
89  {
90  u=ubest;
91  }
92 
93  fmc1.dfn=1.e-2;
94  dvariable pen=0.0;
95  //cout << "starting norm(u) = " << norm(u) << endl;
96  while (fmc1.ireturn>=0)
97  {
98  /*
99  double nu=norm(value(u));
100  if (nu>400)
101  {
102  cout << "U norm(u) = " << nu << endl;
103  }
104  cout << "V norm(u) = " << nu
105  << " f = " << f << endl;
106  */
107  fmc1.fmin(f,u,g);
108  //cout << "W norm(u) = " << norm(value(u)) << endl;
109  if (fmc1.ireturn>0)
110  {
111  dvariable vf=0.0;
114  pfmin->AD_uf_inner();
115  if (saddlepointflag)
116  {
118  }
120  {
122  }
124 
125  /* this is now done in the operator = function
126  if (quadratic_prior::get_num_quadratic_prior()>0)
127  {
128  vf+= quadratic_prior::get_quadratic_priors();
129  }
130  */
131 
133 
134  //cout << " pen = " << pen << endl;
135  if (noboundepen_flag==0)
136  {
137  vf+=pen;
138  }
139  f=value(vf);
140  if (f<fb)
141  {
142  fb=f;
143  ub=u;
144  }
145  gradcalc(usize,g);
146  //cout << " f = " << setprecision(17) << f << " " << norm(g)
147  // << " " << norm(u) << endl;
148  }
149  u=ub;
150  }
151  output_stream << " inner maxg = " << std::scientific << setprecision(10) << fmc1.gmax;
152 #ifdef DIAG
153  if (fabs(fmc1.gmax)>1.e+3)
154  trapper();
155 #endif
156 
157  if (fabs(fmc1.gmax)>1.e-4)
158  {
159  fmc1.itn=0;
160  //fmc1.crit=1.e-9;
161  fmc1.ifn=0;
162  fmc1.ireturn=0;
163  fmc1.ihang=0;
164  fmc1.ihflag=0;
165  fmc1.ialph=0;
166  initial_params::xinit(u); // get the initial values into the
167  //u.initialize();
168  while (fmc1.ireturn>=0)
169  {
170  fmc1.fmin(f,u,g);
171  if (fmc1.ireturn>0)
172  {
173  dvariable vf=0.0;
176  pfmin->AD_uf_inner();
178  {
180  }
183 
184  vf+=pen;
185  f=value(vf);
186  if (f<fb)
187  {
188  fb=f;
189  ub=u;
190  }
191  gradcalc(usize,g);
192  //cout << " f = " << setprecision(15) << f << " " << norm(g) << endl;
193  }
194  }
195  u=ub;
196  output_stream << " Inner second time = " << std::scientific << setprecision(10) << fmc1.gmax;
197  }
198  output_stream << " Inner f = " << std::scientific << setprecision(10) << fb << endl;
199  fmc1.ireturn=0;
200  fmc1.fbest=fb;
203  pfmin->AD_uf_inner();
205  {
207  }
209  pfmin->inner_opt_flag=0;
210  return u;
211 }
212 
218  (const dvector& x,function_minimizer * pfmin)
219 {
220  double f=0.0;
221  //dvector g(1,usize);
222  //g.initialize();
223  independent_variables u(1,usize);
224  if (init_switch)
225  {
226  initial_params::xinit(u); // get the initial values into the
227  }
228  else
229  {
230  u=ubest;
231  }
232 
233  int iprint_save=pfmin->iprint;
234  pfmin->iprint=inner_iprint;
235  pfmin->limited_memory_quasi_newton(f,u,5,1,inner_maxfn,inner_crit);
236  pfmin->iprint=iprint_save;
238  {
239  cout << " Inner f = " << std::scientific << setprecision(10) << f << endl;
240  }
241  //gradient_structure::set_NO_DERIVATIVES();
242  //initial_params::set_inactive_only_random_effects();
243  //initial_params::reset(dvar_vector(x));
244  return u;
245 }
246 
251 void set_partition_sizes(int & num_der_blocks,ivector& minder,
252  ivector& maxder,int nvariables)
253 {
254  if (num_der_blocks==0)
255  {
256  num_der_blocks=1;
257  }
258 #if defined(USE_ADPVM)
260  {
261  switch (ad_comm::pvm_manager->mode)
262  {
263  case 1: //master
264  {
265  int i;
266  int ndb=ad_comm::pvm_manager->num_slave_processes+1;
267 
268  minder.allocate(1,ndb);
269  maxder.allocate(1,ndb);
270 
271  int nd=nvariables/ndb;
272  int r= nvariables - nd * ndb;
273  ivector partition(1,ndb);
274  partition=nd;
275  partition(1,r)+=1;
276  minder(1)=1;
277  maxder(1)=partition(1);
278 
279  for (i=2;i<=ndb;i++)
280  {
281  minder(i)=maxder(i-1)+1;
282  maxder(i)=minder(i)+partition(i)-1;
283  }
284  send_int_to_slaves(minder(2,ndb));
285  send_int_to_slaves(maxder(2,ndb));
286  break;
287  }
288 
289  case 2: //slave
290  minder.allocate(1,1);
291  maxder.allocate(1,1);
292  minder(1)=get_int_from_master();
293  maxder(1)=get_int_from_master();
294  break;
295  }
296  }
297  else
298 #endif // #USE_ADPVM
299  {
300  minder.allocate(1,num_der_blocks);
301  maxder.allocate(1,num_der_blocks);
302 
303  int nd=nvariables/num_der_blocks;
304  int r= nvariables - nd * num_der_blocks;
305  ivector partition(1,num_der_blocks);
306  partition=nd;
307  partition(1,r)+=1;
308 
309  minder(1)=1;
310  maxder(1)=partition(1);
311  for (int i=2;i<=num_der_blocks;i++)
312  {
313  minder(i)=maxder(i-1)+1;
314  maxder(i)=minder(i)+partition(i)-1;
315  }
316  }
317 }
318 
324 (
325  int _xsize,
326  int _usize,
327  int _minder,
328  int _maxder,
329  function_minimizer* _pmin
330 ):
331  init_switch(1),
332  separable_call_level(0),
333  triplet_information(0),
334  compressed_triplet_information(0),
335  nested_separable_calls_counter(1,10),
336  nested_tree_position(1,5),
337  local_dtemp(1,_xsize),
338  nr_debug(0),
339  pmin(_pmin),
340  block_diagonal_flag(0),
341  xsize(_xsize),
342  usize(_usize),
343  bHess_pd_flag(0),
344  sparse_triplet(0),
345  sparse_iterator(0),
346  sparse_count(0),
347  sparse_count_adjoint(0),
348  sparse_triplet2(0),
349  vsparse_triplet(0),
350  vsparse_triplet_adjoint(0),
351  sparse_symbolic(0),
352  sparse_symbolic2(0),
353  fmc1(usize,1),
354  fmc(_xsize),
355  xadjoint(1,_xsize),
356  check_local_uadjoint(1,_usize),
357  check_local_uadjoint2(1,_usize),
358  check_local_xadjoint(1,_xsize),
359  check_local_xadjoint2(1,_xsize),
360  uadjoint(1,_usize),
361  uhat(1,_usize),
362  ubest(1,_usize)
363 {
364  ubest.initialize();
365  nested_shape.allocate(100,100,10);
366  nested_shape.initialize();
367  nested_tree_position.initialize();
368  nested_separable_calls_counter.initialize();
369  calling_set=0;
370  antiepsilon=0;
371  dd_nr_flag=0;
372  no_re_ders_flag =0;
373  importance_sampling_components=0;
374  is_diagnostics_flag=0;
375  importance_sampling_values = 0;
376  importance_sampling_weights = 0;
377  no_function_component_flag=0;
378  uhat.initialize();
379  hesstype=1;
380  use_gauss_hermite=0;
381  in_gauss_hermite_phase=0;
382  multi_random_effects=0;
383  //var_flag=1;
384  var_flag=0;
385  num_separable_calls=0;
386  separable_calls_counter=0;
387  importance_sampling_counter=0;
388  num_local_re_array=0;
389  num_local_fixed_array=0;
390  isfunnel_flag=0;
391  antiflag=0;
392  rseed=3457;
393  nfunnelblocks=0;
394  separable_function_difference=0;
395  gh=0;
396  block_diagonal_hessian=0;
397  block_diagonal_Dux=0;
398  block_diagonal_re_list=0;
399  block_diagonal_fe_list=0;
400  block_diagonal_vhessian=0;
401  block_diagonal_vhessianadjoint=0;
402  block_diagonal_ch=0;
403  block_diagonal_vch=0;
404  have_users_hesstype=0;
405  pHess_non_quadprior_part=0;
406  bw=0;
407  bHess=0;
408  grad_x_u=0;
409  grad_x=0;
410  int ndi=20000;
411  int nopt=0;
412  int on=-1;
413  if ( (inner_lmnflag=option_match(ad_comm::argc,ad_comm::argv,"-ndi",nopt))
414  >-1)
415  {
416  if (!nopt)
417  {
418  cerr << "Usage -ndi option needs integer -- set to default 20000.\n";
419  }
420  else
421  {
422  int jj=atoi(ad_comm::argv[inner_lmnflag+1]);
423  if (jj<=0)
424  {
425  cerr << "Usage -ndi option needs positive integer"
426  " -- set to defalt 20000" << endl;
427  }
428  else
429  {
430  ndi=jj;
431  }
432  }
433  }
434  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-noreders",nopt))>-1)
435  {
436  no_re_ders_flag=1;
437  }
438  derindex=new imatrix(1,ndi);
439  bHessadjoint=0;
441  {
443  {
445  }
447  }
448  inner_lmnflag=0;
449  inner_lmnsteps=10;
450  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-noinit",nopt))>-1)
451  {
452  init_switch=0;
453  }
454  if ( (inner_lmnflag=option_match(ad_comm::argc,ad_comm::argv,"-ilmn",nopt))
455  >-1)
456  {
457  if (!nopt)
458  {
459  cerr << "Usage -ilmn option needs integer -- set to default 10" << endl;
460  }
461  else
462  {
463  int jj=atoi(ad_comm::argv[inner_lmnflag+1]);
464  if (jj<=0)
465  {
466  cerr << "Usage -ilmn option needs positive integer"
467  " -- set to defalt 10" << endl;
468  }
469  else
470  {
471  inner_lmnsteps=jj;
472  }
473  }
474  }
475  else
476  {
477  inner_lmnflag=0;
478  }
479  inner_noprintx=0;
480 
481  num_der_blocks=1; // default value
482  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-ndb",nopt))>-1)
483  {
484  if (!nopt)
485  {
486  cerr << "Usage -ndb option needs non-negative integer -- ignored.\n";
487  }
488  else
489  {
490  int _num_der_blocks=atoi(ad_comm::argv[on+1]);
491  if (_num_der_blocks<=0)
492  {
493  cerr << "Usage -ndb option needs positive integer -- ignored" << endl;
494  }
495  else
496  {
497  num_der_blocks=_num_der_blocks;
498  }
499  }
500  }
501 
502  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-isf",nopt))>-1)
503  {
504  if (!nopt)
505  {
506  cerr << "Usage -isf option needs non-negative integer -- ignored.\n";
507  }
508  else
509  {
510  int _nfunnelblocks=atoi(ad_comm::argv[on+1]);
511  if (_nfunnelblocks<=0)
512  {
513  cerr << "Usage -isf option needs positive integer -- ignored" << endl;
514  }
515  else
516  {
517  nfunnelblocks=_nfunnelblocks;
518  isfunnel_flag=1;
519  }
520  }
521  }
522 
523  antiflag=0;
524  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-anti",nopt))>-1)
525  {
526  antiflag=1;
527  }
528 
529  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-nrdbg",nopt))>-1)
530  {
531  nr_debug=1;
532  }
533 
534  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-ddnr",nopt))>-1)
535  {
536  dd_nr_flag=1;
537  }
538 
539  nvariables=xsize+usize;
540  /*
541  int rem=0;
542  if (nvariables%xsize!=0)
543  rem=1;
544  */
545 
546  if (nvariables/num_der_blocks<xsize)
547  {
548  cerr << " Error -- the number of der blocks (-ndb) can not be larger than "
549  << nvariables/xsize << endl;
550  ad_exit(1);
551  }
552 
553  set_partition_sizes(num_der_blocks,minder,maxder,nvariables);
554 
555  // !!! remove
556  //maxder(1)=5;
557 
558  fmc1.crit=1.e-3;
559 
560  //cout << "Need to fix Hess" << endl;
561  for (int i=1;i<=num_der_blocks;i++)
562  {
563  if (minder(i)<1 || maxder(i) > nvariables || maxder(i) < minder(i))
564  {
565  cerr << " minder or maxder value out of bounds in"
566  "laplace_approximation_calculator constructor "
567  << endl << " values are minder = " << minder
568  << " maxder = " << maxder << endl;
569  ad_exit(1);
570  }
571  }
572 
573  num_nr_iters=2;
574  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-nr",nopt))>-1)
575  {
576  if (!nopt)
577  {
578  cerr << "Usage -nr option needs non-negative integer -- ignored" << endl;
579  }
580  else
581  {
582  int _num_nr_iters=atoi(ad_comm::argv[on+1]);
583  if (_num_nr_iters<0)
584  {
585  cerr << "Usage -nr option needs non-negative integer -- ignored.\n";
586  }
587  else
588  {
589  num_nr_iters=_num_nr_iters;
590  }
591  }
592  }
593 
594  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-nochol",nopt))>-1)
595  {
597  }
598 
600  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-phe",nopt))>-1)
601  {
603  }
604 
605  double _nr_crit=1.e-11;
606  nr_crit=1.e-11;
607  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-nrcrit",nopt))>-1)
608  {
609  if (!nopt)
610  {
611  cerr << "Usage -nrcrit option needs number -- ignored" << endl;
612  }
613  else
614  {
615  istringstream ist(ad_comm::argv[on+1]);
616  ist >> _nr_crit;
617 
618  if (_nr_crit<=0)
619  {
620  cerr << "Usage -nrcrit option needs positive number -- ignored.\n";
621  _nr_crit=0.0;
622  }
623  }
624  if (_nr_crit>0)
625  {
626  nr_crit=_nr_crit;
627  }
628  }
629  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-gh",nopt))>-1)
630  {
631  if (!nopt)
632  {
633  cerr << "Usage -gh option needs positive integer -- ignored" << endl;
634  }
635  else
636  {
637  int _inner_gh=atoi(ad_comm::argv[on+1]);
638  if (_inner_gh<=0)
639  {
640  cerr << "Usage -gh option needs positive integer -- ignored" << endl;
641  }
642  else
643  {
644  use_gauss_hermite=_inner_gh;
645  }
646  }
647  }
648 
649  inner_crit=1.e-3;
650  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-icrit",nopt))>-1)
651  {
652  double _inner_crit=0.0;;
653  if (!nopt)
654  {
655  cerr << "Usage -icrit option needs number -- ignored" << endl;
656  }
657  else
658  {
659  istringstream ist(ad_comm::argv[on+1]);
660  ist >> _inner_crit;
661 
662  if (_inner_crit<=0)
663  {
664  cerr << "Usage -icrit option needs positive number -- ignored" << endl;
665  _inner_crit=0.0;
666  }
667  }
668  if (_inner_crit>0)
669  {
670  inner_crit=_inner_crit;
671  }
672  }
673  fmc1.crit=inner_crit;
674 
675  fmc.dfn=.01;
676 
677  inner_iprint=0;
678  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-iiprint",nopt))>-1)
679  {
680  if (!nopt)
681  {
682  cerr << "Usage -iprint option needs non-negative integer -- ignored.\n";
683  }
684  else
685  {
686  int _inner_iprint=atoi(ad_comm::argv[on+1]);
687  if (_inner_iprint<=0)
688  {
689  cerr << "Usage -iip option needs non-negative integer -- ignored.\n";
690  }
691  else
692  {
693  inner_iprint=_inner_iprint;
694  }
695  }
696  }
697  fmc1.iprint=inner_iprint;
698 
699  inner_maxfn=1000;
700  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-imaxfn",nopt))>-1)
701  {
702  if (!nopt)
703  {
704  cerr << "Usage -maxfn option needs non-negative integer -- ignored.\n";
705  }
706  else
707  {
708  int _inner_maxfn=atoi(ad_comm::argv[on+1]);
709  if (_inner_maxfn<0)
710  {
711  cerr << "Usage -iip option needs non-negative integer -- ignored.\n";
712  }
713  else
714  {
715  inner_maxfn=_inner_maxfn;
716  }
717  }
718  }
719  fmc1.maxfn=inner_maxfn;
720  // what sort of structure on the Hessian do we have
721  nopt=0;
722 
723  rseed=3456;
724  num_importance_samples=0;
725  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-is",nopt))>-1)
726  {
727  if (!nopt)
728  {
729  cerr << "Usage -is option needs positive integer -- ignored" << endl;
730  }
731  else
732  {
733  int tht=atoi(ad_comm::argv[on+1]);
734  if (tht<=0)
735  {
736  cerr << "Usage -is option needs non-negative integer -- ignored.\n";
737  }
738  else
739  {
740  num_importance_samples=tht;
741  }
742  if (nopt==2)
743  {
744  int rseed1=atoi(ad_comm::argv[on+2]);
745  if (rseed1<=0)
746  {
747  cerr << "Usage -is option needs non-negative integer -- ignored.\n";
748  }
749  else
750  {
751  rseed=rseed1;
752  }
753  }
754  }
755  }
756  int use_balanced=0;
757  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-isb",nopt))>-1)
758  {
759  use_balanced=1;
760  if (!nopt)
761  {
762  cerr << "Usage -isb option needs positive integer -- ignored" << endl;
763  }
764  else
765  {
766  int tht=atoi(ad_comm::argv[on+1]);
767  if (tht<=0)
768  {
769  cerr << "Usage -isb option needs non-negative integer -- ignored.\n";
770  }
771  else
772  {
773  num_importance_samples=2*tht;
774  }
775  if (nopt==2)
776  {
777  int rseed1=atoi(ad_comm::argv[on+2]);
778  if (rseed1<=0)
779  {
780  cerr << "Usage -isb option needs non-negative integer -- ignored.\n";
781  }
782  else
783  {
784  rseed=rseed1;
785  }
786  }
787  }
788  }
789  use_outliers=0;
790  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-iso",nopt))>-1)
791  {
792  use_outliers=1;
793  }
794  if (num_importance_samples)
795  {
796  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-isdiag",nopt))>-1)
797  {
798  is_diagnostics_flag=1;
799  }
800  if (importance_sampling_values)
801  {
802  delete importance_sampling_values;
803  importance_sampling_values=0;
804  }
805  importance_sampling_values =
806  new dvector(1,num_importance_samples);
807 
808  if (importance_sampling_weights)
809  {
810  delete importance_sampling_weights;
811  importance_sampling_weights=0;
812  }
813  importance_sampling_weights =
814  new dvector(1,num_importance_samples);
815 
816  random_number_generator rng(rseed);
817  if (allocated(epsilon)) epsilon.deallocate();
818  if (use_balanced)
819  {
820  // check for eveb num_importance samples
821  if (num_importance_samples%2)
822  num_importance_samples+=1;
823  }
824  epsilon.allocate(1,num_importance_samples,1,usize);
825  if (use_balanced)
826  {
827  int n2=num_importance_samples/2;
828  epsilon.sub(1,n2).fill_randn(rng);
829  if (use_outliers)
830  {
831  dmatrix os(1,num_importance_samples,1,usize);
832  os.fill_randu(rng);
833  for (int i=1;i<=num_importance_samples;i++)
834  {
835  for (int j=1;j<=usize;j++)
836  {
837  if (os(i,j)<0.05) epsilon(i,j)*=3.0;
838  }
839  }
840  }
841  for (int i=1;i<=n2;i++)
842  {
843  epsilon(i+n2)=-epsilon(i);
844  }
845  }
846  else
847  {
848  epsilon.fill_randn(rng);
849  if (use_outliers)
850  {
851  dmatrix os(1,num_importance_samples,1,usize);
852  os.fill_randu(rng);
853  for (int i=1;i<=num_importance_samples;i++)
854  {
855  for (int j=1;j<=usize;j++)
856  {
857  if (os(i,j)<0.05) epsilon(i,j)*=3.0;
858  }
859  }
860  }
861  }
862 
863  double eps_mult=1.0;
864  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-epsmult",nopt))>-1)
865  {
866  if (!nopt)
867  {
868  cerr << "Usage -epsmult option needs number -- ignored" << endl;
869  }
870  else
871  {
872  istringstream ist(ad_comm::argv[on+1]);
873  ist >> eps_mult;
874 
875  if (eps_mult<=0.0 || eps_mult>1.0)
876  {
877  cerr << "Usage -epsmult option needs positive number between 0 and 1 "
878  "-- ignored" << endl;
879  eps_mult=1.0;
880  }
881  }
882  for (int i=1;i<=num_importance_samples;i++)
883  {
884  epsilon(i)*=eps_mult;
885  }
886  }
887  }
888 
889  nopt=0;
890  if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-ht",nopt))>-1)
891  {
892  if (!nopt)
893  {
894  cerr << "Usage -ht option needs positive integer -- ignored" << endl;
895  set_default_hessian_type();
896  }
897  else
898  {
899  int tht=atoi(ad_comm::argv[on+1]);
900  if (tht<=0)
901  {
902  cerr << "Usage -ht option needs non-negative integer -- ignored.\n";
903  set_default_hessian_type();
904  }
905  else
906  {
907  have_users_hesstype=1;
908  hesstype=tht;
909  }
910  if (nopt>=2)
911  {
912  int tbw=atoi(ad_comm::argv[on+2]);
913  if (tbw<=0)
914  {
915  cerr << "Usage -ht option needs non-negative bandwidth"
916  " -- ignored" << endl;
917  }
918  else
919  {
920  ad_comm::bandwidth=tbw;
921  }
922  }
923  }
924  }
925  else
926  {
927  set_default_hessian_type();
928  }
929 
930  // !! need to check nvar calculation
931 #ifndef OPT_LIB
932  assert(maxder(1) >= minder(1));
933 #endif
934  nvar = (unsigned int)(maxder(1) - minder(1) + 1);
935 
936  switch (hesstype)
937  {
938  case 0:
939  case 1:
940  case 4:
941  grad.allocate(1,usize);
942  Hess.allocate(1,usize,1,usize);
943  Hessadjoint.allocate(1,usize,1,usize);
944  Dux.allocate(1,usize,1,xsize);
945  break;
946  case 3:
947  {
948  int bw=2;
950  bHess=new banded_symmetric_dmatrix(1,usize,bw);
951  bHessadjoint=new banded_symmetric_dmatrix(1,usize,bw);
952  grad.allocate(1,usize);
953  Dux.allocate(1,usize,1,xsize);
954  }
955  break;
956  default:
957  break;
958  }
959 
960  step.allocate(1,usize);
961  // !!! nov 12
962  f1b2gradlist = new df1b2_gradlist(4000000U,200000U,8000000U,400000U,
963  2000000U,100000U,adstring("f1b2list1"));
964 
965  if (hesstype==2)
966  {
967  ad_dstar::allocate(nvar);
968  global_nvar=nvar;
970  df1b2variable::set_minder(minder(1));
971  df1b2variable::set_maxder(maxder(1));
973  }
974  if (hesstype!=2)
975  {
976  if (sparse_hessian_flag==0)
977  {
978  ad_dstar::allocate(nvar);
979  global_nvar=nvar;
981  df1b2variable::set_minder(minder(1));
982  df1b2variable::set_maxder(maxder(1));
984  y.allocate(1,nvariables);
985  }
986  else
987  {
988  unsigned int nsave=nvar;
989  nvar=1;
990  ad_dstar::allocate(nvar);
991  global_nvar=nvar;
993  df1b2variable::set_minder(minder(1));
994  df1b2variable::set_maxder(maxder(1));
996  y.allocate(1,nvariables);
997  nvar=nsave;
998  global_nvar=nvar;
999  }
1000  }
1002  {
1003  cout << "this can't happen" << endl;
1004  ad_exit(1);
1005  }
1008 #ifndef OPT_LIB
1009  assert(nvariables >= 0);
1010 #endif
1012  (unsigned int)nvariables;
1013  //df1b2variable::adpool_counter++;
1015 }
1016 
1022  int _xsize,
1023  int _usize,
1024  ivector _minder,
1025  ivector _maxder,
1026  function_minimizer* _pmin
1027 ):
1028  init_switch(1),
1029  separable_call_level(1),
1030  triplet_information(0),
1031  compressed_triplet_information(0),
1032  nested_separable_calls_counter(1,10),
1033  nested_tree_position(1,5),
1034  nr_debug(0),
1035  pmin(_pmin),
1036  xsize(_xsize),
1037  usize(_usize),
1038  bHess_pd_flag(0),
1039  sparse_triplet(0),
1040  sparse_iterator(0),
1041  sparse_count(0),
1042  sparse_count_adjoint(0),
1043  sparse_triplet2(0),
1044  vsparse_triplet(0),
1045  vsparse_triplet_adjoint(0),
1046  sparse_symbolic(0),
1047  sparse_symbolic2(0),
1048  fmc1(0),
1049  fmc(_xsize),
1050  xadjoint(1,_xsize),
1051  check_local_uadjoint(1,_usize),
1052  check_local_uadjoint2(1,_usize),
1053  check_local_xadjoint(1,_xsize),
1054  check_local_xadjoint2(1,_xsize),
1055  uadjoint(1,_usize),
1056  uhat(1,_usize)
1057 {
1060  //hesstype=1;
1061  uhat.initialize();
1062  //var_flag=1;
1063  var_flag=0;
1064  calling_set=0;
1065  antiepsilon=0;
1069  isfunnel_flag=0;
1070  antiflag=0;
1071  rseed=3457;
1072  nfunnelblocks=0;
1080  bw=0;
1081  bHess=0;
1082  grad_x_u=0;
1083  grad_x=0;
1085  int mmin=_minder.indexmin();
1086  int mmax=_minder.indexmax();
1087  num_der_blocks= mmax-mmin+1;
1088  minder.allocate(mmin,mmax);
1089  maxder.allocate(mmin,mmax);
1090  minder=_minder;
1091  maxder=_maxder;
1092  inner_iprint = 0;
1094  fmc1.crit=1.e-3;
1096  Hess.allocate(1,usize,1,usize);
1097  for (int i=1;i<=num_der_blocks;i++)
1098  {
1099  if (minder(i)<1 || maxder(i) > nvariables || maxder(i) < minder(i))
1100  {
1101  cerr << " minder or maxder value out of bounds in"
1102  "laplace_approximation_calculator constructor "
1103  << endl << " values are minder = " << minder
1104  << " maxder = " << maxder << endl;
1105  ad_exit(1);
1106  }
1107  }
1109  grad.allocate(1,usize);
1110 
1111 #ifndef OPT_LIB
1112  assert(maxder(1) >= minder(1));
1113 #endif
1114  nvar = (unsigned int)(maxder(1) - minder(1) + 1);
1115  Hessadjoint.allocate(1,usize,1,usize);
1116  Dux.allocate(1,usize,1,xsize);
1117  // !!! nov 12
1118  f1b2gradlist = new df1b2_gradlist(4000000U,200000U,8000000U,400000U,
1119  2000000U,100000U,adstring("f1b2list1"));
1121  global_nvar=nvar;
1126  y.allocate(1,nvariables);
1127 
1128  gh = nullptr;
1129  importance_sampling_weights = nullptr;
1131  importance_sampling_values = nullptr;
1133  Hess_components = nullptr;
1134  derindex = nullptr;
1135  bHessadjoint = nullptr;
1136  block_diagonal_vch = nullptr;
1137  block_diagonal_ch = nullptr;
1138 }
1141 {
1143  {
1146  }
1148  {
1151  }
1153  {
1156  }
1157  if (block_diagonal_vch)
1158  {
1159  delete block_diagonal_vch;
1161  }
1162  if (block_diagonal_ch)
1163  {
1164  delete block_diagonal_ch;
1166  }
1168  {
1169  delete block_diagonal_hessian;
1171  }
1172  if (block_diagonal_Dux )
1173  {
1174  delete block_diagonal_Dux;
1175  block_diagonal_Dux =0;
1176  }
1178  {
1179  delete block_diagonal_re_list;
1181  }
1183  {
1184  delete block_diagonal_fe_list;
1186  }
1188  {
1189  delete block_diagonal_vhessian;
1191  }
1193  {
1196  }
1198  {
1201  }
1202 
1203  if (derindex)
1204  {
1205  delete derindex;
1206  derindex=0;
1207  }
1208 
1210  {
1211  delete pHess_non_quadprior_part;
1213  }
1214 
1215  if (triplet_information)
1216  {
1217  delete triplet_information;
1219  }
1220 
1221  if (bHessadjoint)
1222  {
1223  delete bHessadjoint;
1224  bHessadjoint=0;
1225  }
1226  if (grad_x)
1227  {
1228  delete grad_x;
1229  grad_x=0;
1230  }
1231  if (grad_x_u)
1232  {
1233  delete grad_x_u;
1234  grad_x_u=0;
1235  }
1236  if (bHess)
1237  {
1238  delete bHess;
1239  bHess=0;
1240  }
1241  if (f1b2gradlist)
1242  {
1244  delete f1b2gradlist;
1245  f1b2gradlist=0;
1246  }
1247  if (df1b2variable::pool)
1248  {
1249  //df1b2variable::pool->set_size(-1);
1250  }
1252  {
1253  delete vsparse_triplet_adjoint;
1255  }
1256  if (vsparse_triplet)
1257  {
1258  delete vsparse_triplet;
1259  vsparse_triplet=0;
1260  }
1261  if (sparse_triplet2)
1262  {
1263  delete sparse_triplet2;
1264  sparse_triplet2=0;
1265  }
1266  if (sparse_triplet)
1267  {
1268  delete sparse_triplet;
1269  sparse_triplet=0;
1270  }
1271  if (sparse_symbolic)
1272  {
1273  delete sparse_symbolic;
1274  sparse_symbolic=0;
1275  }
1276  if (sparse_symbolic2)
1277  {
1278  delete sparse_symbolic2;
1279  sparse_symbolic2=0;
1280  }
1281  if (num_local_re_array)
1282  {
1283  delete num_local_re_array;
1284  num_local_re_array = NULL;
1285  }
1287  {
1288  delete num_local_fixed_array;
1289  num_local_fixed_array = NULL;
1290  }
1291 }
1292 
1297 dvector laplace_approximation_calculator::operator()
1298  (const dvector& _x, const double& _f, function_minimizer * pfmin)
1299 {
1300  dvector g;
1301 
1302 #if defined(USE_ADPVM)
1303  if (pfmin->test_trust_flag)
1304  {
1305  return test_trust_region_method(_x,_f,pfmin);
1306  }
1308  {
1309  switch (ad_comm::pvm_manager->mode)
1310  {
1311  case 1:
1312  return default_calculations_parallel_master(_x,_f,pfmin);
1313  break;
1314  case 2:
1315  {
1316  dvector g(1,1);
1317  default_calculations_parallel_slave(_x,_f,pfmin);
1318  return g;
1319  break;
1320  }
1321  default:
1322  cerr << "illegal value for mode " << endl;
1323  ad_exit(1);
1324  }
1325  }
1326  else
1327 #endif //# if defined(USE_ADPVM)
1328 
1329  {
1330  //hesstype=1;
1331  //cout << hesstype << endl;
1332  switch (hesstype)
1333  {
1334  case 1:
1335  {
1336  int check_der_flag=0;
1337  int on=-1;
1338  int nopt=0;
1339  if ((on=option_match(ad_comm::argc,ad_comm::argv,"-chkder",nopt))>-1)
1340  {
1341  check_der_flag=1;
1342  }
1343  if (check_der_flag==1)
1344  {
1345  g = default_calculations_check_derivatives(_x,pfmin,_f);
1346  }
1347  else
1348  {
1349  g = default_calculations(_x,_f,pfmin);
1350  }
1351  break;
1352  }
1353  case 2:
1354  {
1355  // Hessian for the random effects is block diagonal
1356  g = block_diagonal_calculations(_x,_f,pfmin);
1357  break;
1358  }
1359  case 3:
1360  case 4: // not banded but full hessian
1361  {
1362  // Hessian for the random effects is block diagonal
1363  if (laplace_approximation_calculator::variance_components_vector)
1364  {
1365  g = banded_calculations_lme(_x,_f,pfmin);
1366  }
1367  else
1368  {
1369  g = banded_calculations(_x,_f,pfmin);
1370  }
1371  break;
1372  }
1373  default:
1374  {
1375  cerr << "illegal value for hesstype " << endl;
1376  ad_exit(1);
1377  }
1378  }
1379  }
1380 
1381  return g;
1382 }
1383 
1384 void random_effects_userfunction(double f,const dvector& x, const dvector& g);
1385 
1390 void get_second_ders(int xs,int us,const init_df1b2vector _y,dmatrix& Hess,
1393 {
1394  // Note that xs is the number of active non random effects
1395  // parameters
1396  // us is the number of random effects parameters
1397  int j;
1398  int i;
1400  int num_der_blocks=lpc->num_der_blocks;
1401  int xsize=lpc->xsize;
1402  int usize=lpc->usize;
1403 
1404  for (int ip=1;ip<=num_der_blocks;ip++)
1405  {
1406  df1b2variable::minder=lpc->minder(ip);
1407  df1b2variable::maxder=lpc->maxder(ip);
1408  lpc->set_u_dot(ip);
1410  (*re_objective_function_value::pobjfun)=0;
1411  df1b2variable pen=0.0;
1412  df1b2variable zz=0.0;
1414  //initial_params::straight_through_flag=1;
1418  {
1420  Hess.initialize();
1421  Dux.initialize();
1422  }
1423 
1424  //cout << "2D" << endl;
1425  pfmin->user_function();
1426 
1427  //pfmin->user_function(y,zz);
1428  (*re_objective_function_value::pobjfun)+=pen;
1429  (*re_objective_function_value::pobjfun)+=zz;
1431  {
1435  df1b2_gradcalc1();
1436  int mind=y(1).minder;
1437  int jmin=max(mind,xsize+1);
1438  int jmax=min(y(1).maxder,xsize+usize);
1439  for (i=1;i<=usize;i++)
1440  for (j=jmin;j<=jmax;j++){
1441  //cout<<i<<" "<<j<<" "<<mind<<" "<<xsize<<" "<<usize<<endl;
1442  Hess(i,j-xsize)=y(i+xsize).u_bar[j-mind];
1443  }
1444  jmin=max(mind,1);
1445  jmax=min(y(1).maxder,xsize);
1446  for (i=1;i<=usize;i++)
1447  for (j=jmin;j<=jmax;j++)
1448  Dux(i,j)=y(i+xsize).u_bar[j-1];
1449  }
1450  if (ip<num_der_blocks)
1451  f1b2gradlist->reset();
1452  }
1453 
1455  {
1456  for (i=1;i<=usize;i++)Hess(i,i) -= 2*CHECK_HESSIAN_PENALTY;
1457  cout<< "--------------------------------------------------------------------"<<endl;
1458  cout << "param:" << endl;
1459  for (i=1;i<=usize;i++)cout << y(i+xsize)<<" ";
1460  cout << endl;
1461  cout<< "--------------------------------------------------------------------"<<endl;
1462  cout << "Hess:" << endl;
1463  cout << Hess << endl;
1464  cout<< "--------------------------------------------------------------------"<<endl;
1465  ad_exit(0);
1466  }
1467 
1468  //cout << "Dux" << endl;
1469  //cout << setprecision(16) << Dux << endl;
1470 }
1471 
1477  const dmatrix& _Hess,const dvector& _xadjoint,const dvector& _uadjoint,
1478  const dmatrix& _Hessadjoint,function_minimizer * pmin)
1479 {
1480  ADUNCONST(dvector,xadjoint)
1481  ADUNCONST(dvector,uadjoint)
1482  ADUNCONST(dmatrix,Hessadjoint)
1483  ADUNCONST(dmatrix,Hess)
1484 #if !defined(OPT_LIB) && (__cplusplus >= 201103L)
1485  const int xs = [](unsigned int size)->int
1486  {
1487  assert(size <= INT_MAX);
1488  return static_cast<int>(size);
1489  }(x.size());
1490  const int us = [](unsigned int size)->int
1491  {
1492  assert(size <= INT_MAX);
1493  return static_cast<int>(size);
1494  }(u0.size());
1495 #else
1496  const int xs = static_cast<int>(x.size());
1497  const int us = static_cast<int>(u0.size());
1498 #endif
1500  int nvar;
1501  //dcompressed_triplet & lst2 = *(pmin->lapprox->sparse_triplet);
1502  //hs_symbolic & ssymb=*(pmin->lapprox->sparse_symbolic);
1503  //dcompressed_triplet & xxxt = *(pmin->lapprox->sparse_triplet);
1504  dcompressed_triplet & lst = *(pmin->lapprox->sparse_triplet2);
1505  hs_symbolic & ssymb=*(pmin->lapprox->sparse_symbolic2);
1506  {
1507  /*
1508  cout << norm2(make_dmatrix(lst)-make_dmatrix(xxxt)) << endl;
1509  ofstream ofs1("m1");
1510  ofs1 << setfixed() << setprecision(3) << setw(10)
1511  << make_dmatrix(lst) << endl;
1512  ofstream ofs2("m2");
1513  ofs2 << setfixed() << setprecision(3) << setw(10)
1514  << make_dmatrix(xxxt) << endl;
1515  */
1516  }
1517 
1519  {
1520  if (pmin->lapprox->sparse_hessian_flag==0)
1521  {
1522  nvar = xs + us + us * us;
1523  }
1524  else
1525  {
1526  int sz= lst.indexmax()-lst.indexmin()+1;
1527  nvar = xs + us + sz;
1528  }
1529  }
1530  else
1531  {
1532  nvar = xs + us;
1533  }
1534  independent_variables y(1,nvar);
1535 
1536  // need to set random effects active together with whatever
1537  // init parameters should be active in this phase
1540  /*int onvar=*/initial_params::nvarcalc();
1541  initial_params::xinit(y); // get the initial values into the
1542  y(1,xs)=x;
1543 
1544  int i,j;
1545  dvar_vector d(1,xs+us);
1546 
1547  dmatrix Hess_save;
1548  // contribution for quadratic prior
1550  {
1551  if (allocated(Hess_save)) Hess_save.deallocate();
1552  int mmin=Hess.indexmin();
1553  int mmax=Hess.indexmax();
1554  Hess_save.allocate(mmin,mmax,mmin,mmax);
1555  Hess_save=Hess;
1556  int & vxs = (int&)(xs);
1558  }
1559 
1560  dvar_matrix vHess;
1561  int ii=xs+us+1;
1563  {
1564  dvar_vector vy=dvar_vector(y);
1565  initial_params::reset(vy); // get the initial values into the
1566  }
1567  else
1568  {
1569  // Here need hooks for sparse matrix structures
1570  if (pmin->lapprox->sparse_hessian_flag==0)
1571  {
1572  for (i=1;i<=us;i++)
1573  for (j=1;j<=us;j++)
1574  y(ii++)=Hess(i,j);
1575  }
1576  else
1577  {
1578  int smin=lst.indexmin();
1579  int smax=lst.indexmax();
1580  for (i=smin;i<=smax;i++)
1581  y(ii++)=lst(i);
1582  }
1583 
1585  {
1586  Hess=Hess_save;
1587  }
1588 
1589  dvar_vector vy=dvar_vector(y);
1591  if (pmin->lapprox->sparse_hessian_flag==0)
1592  {
1593  vHess.allocate(1,us,1,us);
1594  }
1595  initial_params::reset(vy); // get the initial values into the
1596  ii=xs+us+1;
1598  {
1599  for (i=1;i<=us;i++)
1600  {
1601  if (d(i+xs)>0.0)
1602  {
1603  for (j=1;j<=us;j++)
1604  {
1605  if (d(j+xs)>0.0)
1606  vHess(i,j)=vy(ii++)/(d(i+xs)*d(j+xs));
1607  else
1608  vHess(i,j)=vy(ii++)/d(i+xs);
1609  }
1610  }
1611  else
1612  {
1613  for (j=1;j<=us;j++)
1614  {
1615  if (d(j+xs)>0.0)
1616  vHess(i,j)=vy(ii++)/d(j+xs);
1617  else
1618  vHess(i,j)=vy(ii++);
1619  }
1620  }
1621  }
1622  }
1623  else
1624  {
1626  {
1627  for (i=1;i<=us;i++)
1628  {
1629  for (j=1;j<=us;j++)
1630  {
1631  vHess(i,j)=vy(ii++);
1632  }
1633  }
1634  }
1635  else
1636  {
1637  int mmin=lst.indexmin();
1638  int mmax=lst.indexmax();
1639  dvar_compressed_triplet * vsparse_triplet =
1640  pmin->lapprox->vsparse_triplet;
1641 
1642  if (vsparse_triplet==0)
1643  {
1644  pmin->lapprox->vsparse_triplet=
1645  new dvar_compressed_triplet(mmin,mmax,us,us);
1646  vsparse_triplet = pmin->lapprox->vsparse_triplet;
1647  for (i=mmin;i<=mmax;i++)
1648  {
1649  (*vsparse_triplet)(1,i)=lst(1,i);
1650  (*vsparse_triplet)(2,i)=lst(2,i);
1651  }
1652  }
1653  else
1654  {
1655  if (!allocated(*vsparse_triplet))
1656  {
1657  (*vsparse_triplet).allocate(mmin,mmax,us,us);
1658  for (i=mmin;i<=mmax;i++)
1659  {
1660  (*vsparse_triplet)(1,i)=lst(1,i);
1661  (*vsparse_triplet)(2,i)=lst(2,i);
1662  }
1663  }
1664  }
1665  dcompressed_triplet * vsparse_triplet_adjoint =
1667 
1668  if (vsparse_triplet_adjoint==0)
1669  {
1671  new dcompressed_triplet(mmin,mmax,us,us);
1672  vsparse_triplet_adjoint = pmin->lapprox->vsparse_triplet_adjoint;
1673  for (i=mmin;i<=mmax;i++)
1674  {
1675  (*vsparse_triplet_adjoint)(1,i)=lst(1,i);
1676  (*vsparse_triplet_adjoint)(2,i)=lst(2,i);
1677  }
1678  }
1679  else
1680  {
1681  if (!allocated(*vsparse_triplet_adjoint))
1682  {
1683  (*vsparse_triplet_adjoint).allocate(mmin,mmax,us,us);
1684  for (i=mmin;i<=mmax;i++)
1685  {
1686  (*vsparse_triplet_adjoint)(1,i)=lst(1,i);
1687  (*vsparse_triplet_adjoint)(2,i)=lst(2,i);
1688  }
1689  }
1690  }
1691  vsparse_triplet->get_x()=vy(ii,ii+mmax-mmin).shift(1);
1692  }
1693  }
1694  }
1695 
1696  dvariable vf=0.0;
1697 
1700  {
1702  }
1703  pmin->AD_uf_outer();
1705  {
1707  }
1709  {
1711  }
1713  //cout << setprecision(15) << vf << endl;
1714  // *********************************************
1715  // *********************************************
1716  // *********************************************
1717  // dvector tmpg(1,nvar);
1718  // tmpg.initialize();
1719  // gradcalc(nvar,tmpg);
1720  // *********************************************
1721  // *********************************************
1722  // *********************************************
1723 
1724  int sgn=0;
1725  dvariable ld = 0;
1726  double f=0.0;
1728  {
1730  {
1732  {
1733  ld=0.5*ln_det(vHess,sgn);
1734  }
1735  else
1736  {
1737  ld=0.5*ln_det(-vHess,sgn);
1738  }
1739  }
1740  else
1741  {
1743  {
1744  int ierr=0;
1745  if (pmin->lapprox->sparse_hessian_flag==0)
1746  {
1747  ld=0.5*ln_det_choleski_error(vHess,ierr);
1748  if (ierr==1)
1749  {
1750  ofstream ofs("hessian.diag");
1751  ofs << vHess << endl;
1752  ofs << x << endl;
1753  ofs << u0 << endl;
1754  ofs << "Matrix not positive definite in Ln_det_choleski"
1755  << endl;
1756  cerr << "Matrix not positive definite in Ln_det_choleski\n"
1757  << "see file hessian.diag for details"
1758  << endl;
1759  ad_exit(1);
1760  }
1761  }
1762  else
1763  {
1764  //double ld1=0.5*ln_det(*(pmin->lapprox->sparse_triplet),
1765  // *(pmin->lapprox->sparse_symbolic));
1766 
1767 #ifdef DIAG
1768  if (write_sparse_flag)
1769  {
1770  //ofstream ofs("sparse");
1771  //ofs << *(pmin->lapprox->vsparse_triplet) << endl;
1772  }
1773 #endif
1774  ld=0.5*ln_det(*(pmin->lapprox->vsparse_triplet),
1775  ssymb,*(pmin->lapprox->sparse_triplet2));
1776  //*(pmin->lapprox->sparse_symbolic2),pmin->lapprox);
1777  //cout << ld-ld1 << endl;
1778  }
1779  }
1780  else
1781  {
1782  if (pmin->lapprox->sparse_hessian_flag==0)
1783  {
1784  ld=0.5*ln_det_choleski(-vHess);
1785  }
1786  else
1787  {
1788  cerr << "need to fix this" << endl;
1789  ad_exit(1);
1790  //ld+=ln_det(-*(pmin->lapprox->vsparse_triplet));
1791  }
1792  }
1793  }
1794 
1795 #ifdef DIAG
1796  int ps1=0;
1797  if (ps1)
1798  {
1799  dmatrix cHess=value(vHess);
1800  cout << " ln_det = " << ld << " ";
1801  dvector eig=eigenvalues(cHess);
1802  dmatrix r(1,2,1,eig.indexmax());
1803  dvector num(1,eig.indexmax());
1804  num.fill_seqadd(1,1);
1805  r(1)=eig;
1806  r(2)=num;
1807  dmatrix tt=trans(r);
1808  dmatrix ss=trans(sort(tt,1));
1809  cout << ss(1,1) << " " << ss(1,eig.indexmax()) << " " ;
1810  int nnn=(int)ss(2,1);
1811  dmatrix d=eigenvectors(cHess);
1812  //cout << " d " << d(nnn) << " " << d(nnn)*cHess*d(nnn) << endl;
1813  dmatrix t=trans(d);
1814  r(1)=t(nnn);
1815  r(2)=num;
1816  dmatrix tt2=trans(r);
1817  dmatrix ss2=trans(sort(tt2,1));
1818 
1819  cout << " t " << setprecision(3) << ss2(1)(1,5) << " --- "
1820  << t(nnn)*cHess*t(nnn) << endl;
1821  cout << " " << setprecision(3) << ss2(2)(1,5) << endl;
1822  //cout << " t " << t(1) << " " << t(1)*cHess*t(2) << endl;
1823  }
1824 #endif
1825 
1826  int nx=0;
1827  if (nx==0)
1828  {
1829  vf+=ld;
1830  }
1831  f=value(vf);
1832  f-=us*0.5*log(2.0*PI);
1833  }
1834  else
1835  {
1836  f=value(vf);
1837  }
1838 
1839  dvector g(1,nvar);
1840  gradcalc(nvar,g);
1841 
1842  ii=1;
1843  for (i=1;i<=xs;i++)
1844  xadjoint(i)=g(ii++);
1845  for (i=1;i<=us;i++)
1846  uadjoint(i)=g(ii++);
1847 
1849  {
1850  if (pmin->lapprox->sparse_hessian_flag==0)
1851  {
1852  if (allocated(Hessadjoint))
1853  {
1854  for (i=1;i<=us;i++)
1855  for (j=1;j<=us;j++)
1856  Hessadjoint(i,j)=g(ii++);
1857  }
1858  }
1859  else
1860  {
1861  dcompressed_triplet * vsparse_triplet_adjoint =
1863 
1864  int smin=lst.indexmin();
1865  int smax=lst.indexmax();
1866  for (i=smin;i<=smax;i++)
1867  {
1868  (*vsparse_triplet_adjoint)(i)=g(ii);
1869  ii++;
1870  }
1871  }
1872  }
1873 
1874  /*
1875  if (quadratic_prior::get_num_quadratic_prior()>0)
1876  {
1877  // *******************************************************
1878  // *******************************************************
1879  // for quadratic prior option
1880  // temporarily get der wrt x of x ---> ln_det(F_uu + inv(A(x)))
1881  int nvar=x.size()+u0.size();
1882  independent_variables y(1,nvar);
1883  initial_params::xinit(y); // get the initial values into the
1884  y(1,xs)=x;
1885  dvar_vector vy=dvar_vector(y);
1886  initial_params::reset(vy); // get the initial values into the
1887 
1888  pmin->AD_uf_outer();
1889 
1890  dvar_matrix tmpH=quadratic_prior::get_vHessian_contribution();
1891 
1892  tmpH+=Hess;
1893  dvariable ld;
1894  ld=0.5*ln_det(tmpH,sgn);
1895  dvector g(1,nvar);
1896  gradcalc(nvar,g);
1897  int ii=1;
1898  double checksum=0.0;
1899  for (i=1;i<=xs;i++)
1900  xadjoint(i)+=g(ii++);
1901  for (i=1;i<=us;i++)
1902  checksum+=g(ii++);
1903 
1904  if (fabs(checksum)> 1.e-10)
1905  {
1906  cerr << "checksum too big " << checksum << endl;
1907  ad_exit(1);
1908  }
1909 
1910  }
1911  */
1912 
1913 
1914  // *******************************************************
1915  // *******************************************************
1916 
1917  return f;
1918 }
1919 
1924 void get_newton_raphson_info(int xs,int us,const init_df1b2vector _y,
1925  dmatrix& Hess, dvector& grad, df1b2_gradlist * f1b2gradlist,
1926  function_minimizer * pfmin)
1927 {
1928  // Note that xs is the number of active non random effects
1929  // parameters
1930  // us is the number of random effects parameters
1931  int j;
1932  int i;
1934  {
1936  //cout << re_objective_function_value::pobjfun << endl;
1937  //cout << re_objective_function_value::pobjfun->ptr << endl;
1938  (*re_objective_function_value::pobjfun)=0;
1939  df1b2variable pen=0.0;
1940  df1b2variable zz=0.0;
1942  pfmin->user_function();
1943  //pfmin->user_function(y,zz);
1944  (*re_objective_function_value::pobjfun)+=pen;
1945  (*re_objective_function_value::pobjfun)+=zz;
1949  }
1950 
1951  int mind=y(1).minder;
1952  int jmin=max(mind,xs+1);
1953  int jmax=min(y(1).maxder,xs+us);
1955  {
1956  df1b2_gradcalc1();
1957  for (i=1;i<=us;i++)
1958  for (j=jmin;j<=jmax;j++)
1959  Hess(i,j-xs)=y(i+xs).u_bar[j-mind];
1960  for (j=jmin;j<=jmax;j++)
1961  grad(j-xs)= re_objective_function_value::pobjfun->u_dot[j-mind];
1962  }
1963 }
1964 
1969 {
1970 }
1971  //void laplace_approximation_calculator::get_newton_raphson_info
1972  // (function_minimizer * pfmin)
1973  //{
1974  // int i,j,ip;
1975  //
1976  // for (ip=1;ip<=num_der_blocks;ip++)
1977  // {
1978  // // do we need to reallocate memory for df1b2variables?
1979  // check_for_need_to_reallocate(ip);
1980  // df1b2_gradlist::set_yes_derivatives();
1981  // //cout << re_objective_function_value::pobjfun << endl;
1982  // //cout << re_objective_function_value::pobjfun->ptr << endl;
1983  // (*re_objective_function_value::pobjfun)=0;
1984  // df1b2variable pen=0.0;
1985  // df1b2variable zz=0.0;
1986  // initial_df1b2params::reset(y,pen);
1987  // if (initial_df1b2params::separable_flag)
1988  // {
1989  // Hess.initialize();
1990  // grad.initialize();
1991  // }
1992  // pfmin->user_function();
1993  // //pfmin->user_function(y,zz);
1994  // (*re_objective_function_value::pobjfun)+=pen;
1995  // (*re_objective_function_value::pobjfun)+=zz;
1996  //
1997  // if (!initial_df1b2params::separable_flag)
1998  // {
1999  // set_dependent_variable(*re_objective_function_value::pobjfun);
2000  // df1b2_gradlist::set_no_derivatives();
2001  // df1b2variable::passnumber=1;
2002  // df1b2_gradcalc1();
2003  // int mind=y(1).minder;
2004  // int jmin=max(mind,xsize+1);
2005  // int jmax=min(y(1).maxder,xsize+usize);
2006  // for (i=1;i<=usize;i++)
2007  // for (j=jmin;j<=jmax;j++)
2008  // Hess(i,j-xsize)=y(i+xsize).u_bar[j-mind];
2009  //
2010  // for (j=jmin;j<=jmax;j++)
2011  // grad(j-xsize)= re_objective_function_value::pobjfun->u_dot[j-mind];
2012  // }
2013  // }
2014  // {
2015  // //ofstream ofs("hess.tmp");
2016  // //ofs << Hess << endl;
2017  // //ofs << grad << endl;
2018  // }
2019  //}
2020  //
2021 
2027 {
2028  int usize=initial_params::nvarcalc();
2029  //double f=0.0;
2030  dvector g(1,usize);
2031  independent_variables u(1,usize);
2032  u=x;
2033  dvariable vf=0.0;
2035  //vf=0.0;
2038  pfmin->AD_uf_inner();
2040  {
2042  }
2046  gradcalc(usize,g);
2047  double maxg=max(fabs(g));
2049  {
2050  std::ostream& output_stream = get_output_stream();
2051  std::streamsize save = output_stream.precision();
2052  output_stream << std::fixed << setprecision(10) << " f = " << vf
2053  << " max g = " << std::scientific << maxg << endl;
2054  output_stream.precision(save);
2055  }
2056  return maxg;
2057 }
2058 
2063 double evaluate_function(double& fval,const dvector& x,
2064  function_minimizer* pfmin)
2065 {
2066  int usize=initial_params::nvarcalc();
2067  //double f=0.0;
2068  dvector g(1,usize);
2069  independent_variables u(1,usize);
2070  u=x;
2071  dvariable vf=0.0;
2073  //vf=0.0;
2076  pfmin->AD_uf_inner();
2078  {
2080  }
2084  gradcalc(usize,g);
2085  double maxg=max(fabs(g));
2086  fval=value(vf);
2088  {
2089  std::ostream& output_stream = get_output_stream();
2090  output_stream << std::fixed << setprecision(10) << " f = " << vf
2091  << " max g = " << std::scientific << maxg << endl;
2092  }
2093  return maxg;
2094 }
2095 
2100 double evaluate_function(double& fval,const dvector& x,const dvector& g,
2101  function_minimizer * pfmin)
2102 {
2103  int usize=initial_params::nvarcalc();
2104  //double f=0.0;
2105  //dvector g(1,usize);
2106  independent_variables u(1,usize);
2107  u=x;
2108  dvariable vf=0.0;
2110  //vf=0.0;
2113  pfmin->AD_uf_inner();
2115  {
2117  }
2121  gradcalc(usize,g);
2122  double maxg=max(fabs(g));
2123  fval=value(vf);
2125  {
2126  std::ostream& output_stream = get_output_stream();
2127  output_stream << std::fixed << setprecision(10) << " f = " << vf
2128  << " max g = " << std::scientific << maxg << endl;
2129  }
2130  return maxg;
2131 }
2132 
2138 {
2139  int usize=initial_params::nvarcalc();
2140  //double f=0.0;
2141  dvector g(1,usize);
2142  independent_variables u(1,usize);
2143  u=x;
2144  dvariable vf=0.0;
2146  //vf=0.0;
2149  pfmin->AD_uf_inner();
2151  {
2153  }
2157  gradcalc(usize,g);
2158  double maxg=max(fabs(g));
2159  return maxg;
2160 }
2161 
2166 void evaluate_function_gradient(double& f,const dvector& x,
2167  function_minimizer * pfmin,dvector& xadjoint,dvector& uadjoint)
2168 {
2169  int usize=initial_params::nvarcalc();
2170  dvector g(1,usize);
2171  independent_variables u(1,usize);
2172  u=x;
2173  dvariable vf=0.0;
2175  //vf=0.0;
2177  pfmin->AD_uf_inner();
2179  {
2181  }
2184  f=value(vf);
2185  gradcalc(usize,g);
2186  int xsize=xadjoint.indexmax();
2187  int us=uadjoint.indexmax();
2188  xadjoint=g(1,xsize);
2189  uadjoint=g(xsize+1,xsize+us).shift(1);
2190 }
2191 
2197  function_minimizer* pfmin)
2198 {
2199  double fval;
2201  int usize=initial_params::nvarcalc();
2202  //double f=0.0;
2203  dvector g(1,usize);
2204  independent_variables u(1,usize);
2205  u=x;
2206  dvariable vf=0.0;
2208  //vf=0.0;
2210  pfmin->AD_uf_inner();
2212  {
2214  }
2217  fval=value(vf);
2219  return fval;
2220 }
2221 
2227 {
2228  if (l)
2229  {
2230  delete l;
2231  l=0;
2233 
2234  for (int i=0;i<df1b2variable::adpool_counter;i++)
2235  {
2236  delete df1b2variable::adpool_vector[i];
2239  }
2240 
2241  df1b2variable::adpool_counter=0;
2242  }
2243 }
2244 
2251  int use_gauss_hermite,int num_separable_calls ,const ivector& itmp) :
2252  x(1,use_gauss_hermite), w(1,use_gauss_hermite), mi(0)
2253 {
2254  // for now we must have the same number of random effects in each
2255  // separable function call"
2256  int i;
2257  for (i=2;i<=num_separable_calls;i++)
2258  {
2259  if (itmp(i)!=itmp(i-1))
2260  {
2261  cerr << " At present for the adaptive gauss_hermite must have the same"
2262  << endl
2263  << " number of random effects in each separable function call"
2264  << endl;
2265  ad_exit(1);
2266  }
2267  }
2268  for (i=1;i<=num_separable_calls;i++)
2269  {
2270  if (itmp(i)>1)
2271  {
2272  lapprox->multi_random_effects=
2273  max(lapprox->multi_random_effects,itmp(i));
2274  /*
2275  cerr << "At present gauss-hermite is only implemented for"
2276  " one random effect per separable function call "
2277  << endl;
2278  ad_exit(1);
2279  */
2280  }
2281  }
2282  if (allocated(gauss_hermite_values))
2283  gauss_hermite_values.deallocate();
2284  if (lapprox->multi_random_effects==0)
2285  {
2286  gauss_hermite_values.allocate(1,num_separable_calls,1,use_gauss_hermite);
2287  }
2288  else
2289  {
2290  ivector indx=pow(use_gauss_hermite,itmp);
2291  gauss_hermite_values.allocate(1,num_separable_calls,1,indx);
2292  }
2293 
2295  is=0;
2296 }
2297 
2303 {
2304  //if (function_minimizer::hesstype==0)
2305  {
2307  {
2308  // have SEPARABLE_FUNCTION but no object of type quadratic_penalty or
2309  // normal prior -- can't tell if we should use hesstype 2 or 3
2310  // until we run
2311  hesstype=3;
2312  }
2313  else
2314  {
2315  hesstype=1; // assume block diagonal
2316  }
2317  }
2318  /*
2319  else
2320  {
2321  hesstype=function_minimizer::hesstype;
2322  }
2323  */
2324 }
2325 
2331 {
2334 
2335  if (grad_x==0)
2336  {
2337  grad_x=new dvector(1,xsize);
2338  }
2339  grad_x->initialize();
2340 
2341  if (grad_x_u==0)
2342  {
2343  grad_x_u=new dvector(1,xsize+usize);
2344  }
2345  pfmin->inner_opt_flag=0;
2347  //gradcalc(0,*grad_x_u);
2348  initial_params::xinit(u); // get the initial values into the
2349 
2350  dvariable pen=0.0;
2351  dvariable vf=0.0;
2354  pfmin->AD_uf_inner();
2357  vf+=pen;
2359  double f=value(vf);
2360  return f;
2361 }
2362 
2368  {
2370  //build_up_nested_shape();
2371  //clean(nested_tree_position,separable_call_level);
2372  //nested_separable_calls_counter(separable_call_level)++;
2373  //nested_tree_position(separable_call_level)++;
2374  }
2375 
2381  {
2382  //build_up_nested_shape();
2383  //clean(nested_tree_position,separable_call_level);
2385  }
2386 
2392 {
2393  int ll;
2395  int clean_level=0;
2396  switch(separable_call_level)
2397  {
2398  case 1:
2399  ll=1;
2401  {
2402  if (clean_level==0) clean_level=ll+1;
2403  if (nested_separable_calls_counter(ll+1)>0)
2404  {
2405  nested_shape(a(ll))=a(ll+1);
2406  }
2407  else
2408  {
2409  break;
2410  }
2411  }
2412  else
2413  {
2414  break;
2415  }
2416  case 2:
2417  ll=2;
2419  {
2420  if (clean_level==0) clean_level=ll+1;
2421  if (nested_separable_calls_counter(ll+1)>0)
2422  {
2423  nested_shape(a(1),a(2))=a(ll+1);
2424  }
2425  else
2426  {
2427  break;
2428  }
2429  }
2430  else
2431  {
2432  break;
2433  }
2434  case 3:
2435  ll=3;
2437  {
2438  if (clean_level==0) clean_level=ll+1;
2439  if (nested_separable_calls_counter(ll+1)>0)
2440  {
2441  nested_shape(a(1),a(2),a(3))=a(ll+1);
2442  }
2443  else
2444  {
2445  break;
2446  }
2447  }
2448  else
2449  {
2450  break;
2451  }
2452  case 4:
2453  ll=4;
2455  {
2456  if (clean_level==0) clean_level=ll+1;
2457  if (nested_separable_calls_counter(ll+1)>0)
2458  {
2459  nested_shape(a(1),a(2),a(3),a(4))=a(ll+1);
2460  }
2461  else
2462  {
2463  break;
2464  }
2465  }
2466  else
2467  {
2468  break;
2469  }
2470  default:
2471  cerr << "illegal value in " <<
2472  "laplace_approximation_calculator::build_up_nested_shape"
2473  << endl;
2474  }
2475  if (clean_level>0)
2476  {
2477  int mmax=a.indexmax();
2478  for (int i=clean_level;i<=mmax;i++)
2479  {
2480  a(i)=0;
2481  }
2482  }
2483 }
2484 
2490 {
2491  return ptr1;
2492 }
2493 
2499 {
2500  return ptr2;
2501 }
2502 
2508 {
2509  return ptr3;
2510 }
2511 
2517 {
2518  return ptr4;
2519 }
2520 
2525 ostream & operator << (const ostream& _s,const nested_calls_shape& _m)
2526 {
2528  ADUNCONST(ofstream,s)
2529  if (m.get_ptr1())
2530  s<< *(m.get_ptr1()) << endl << endl;
2531 
2532  if (m.get_ptr2())
2533  s<< *(m.get_ptr2()) << endl << endl;
2534 
2535  if (m.get_ptr3())
2536  s<< *(m.get_ptr3()) << endl << endl;
2537 
2538  if (m.get_ptr4())
2539  s<< *(m.get_ptr4()) << endl << endl;
2540 
2541  return s;
2542 }
2543 
2549 {
2550  int mmax1=0;
2551  if (ptr1)
2552  {
2553  int imin=ptr1->indexmin();
2554  int imax=ptr1->indexmax();
2555  for (int i=imin;i<=imax;i++)
2556  {
2557  if ( (*ptr1)(i)==0) break;
2558  mmax1++;
2559  }
2560  }
2561  if (mmax1==0)
2562  {
2563  delete ptr1;
2564  ptr1=0;
2565  }
2566  else
2567  {
2568  ivector * tmp=new ivector(1,mmax1);
2569  (*tmp)=(*ptr1)(1,mmax1);
2570  delete ptr1;
2571  ptr1=tmp;
2572  }
2573  if (ptr2)
2574  {
2575  if (!ptr1)
2576  {
2577  delete ptr2;
2578  ptr2=0;
2579  }
2580  else
2581  {
2582  imatrix * tmp=new imatrix(1,mmax1);
2583  int imin=tmp->indexmin();
2584  int imax=tmp->indexmax();
2585  for (int i=imin;i<=imax;i++)
2586  {
2587  int jmin=(*ptr2)(i).indexmin();
2588  int jmax=(*ptr2)(i).indexmax();
2589  int mmax2=0;
2590  for (int j=jmin;j<=jmax;j++)
2591  {
2592  if ((*ptr2)(i,j)==0) break;
2593  mmax2++;
2594  }
2595  if (mmax2>0)
2596  {
2597  (*tmp)(i).allocate(1,mmax2);
2598  (*tmp)(i)=(*ptr2)(i)(1,mmax2);
2599  }
2600  }
2601  delete ptr2;
2602  ptr2=tmp;
2603  }
2604  }
2605  if (ptr3)
2606  {
2607  cerr << "warning not dealitn with prt3" << endl;
2608  delete ptr3;
2609  ptr3=0;
2610  }
2611 }
2612 
2618 {
2619  if (ptr1)
2620  {
2621  delete ptr1;
2622  ptr1=0;
2623  }
2624  if (ptr2)
2625  {
2626  delete ptr2;
2627  ptr2=0;
2628  }
2629  if (ptr3)
2630  {
2631  delete ptr3;
2632  ptr3=0;
2633  }
2634  if (ptr4)
2635  {
2636  delete ptr4;
2637  ptr4=0;
2638  }
2639 }
2640 
2645 void nested_calls_shape::allocate(int n,int m,int p )
2646 {
2647  if (ptr1)
2648  {
2649  delete ptr1;
2650  ptr1=0;
2651  }
2652  ptr1=new ivector(1,n);
2653 
2654  if (ptr2)
2655  {
2656  delete ptr2;
2657  ptr2=0;
2658  }
2659  ptr2=new imatrix(1,n,1,m);
2660  if (ptr3)
2661  {
2662  delete ptr3;
2663  ptr3=0;
2664  }
2665  ptr3=new i3_array(1,n,1,m,1,p);
2666  /*
2667  if (ptr4)
2668  {
2669  delete ptr4;
2670  ptr4=0;
2671  }
2672  */
2673 }
2674 
2680 {
2681  if (ptr1)
2682  {
2683  delete ptr1;
2684  ptr1=0;
2685  }
2686  ptr1=new ivector(1,n);
2687 
2688  if (ptr2)
2689  {
2690  delete ptr2;
2691  ptr2=0;
2692  }
2693  ptr2=new imatrix(1,n,1,m);
2694  /*
2695  if (ptr3)
2696  {
2697  delete ptr3;
2698  ptr3=0;
2699  }
2700  ptr=new i3_array(1,n,1,m,1,p);
2701  if (ptr4)
2702  {
2703  delete ptr4;
2704  ptr4=0;
2705  }
2706  */
2707 }
2708 
2714 {
2715  if (ptr1)
2716  {
2717  ptr1->initialize();
2718  }
2719 
2720  if (ptr2)
2721  {
2722  ptr2->initialize();
2723  }
2724 
2725  if (ptr3)
2726  {
2727  ptr3->initialize();
2728  }
2729 
2730  if (ptr4)
2731  {
2732  ptr4->initialize();
2733  }
2734 }
2735 
2741 {
2742  if (ptr1)
2743  {
2744  delete ptr1;
2745  ptr1=0;
2746  }
2747  ptr1=new ivector(1,n);
2748 
2749  /*
2750  if (ptr2)
2751  {
2752  delete ptr2;
2753  ptr2=0;
2754  }
2755  ptr=new imatrix(1,n,1,m);
2756  if (ptr3)
2757  {
2758  delete ptr3;
2759  ptr3=0;
2760  }
2761  ptr=new i3_array(1,n,1,m,1,p);
2762  if (ptr4)
2763  {
2764  delete ptr4;
2765  ptr4=0;
2766  }
2767  */
2768 }
2769 
2775 {
2776  int mmax=0;
2778  mmax=nsc.get_ptr1()->indexmax();
2779  if (ptr1)
2780  {
2781  delete ptr1;
2782  ptr1=0;
2783  }
2784  if (nsc.get_ptr1())
2785  {
2786  ptr1=new imatrix(1,mmax,1,*nsc.get_ptr1());
2787  }
2788  if (ptr2)
2789  {
2790  delete ptr2;
2791  ptr2=0;
2792  }
2793  if (nsc.get_ptr2())
2794  {
2795  ptr2=new i3_array(1,mmax,1,*nsc.get_ptr1(),1,*nsc.get_ptr2());
2796  }
2797  if (ptr3)
2798  {
2799  delete ptr3;
2800  ptr3=0;
2801  }
2802  if (nsc.get_ptr3())
2803  {
2804  ptr3=new i4_array(1,mmax,1,*nsc.get_ptr1(),1,*nsc.get_ptr2(),
2805  1,*nsc.get_ptr3());
2806  }
2807 }
2808 
2814  (const dvector& x,function_minimizer * pfmin)
2815 {
2816  std::ostream& output_stream = get_output_stream();
2817 
2818  //int on,nopt;
2819  pfmin->inner_opt_flag=1;
2820  double f=0.0;
2821  double fb=1.e+100;
2822  dvector g(1,usize);
2823  dvector ub(1,usize);
2824  independent_variables u(1,usize);
2825  gradcalc(0,g);
2826  fmc1.itn=0;
2827  fmc1.ifn=0;
2828  fmc1.ireturn=0;
2829  initial_params::xinit(u); // get the initial values into the
2830  initial_params::xinit(ubest); // get the initial values into the
2831  fmc1.ialph=0;
2832  fmc1.ihang=0;
2833  fmc1.ihflag=0;
2834 
2835  if (init_switch)
2836  {
2837  u.initialize();
2838  }
2839  else
2840  {
2841  u=ubest;
2842  }
2843 
2844  fmc1.dfn=1.e-2;
2845  dvariable pen=0.0;
2846  //cout << "starting norm(u) = " << norm(u) << endl;
2847  while (fmc1.ireturn>=0)
2848  {
2849  /*
2850  double nu=norm(value(u));
2851  if (nu>400)
2852  {
2853  cout << "U norm(u) = " << nu << endl;
2854  }
2855  cout << "V norm(u) = " << nu
2856  << " f = " << f << endl;
2857  */
2858  fmc1.fmin(f,u,g);
2859  //cout << "W norm(u) = " << norm(value(u)) << endl;
2860  if (fmc1.ireturn>0)
2861  {
2862  dvariable vf=0.0;
2865  pfmin->AD_uf_inner();
2866  if (saddlepointflag)
2867  {
2869  }
2871  {
2873  }
2875 
2876  /* this is now done in the operator = function
2877  if (quadratic_prior::get_num_quadratic_prior()>0)
2878  {
2879  vf+= quadratic_prior::get_quadratic_priors();
2880  }
2881  */
2882 
2884 
2885  //cout << " pen = " << pen << endl;
2886  if (noboundepen_flag==0)
2887  {
2888  vf+=pen;
2889  }
2890  f=value(vf);
2891  if (f<fb)
2892  {
2893  fb=f;
2894  ub=u;
2895  }
2896  gradcalc(usize,g);
2897  //cout << " f = " << setprecision(17) << f << " " << norm(g)
2898  // << " " << norm(u) << endl;
2899  }
2900  u=ub;
2901  }
2902  output_stream << " inner maxg = " << std::scientific << setprecision(10) << fmc1.gmax;
2903 #ifdef DIAG
2904  if (fabs(fmc1.gmax)>1.e+3)
2905  trapper();
2906 #endif
2907 
2908  if (fabs(fmc1.gmax)>1.e-4)
2909  {
2910  fmc1.itn=0;
2911  //fmc1.crit=1.e-9;
2912  fmc1.ifn=0;
2913  fmc1.ireturn=0;
2914  fmc1.ihang=0;
2915  fmc1.ihflag=0;
2916  fmc1.ialph=0;
2917  initial_params::xinit(u); // get the initial values into the
2918  //u.initialize();
2919  while (fmc1.ireturn>=0)
2920  {
2921  fmc1.fmin(f,u,g);
2922  if (fmc1.ireturn>0)
2923  {
2924  dvariable vf=0.0;
2927  pfmin->AD_uf_inner();
2929  {
2931  }
2934 
2935  vf+=pen;
2936  f=value(vf);
2937  if (f<fb)
2938  {
2939  fb=f;
2940  ub=u;
2941  }
2942  gradcalc(usize,g);
2943  //cout << " f = " << setprecision(15) << f << " " << norm(g) << endl;
2944  }
2945  }
2946  u=ub;
2947  output_stream << " Inner second time = " << std::scientific << setprecision(10) << fmc1.gmax;
2948  }
2949  output_stream << " Inner f = " << std::scientific << setprecision(10) << fb << endl;
2950  fmc1.ireturn=0;
2951  fmc1.fbest=fb;
2954  pfmin->AD_uf_inner();
2956  {
2958  }
2960  pfmin->inner_opt_flag=0;
2961  return u;
2962 }
static int straight_through_flag
Definition: admodel.h:839
Description not yet available.
Definition: fvar.hpp:7981
void cleanup_laplace_stuff(laplace_approximation_calculator *)
Description not yet available.
Definition: df1b2lap.cpp:2226
unsigned int global_nvar
Definition: df1b2lap.cpp:32
static adpool * pool
Definition: df1b2fun.h:273
i3_array * ptr3
Definition: adrndeff.h:154
laplace_approximation_calculator * lapprox
Definition: admodel.h:1862
static dvar_vector * variance_components_vector
Definition: adrndeff.h:204
static void reset(const init_df1b2vector &, const df1b2variable &)
Description not yet available.
Definition: df1b2f15.cpp:51
dcompressed_triplet * sparse_triplet2
Definition: adrndeff.h:271
int noboundepen_flag
Definition: df1b2lap.cpp:31
static adpvm_manager * pvm_manager
Definition: fvar.hpp:8849
virtual void AD_uf_inner()
Definition: xmodelm4.cpp:40
static void set_yes_derivatives(void)
Definition: df1b2fun.h:761
df1b2_gradlist * f1b2gradlist
Definition: df1b2glo.cpp:49
void set_dependent_variable(const df1b2variable &_x)
Description not yet available.
Definition: df1b2fn2.cpp:699
static adpool * adpool_vector[]
Definition: df1b2fun.h:282
double crit
Definition: fvar.hpp:3184
int indexmin(void)
Definition: fvar.hpp:9354
static void set_no_derivatives(void)
Definition: df1b2fun.h:760
Description not yet available.
Definition: adrndeff.h:182
int indexmax() const
Definition: imatrix.h:142
i4_array * ptr3
Definition: adrndeff.h:136
i3_array * ptr2
Definition: adrndeff.h:135
void set_partition_sizes(int &num_der_blocks, ivector &minder, ivector &maxder, int nvariables)
Description not yet available.
Definition: df1b2lap.cpp:251
Description not yet available.
Definition: imatrix.h:69
static void set_NO_DERIVATIVES(void)
Disable accumulation of derivative information.
Definition: gradstrc.cpp:641
dcompressed_triplet * sparse_triplet
Definition: adrndeff.h:267
Description not yet available.
dvar_vector & shift(int min)
Description not yet available.
Definition: fvar_arr.cpp:28
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
static double cobjfun
Definition: df1b2fun.h:1360
void deallocate()
Deallocate dmatrix memory.
Definition: dmat.cpp:363
ivector * ptr1
Definition: adrndeff.h:152
double get_fx_fu(function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:2330
int indexmax(void)
Definition: fvar.hpp:9358
int indexmin() const
Definition: imatrix.h:138
static int no_stuff
Definition: df1b2lap.cpp:21
void begin_separable_call_stuff(void)
Description not yet available.
Definition: df1b2lap.cpp:2367
i4_array * ptr4
Definition: adrndeff.h:155
i3_array * get_ptr3(void)
Description not yet available.
Definition: df1b2lap.cpp:2507
dmatrix trans(const dmatrix &m1)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: dmat2.cpp:13
dvector eigenvalues(const banded_symmetric_dmatrix &_SS)
Description not yet available.
Definition: dmat28.cpp:411
#define x
int allocated(const ivector &v)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: fvar_a59.cpp:13
#define ADUNCONST(type, obj)
Creates a shallow copy of obj that is not CONST.
Definition: fvar.hpp:140
Vector of double precision numbers.
Definition: dvector.h:50
static int inner_opt_flag
Definition: admodel.h:1860
dvar_matrix * importance_sampling_components
Definition: adrndeff.h:224
static void set_minder(int n)
Definition: df1b2fun.h:299
static int separable_flag
Definition: df1b2fun.h:1351
static unsigned int nvar_vector[]
Definition: df1b2fun.h:288
static char ** argv
Definition: fvar.hpp:8866
static void set_blocksize(void)
Description not yet available.
Definition: df1b2fn2.cpp:89
static int passnumber
Definition: df1b2fun.h:289
void allocate(int ncl, int ncu)
Allocate memory for a dvector.
Definition: dvector.cpp:409
void fill_seqadd(double, double)
Fills dvector elements with values starting from base and incremented by offset.
Definition: cranfill.cpp:58
static re_objective_function_value * pobjfun
Definition: df1b2fun.h:1693
gauss_hermite_stuff * gh
Definition: adrndeff.h:223
void build_up_nested_shape(void)
Description not yet available.
Definition: df1b2lap.cpp:2391
exitptr ad_exit
Definition: gradstrc.cpp:53
df1_two_variable fabs(const df1_two_variable &x)
Definition: df12fun.cpp:891
gauss_hermite_stuff(laplace_approximation_calculator *lapprox, int use_gauss_hermite, int num_separable_calls, const ivector &itmp)
Description not yet available.
Definition: df1b2lap.cpp:2250
void allocate(int)
Description not yet available.
Definition: df1b2lap.cpp:2740
double calculate_laplace_approximation(const dvector &x, const dvector &u0, const dmatrix &Hess, const dvector &_xadjoint, const dvector &_uadjoint, const dmatrix &_Hessadjoint, function_minimizer *pmin)
Description not yet available.
Definition: df1b2lap.cpp:1476
Description not yet available.
Definition: df1b2fun.h:745
static dvariable reset(const dvar_vector &x)
Definition: model.cpp:345
imatrix * ptr2
Definition: adrndeff.h:153
ADMB variable vector.
Definition: fvar.hpp:2172
void trim(void)
Description not yet available.
Definition: df1b2lap.cpp:2548
laplace_approximation_calculator(int _xsize, int _usize, int _minder, int _maxder, function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:324
void initialize(void)
Description not yet available.
Definition: df1b2lap.cpp:2713
static int maxder
Definition: df1b2fun.h:292
Description not yet available.
Definition: df1b2fun.h:266
void gradcalc(int nvar, const dvector &g)
Definition: sgradclc.cpp:77
int atoi(adstring &s)
Returns a integer converted from input s.
Definition: atoi.cpp:20
static int nvarcalc()
Definition: model.cpp:152
void deallocate(void)
Description not yet available.
Definition: adpool.cpp:319
ivector sgn(const dvector &v)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: dvect24.cpp:11
~laplace_approximation_calculator()
Destructor.
Definition: df1b2lap.cpp:1140
static int no_ln_det_choleski_flag
Definition: fvar.hpp:8841
void set_u_dot(int i)
Description not yet available.
Definition: df1b2lp1.cpp:856
Description not yet available.
Definition: fvar.hpp:5769
static int stddev_vscale(const dvar_vector &d, const dvar_vector &x)
Definition: model.cpp:191
static void set_maxder(int n)
Definition: df1b2fun.h:300
dvariable ln_det_choleski_error(const dvar_matrix &, int &ierr)
Description not yet available.
Definition: fvar_m51.cpp:304
banded_symmetric_dmatrix * bHess
Definition: adrndeff.h:265
dcompressed_triplet * vsparse_triplet_adjoint
Definition: adrndeff.h:273
double evaluate_function(const dvector &x, function_minimizer *pfmin)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: df1b2lap.cpp:2026
void check_for_need_to_reallocate(int ip)
Does Nothing.
Definition: df1b2lap.cpp:1968
void set_default_hessian_type(void)
Description not yet available.
Definition: df1b2lap.cpp:2302
void limited_memory_quasi_newton(const independent_variables &, int)
Definition: lmnewton.cpp:32
static int print_hess_and_exit_flag
Definition: fvar.hpp:8838
dvar3_array * block_diagonal_vch
Definition: adrndeff.h:228
Description not yet available.
Definition: fvar.hpp:7951
dmatrix sort(const dmatrix &m, int column, int NSTACK)
Description not yet available.
Definition: dmsort.cpp:17
dvector get_uhat_quasi_newton(const dvector &x, function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:63
prnstream & endl(prnstream &)
dmatrix eigenvectors(const banded_symmetric_dmatrix &_SS, const dvector &_e)
Description not yet available.
Definition: dmat28.cpp:442
Description not yet available.
Definition: fvar.hpp:1937
Array of integers(int) with indexes from index_min to indexmax.
Definition: ivector.h:50
dvar_vector & get_x(void)
Definition: fvar.hpp:9327
void fill_randu(long int &n)
Fill matrix with random numbers.
Definition: ranfill.cpp:286
void allocate(int lib, int ub)
Definition: df1b2fn2.cpp:351
#define min(a, b)
Definition: cbivnorm.cpp:188
void initialize(void)
Description not yet available.
Definition: imat3.cpp:20
static int argc
Definition: fvar.hpp:8863
int indexmax() const
Get maximum valid index.
Definition: dvector.h:204
void normalized_gauss_hermite(const dvector &_x, const dvector &_w)
Gauss-Hermite quadature.
Definition: gaussher.cpp:343
static objective_function_value * pobjfun
Definition: admodel.h:2394
Description not yet available.
static void xinit(const dvector &x)
Definition: model.cpp:226
int indexmin() const
Definition: ivector.h:99
Description not yet available.
Definition: fvar.hpp:9345
int indexmax() const
Definition: ivector.h:104
Description not yet available.
Definition: fvar.hpp:2819
void initialize(void)
Initialze all elements of dvector to zero.
Definition: dvect5.cpp:10
dvector & shift(int min)
Shift valid range of subscripts.
Definition: dvector.cpp:52
void evaluate_function_gradient(double &f, const dvector &x, function_minimizer *pfmin, dvector &xadjoint, dvector &uadjoint)
Description not yet available.
Definition: df1b2lap.cpp:2166
d3_array * block_diagonal_vhessianadjoint
Definition: adrndeff.h:232
Description not yet available.
Definition: fvar.hpp:9293
double ln_det(const dmatrix &m1, int &sgn)
Compute log determinant of a constant matrix.
Definition: dmat3.cpp:536
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
std::ostream & get_output_stream()
Definition: adglobl.cpp:45
void get_newton_raphson_info(int xs, int us, const init_df1b2vector _y, dmatrix &Hess, dvector &grad, df1b2_gradlist *f1b2gradlist, function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:1924
static void set_nvar(unsigned int n)
Definition: df1b2fun.h:297
void initialize(void)
Description not yet available.
Definition: ivec2.cpp:17
dvar3_array * block_diagonal_vhessian
Definition: adrndeff.h:233
void allocate(const nested_calls_shape &nsc)
Description not yet available.
Definition: df1b2lap.cpp:2774
static void allocate(const unsigned int _n)
Definition: df1b2fn2.cpp:676
Description not yet available.
Definition: df1b2fun.h:373
ostream & operator<<(const ostream &_s, preshowpoint p)
Description not yet available.
Definition: admanip.cpp:48
double CHECK_HESSIAN_PENALTY
Definition: df1b2lap.cpp:34
static int have_bounded_random_effects
Definition: df1b2fun.h:1353
Description not yet available.
Definition: fvar.hpp:9410
Class definition of matrix with derivitive information .
Definition: fvar.hpp:2480
#define w
virtual void AD_uf_outer()
Definition: xmodelm4.cpp:39
~nested_calls_shape()
Description not yet available.
Definition: df1b2lap.cpp:2617
ivector * get_ptr1(void)
Description not yet available.
Definition: df1b2lap.cpp:2489
Description not yet available.
Definition: adrndeff.h:150
double ln_det_choleski(const banded_symmetric_dmatrix &MM, int &ierr)
Definition: dmat38.cpp:218
static int minder
Definition: df1b2fun.h:291
i4_array * get_ptr4(void)
Description not yet available.
Definition: df1b2lap.cpp:2516
double evaluate_function_quiet(const dvector &x, function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:2137
static int bandwidth
Definition: fvar.hpp:8837
dvar_compressed_triplet * vsparse_triplet
Definition: adrndeff.h:272
banded_symmetric_dmatrix * bHessadjoint
Definition: adrndeff.h:266
Description not yet available.
static int test_trust_flag
Definition: admodel.h:1856
Description not yet available.
Definition: admodel.h:1850
static void increment_adpool_counter(void)
Description not yet available.
Definition: df1b2glo.cpp:24
static void set_YES_DERIVATIVES(void)
Enable accumulation of derivative information.
Definition: gradstrc.cpp:650
void allocate(int nrl, int nrh, int ncl, int nch)
Allocates AD variable matrix with dimensions nrl to nrh by ncl to nch.
Definition: fvar_mat.cpp:216
virtual void user_function()
Definition: xmodelm4.cpp:43
void initialize(int sl, int sh, int nrl, const ivector &nrh, int ncl, const ivector &nch)
void get_second_ders(int xs, int us, const init_df1b2vector y, dmatrix &Hess, dmatrix &Dux, df1b2_gradlist *f1b2gradlist, function_minimizer *pfmin)
static void get_cHessian_contribution(dmatrix, int)
Description not yet available.
Definition: quadpri.cpp:591
static lmatrix * pointer_table
Definition: df1b2fun.h:1365
void random_effects_userfunction(double f, const dvector &x, const dvector &g)
void initialize()
Initializes elements of i4_array to zero.
Definition: i4arr.cpp:15
dvector value(const df1_one_vector &v)
Definition: df11fun.cpp:69
#define PI
Definition: fvar.hpp:95
Description not yet available.
Definition: fvar.hpp:3944
long iprint
Definition: fvar.hpp:3183
#define max(a, b)
Definition: cbivnorm.cpp:189
static int mc_phase
Definition: admodel.h:846
dvector get_uhat_lm_newton2(const dvector &x, function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:2814
void end_separable_call_stuff(void)
Description not yet available.
Definition: df1b2lap.cpp:2380
static int adpool_counter
Definition: df1b2fun.h:279
dvector get_uhat_lm_newton(const dvector &x, function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:218
void initialize(void)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: dmat7.cpp:12
static int get_num_quadratic_prior(void)
Definition: df1b2fun.h:1916
static double fun_without_pen
Definition: df1b2fun.h:1694
static double fun_without_pen
Definition: admodel.h:2395
nested_calls_shape nested_shape
Definition: adrndeff.h:187
double evaluate_function_no_derivatives(const dvector &x, function_minimizer *pfmin)
Description not yet available.
Definition: df1b2lap.cpp:2196
void allocate(const ad_integer &ncl, const index_type &ncu)
Allocate vector of integers with dimension [_ncl to _nch].
Definition: indextyp.cpp:488
void reset(void)
Description not yet available.
Definition: df1b2fn2.cpp:581
static int alternative_user_function_flag
Definition: adrndeff.h:199
Fundamental data type for reverse mode automatic differentiation.
Definition: fvar.hpp:1518
double * u_dot
Definition: df1b2fun.h:196
d3_array log(const d3_array &arr3)
Author: David Fournier Copyright (c) 2008-2012 Regents of the University of California.
Definition: d3arr2a.cpp:13
void df1b2_gradcalc1(void)
Description not yet available.
Definition: df1b2fun.cpp:145
static int print_importance_sampling_weights_flag
Definition: adrndeff.h:203
static void set_inactive_only_random_effects(void)
Definition: model.cpp:259
imatrix * get_ptr2(void)
Description not yet available.
Definition: df1b2lap.cpp:2498
d3_array pow(const d3_array &m, int e)
Description not yet available.
Definition: d3arr6.cpp:17
static int separable_calculation_type
Definition: df1b2fun.h:1352
static void get_M_calculations(void)
Description not yet available.
Definition: quadpri.cpp:30