ADMB Documentation  11.1.1020
 All Classes Files Functions Variables Typedefs Friends Defines
df1b2lp6.cpp
Go to the documentation of this file.
00001 /*
00002  * $Id: df1b2lp6.cpp 542 2012-07-10 21:04:06Z johnoel $
00003  *
00004  * Author: David Fournier
00005  * Copyright (c) 2008-2012 Regents of the University of California 
00006  */
00011 //#define USE_DD_STUFF
00012 //#define USE_DD
00013 #if defined(USE_LAPLACE)
00014 #  include <admodel.h>
00015 #  include <df1b2fun.h>
00016 #  include <adrndeff.h>
00017 double calculate_laplace_approximation(const dvector& x,const dvector& u0,
00018   const banded_symmetric_dmatrix& bHess,const dvector& _xadjoint,
00019   const dvector& _uadjoint,
00020   const banded_symmetric_dmatrix& _bHessadjoint,function_minimizer * pmin);
00021 double calculate_laplace_approximation(const dvector& x,const dvector& u0,
00022   const dmatrix& Hess,const dvector& _xadjoint,const dvector& _uadjoint,
00023   const dmatrix& _Hessadjoint,function_minimizer * pmin);
00024 
00025 double calculate_importance_sample_shess(const dvector& x,const dvector& u0,
00026   const dmatrix& Hess,const dvector& _xadjoint,const dvector& _uadjoint,
00027   const dmatrix& _Hessadjoint,function_minimizer * pmin);
00028 
00029 int use_dd_nr=0;
00030 int admb_ssflag=0;
00031 dvector solve(const dmatrix & st,const dmatrix & Hess,
00032   const dvector& grad);
00033 
00034 #if defined(USE_DD_STUFF)
00035 #  if defined(__MSVC32__)
00036     extern "C" _export  void dd_newton_raphson(int n,double * v,double * diag,
00037       double * ldiag, double * yy);
00038 #  else
00039     extern "C" void dd_newton_raphson(int n,double * v,double * diag,
00040       double * ldiag, double * yy);
00041 #  endif
00042 #endif
00043 
00048 void laplace_approximation_calculator::
00049   do_newton_raphson_banded(function_minimizer * pfmin,double f_from_1,
00050   int& no_converge_flag)
00051 {
00052   //quadratic_prior * tmpptr=quadratic_prior::ptr[0];
00053   //cout << tmpptr << endl;
00054 
00055 
00056   laplace_approximation_calculator::where_are_we_flag=2; 
00057   double maxg=fabs(evaluate_function(uhat,pfmin));
00058 
00059 
00060   laplace_approximation_calculator::where_are_we_flag=0; 
00061   dvector uhat_old(1,usize);
00062   for(int ii=1;ii<=num_nr_iters;ii++)
00063   {  
00064     // test newton raphson
00065     switch(hesstype)
00066     {
00067     case 3:
00068       bHess->initialize();
00069       break;
00070     case 4:
00071       Hess.initialize();
00072       break;
00073     default:
00074       cerr << "Illegal value for hesstype here" << endl;
00075       ad_exit(1);
00076     }
00077 
00078 
00079 
00080     grad.initialize();
00081     //int check=initial_params::stddev_scale(scale,uhat);
00082     //check=initial_params::stddev_curvscale(curv,uhat);
00083     //max_separable_g=0.0;
00084     sparse_count=0.0;
00085 
00086 
00087     step=get_newton_raphson_info_banded(pfmin);
00088     //if (bHess)
00089      // cout << "norm(*bHess) = " << norm(*bHess) << endl;
00090     //cout << "norm(Hess) = " << norm(Hess) << endl;
00091     //cout << grad << endl;
00092     //check_pool_depths();
00093     if (!initial_params::mc_phase)
00094       cout << "Newton raphson " << ii << "  ";
00095     if (quadratic_prior::get_num_quadratic_prior()>0)
00096     {
00097       quadratic_prior::get_cHessian_contribution(Hess,xsize);
00098       quadratic_prior::get_cgradient_contribution(grad,xsize);
00099     }
00100 
00101     int ierr=0;
00102     if (hesstype==3)
00103     {
00104       if (use_dd_nr==0)
00105       {
00106         banded_lower_triangular_dmatrix bltd=choleski_decomp(*bHess,ierr);
00107         if (ierr && no_converge_flag ==0)
00108         {
00109           no_converge_flag=1;
00110           //break;
00111         }
00112         if (ierr)
00113         {
00114           double oldval;
00115           evaluate_function(oldval,uhat,pfmin);
00116           uhat=banded_calculations_trust_region_approach(uhat,pfmin);
00117         }
00118         else
00119         {
00120           if (dd_nr_flag==0)
00121           {
00122             dvector v=solve(bltd,grad);
00123             step=-solve_trans(bltd,v);
00124             //uhat_old=uhat;
00125             uhat+=step;
00126           }
00127           else
00128           {
00129 #if defined(USE_DD_STUFF)
00130             int n=grad.indexmax();
00131             maxg=fabs(evaluate_function(uhat,pfmin));
00132             uhat=dd_newton_raphson2(grad,*bHess,uhat);
00133 #else
00134             cerr << "high precision Newton Raphson not implemented" << endl;
00135             ad_exit(1);
00136 #endif
00137           }
00138           maxg=fabs(evaluate_function(uhat,pfmin));
00139           if (f_from_1< pfmin->lapprox->fmc1.fbest)
00140           {
00141             uhat=banded_calculations_trust_region_approach(uhat,pfmin);
00142             maxg=fabs(evaluate_function(uhat,pfmin));
00143           }
00144         }
00145       }
00146       else
00147       {
00148         cout << "error not used" << endl;
00149         ad_exit(1);
00150        /*
00151         banded_symmetric_ddmatrix bHessdd=banded_symmetric_ddmatrix(*bHess);
00152         ddvector gradd=ddvector(grad);
00153         //banded_lower_triangular_ddmatrix bltdd=choleski_decomp(bHessdd,ierr);
00154         if (ierr && no_converge_flag ==0)
00155         {
00156           no_converge_flag=1;
00157           break;
00158         }
00159         if (ierr)
00160         {
00161           double oldval;
00162           evaluate_function(oldval,uhat,pfmin);
00163           uhat=banded_calculations_trust_region_approach(uhat,pfmin);
00164           maxg=fabs(evaluate_function(uhat,pfmin));
00165         }
00166         else
00167         {
00168           ddvector v=solve(bHessdd,gradd);
00169           step=-make_dvector(v);
00170           //uhat_old=uhat;
00171           uhat=make_dvector(ddvector(uhat)+step);
00172           maxg=fabs(evaluate_function(uhat,pfmin));
00173           if (f_from_1< pfmin->lapprox->fmc1.fbest)
00174           {
00175             uhat=banded_calculations_trust_region_approach(uhat,pfmin);
00176             maxg=fabs(evaluate_function(uhat,pfmin));
00177           }
00178         }
00179         */
00180       }
00181     
00182       if (maxg < 1.e-13) 
00183       {
00184         break;
00185       }
00186     } 
00187     else if (hesstype==4)
00188     {
00189       dvector step;
00190      
00191 #     if defined(USE_ATLAS)   
00192         if (!ad_comm::no_atlas_flag) 
00193         {
00194           step=-atlas_solve_spd(Hess,grad,ierr);
00195         }
00196         else
00197         {
00198           dmatrix A=choleski_decomp_positive(Hess,ierr);
00199           if (!ierr)
00200           {
00201             step=-solve(Hess,grad);
00202             //step=-solve(A*trans(A),grad);
00203           }
00204         }
00205         if (!ierr) break;
00206 #     else
00207         if (sparse_hessian_flag)
00208         {
00209           //step=-solve(*sparse_triplet,Hess,grad,*sparse_symbolic);
00210           dvector temp=solve(*sparse_triplet2,grad,*sparse_symbolic2,ierr);
00211           if (ierr)
00212           {
00213             step=-temp;
00214           }
00215           else
00216           {
00217             cerr << "matrix not pos definite in sparse choleski"  << endl;
00218             pfmin->bad_step_flag=1;
00219 
00220             int on;
00221             int nopt;
00222             if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-ieigvec",nopt))>-1)
00223             {
00224               dmatrix M=make_dmatrix(*sparse_triplet2);
00225 
00226               ofstream ofs3("inner-eigvectors");
00227               ofs3 << "eigenvalues and eigenvectors " << endl;
00228               dvector v=eigenvalues(M);
00229               dmatrix ev=trans(eigenvectors(M));
00230               ofs3 << "eigenvectors" << endl;
00231               int i;
00232               for (i=1;i<=ev.indexmax();i++)
00233                {
00234                   ofs3 << setw(4) << i  << " " 
00235                    << setshowpoint() << setw(14) << setprecision(10) << v(i) << " "
00236                    << setshowpoint() << setw(14) << setprecision(10)
00237                    << ev(i) << endl;
00238                }
00239             }
00240           }
00241           //cout << norm2(step-tmpstep) << endl;
00242           //dvector step1=-solve(Hess,grad);
00243           //cout << norm2(step-step1) << endl;
00244         }
00245         else
00246         {
00247           step=-solve(Hess,grad);
00248         }
00249 #     endif
00250       if (pmin->bad_step_flag)
00251         break;
00252       uhat_old=uhat;
00253       uhat+=step;
00254     
00255       double maxg_old=maxg;
00256       maxg=fabs(evaluate_function(uhat,pfmin));
00257       if (maxg>maxg_old)
00258       {
00259         uhat=uhat_old;
00260         evaluate_function(uhat,pfmin);
00261         break;
00262       }
00263       if (maxg < 1.e-13) 
00264       {
00265         break;
00266       }
00267     }
00268     
00269     if (sparse_hessian_flag==0)
00270     {
00271       for (int i=1;i<=usize;i++)
00272       {
00273         y(i+xsize)=uhat(i);
00274       }
00275     }
00276     else  
00277     {
00278       for (int i=1;i<=usize;i++)
00279       {
00280         value(y(i+xsize))=uhat(i);
00281       }
00282     }
00283   }
00284 }
00285 
00290 double laplace_approximation_calculator::
00291   inner_optimization_banded(/*dvector& uhat,*/ dvector& x,
00292   function_minimizer * pfmin,int& no_converge_flag)
00293 {
00294   int reset_flag=0;
00295   if (no_converge_flag) 
00296   {
00297     reset_flag=1;
00298     no_converge_flag=0;
00299   } 
00300   
00301   if (!inner_lmnflag)
00302   {
00303     if (!ADqd_flag)
00304     {
00305       uhat=get_uhat_quasi_newton(x,pfmin);
00306       double maxg=fabs(fmc1.gmax);
00307       if (maxg>1.0)
00308       {
00309         uhat=get_uhat_quasi_newton(x,pfmin);
00310       }
00311     }
00312     else
00313     {
00314       uhat=get_uhat_quasi_newton_qd(x,pfmin);
00315     }
00316   }
00317   else
00318   {
00319     uhat=get_uhat_lm_newton(x,pfmin);
00320     //uhat=get_uhat_lm_newton2(x,pfmin);
00321     //maxg=objective_function_value::gmax;
00322   }
00323   return fmc1.fbest;
00324 }
00325 
00330 dvector laplace_approximation_calculator::banded_calculations
00331   (const dvector& _x,const double& _f,function_minimizer * pfmin)
00332 {
00333   // for use when there is no separability
00334   ADUNCONST(dvector,x)
00335   ADUNCONST(double,f)
00336   //int i,j;
00337   int i;
00338 
00339   initial_params::set_inactive_only_random_effects(); 
00340   gradient_structure::set_NO_DERIVATIVES();
00341   initial_params::reset(x);    // get current x values into the model
00342   gradient_structure::set_YES_DERIVATIVES();
00343  
00344   initial_params::set_active_only_random_effects(); 
00345   if (init_switch==0)
00346   {
00347     gradient_structure::set_NO_DERIVATIVES();
00348     initial_params::xinit(ubest); 
00349     gradient_structure::set_YES_DERIVATIVES();
00350   }
00351   //int lmn_flag=0;
00352   //double maxg;
00353   //double maxg_save;
00354   //dvector uhat(1,usize);
00355   double f_from_1=0.0;
00356 
00357   int no_converge_flag=0;
00358  
00359   // this is the main loop to do inner optimization
00360   do
00361   {
00362     int icount=0;
00363     do
00364     {
00365       icount++;
00366       // do the inner optimziation
00367       if (inner_maxfn>0)
00368       {
00369         f_from_1=inner_optimization_banded(/*uhat,*/ x,pfmin,
00370           no_converge_flag);
00371       }
00372   
00373       if (sparse_hessian_flag==0)
00374       {
00375         for (i=1;i<=xsize;i++) { y(i)=x(i); }
00376         for (i=1;i<=usize;i++) { y(i+xsize)=uhat(i); }
00377       }
00378       else
00379       {
00380         for (i=1;i<=xsize;i++) { value(y(i))=x(i); }
00381         for (i=1;i<=usize;i++) { value(y(i+xsize))=uhat(i); }
00382       }
00383           
00384       laplace_approximation_calculator::where_are_we_flag=2; 
00385       if (admb_ssflag==0)
00386       {
00387         do_newton_raphson_banded(pfmin,f_from_1,no_converge_flag);
00388       }
00389       else
00390       {
00391         do_newton_raphson_state_space(pfmin,f_from_1,no_converge_flag);
00392       }
00393       laplace_approximation_calculator::where_are_we_flag=0; 
00394    
00395       if (num_nr_iters<=0) { evaluate_function(uhat,pfmin); }
00396   
00397       if (sparse_hessian_flag==0)
00398       {
00399         for (i=1;i<=usize;i++) { y(i+xsize)=uhat(i); }
00400       }
00401       else
00402       {
00403         for (i=1;i<=usize;i++) { value(y(i+xsize))=uhat(i); }
00404       }
00405       if (icount>2) pfmin->bad_step_flag=1;
00406       if (pfmin->bad_step_flag)
00407         return xadjoint;
00408     }
00409     while(no_converge_flag);
00410 
00411     /* If we are in mcmc phase we just need to calcualte the
00412        ln_det(Hess) and return
00413     */
00414     hs_symbolic & ssymb=*(pmin->lapprox->sparse_symbolic2);
00415     if (initial_params::mc_phase)
00416     {
00417       do_newton_raphson_banded(pfmin,f_from_1,no_converge_flag);
00418       int sgn=0;
00419       double& f = (double&) _f;
00420       f=initial_df1b2params::cobjfun;
00421      if (pmin->lapprox->sparse_hessian_flag==0)
00422      {
00423         if (bHess==0)
00424         {
00425           cerr << "Block diagonal Hessian is unallocated" << endl;
00426           ad_exit(1);
00427         }
00428         f+=0.5*ln_det_choleski(*bHess,sgn);
00429       }
00430       else
00431       {
00432         //hs_symbolic & ssymb=*(pmin->lapprox->sparse_symbolic2);
00433         //dvariable tmp=0.5*ln_det(*(pmin->lapprox->vsparse_triplet),
00434         //  ssymb,*(pmin->lapprox->sparse_triplet2));
00435         f+=0.5*ln_det(*(pmin->lapprox->sparse_triplet2),ssymb);
00436       }
00437     }
00438     else
00439     {
00440       xadjoint.initialize();
00441       uadjoint.initialize();
00442       Dux.initialize();
00443     
00444       if (hesstype==3)
00445         bHess->initialize();
00446       else if (hesstype==4)
00447         Hess.initialize();
00448     
00449       block_diagonal_flag=2;
00450       used_flags.initialize();
00451       funnel_init_var::lapprox=this;
00452       sparse_count=0.0;
00453     
00454       initial_params::straight_through_flag=1;
00455 
00456       if (sparse_triplet2)
00457         sparse_triplet2->initialize();
00458 
00459       pfmin->user_function();
00460       initial_params::straight_through_flag=0;
00461     
00462       int ierr=0;
00463 
00464       laplace_approximation_calculator::where_are_we_flag=3; 
00465       if (!ierr)
00466       {
00467         if (num_importance_samples==0)
00468         {
00469           if (hesstype==3)
00470           {
00471             f=calculate_laplace_approximation(x,uhat,*bHess,xadjoint,uadjoint,
00472               *bHessadjoint,pfmin);
00473           }
00474           else if (hesstype==4)
00475           {
00476             //cout << "Hess" << endl << Hess << endl;
00477             f=calculate_laplace_approximation(x,uhat,Hess,xadjoint,uadjoint,
00478               Hessadjoint,pfmin);
00479           }
00480           else
00481           {
00482             cerr << "Error in hesstype" << endl;
00483             ad_exit(1);
00484           }
00485         }
00486         else
00487         {
00488           if (hesstype==3)
00489           {
00490             //cerr << "Error not implemented yet" << endl;
00491             //ad_exit(1);
00492             f=calculate_importance_sample(x,uhat,*bHess,xadjoint,uadjoint,
00493               *bHessadjoint,pfmin);
00494           }
00495           else if (hesstype==4)
00496           {
00497             if (pmin->lapprox->sparse_hessian_flag==0)
00498               f=calculate_importance_sample(x,uhat,Hess,xadjoint,uadjoint,
00499                 Hessadjoint,pfmin);
00500             else
00501               f=calculate_importance_sample_shess(x,uhat,Hess,xadjoint,uadjoint,
00502                 Hessadjoint,pfmin);
00503           }
00504           else
00505           {
00506             cerr << "Error in hesstype" << endl;
00507             ad_exit(1);
00508           }
00509         }
00510       }
00511       else
00512       {
00513         f=1.e+30;
00514       }
00515 
00516       // set flag for thrid erivatvies and call function again because
00517       // stack is wiped out
00518       
00519     
00520       if (hesstype==3)
00521       {
00522         bHess->initialize();
00523       }
00524       else if (hesstype==4)
00525       {
00526         if (sparse_hessian_flag==0)
00527         {
00528           Hess.initialize();
00529         }
00530         else
00531         {
00532           sparse_triplet2->initialize();
00533         }
00534       }
00535       else
00536       {
00537         cerr << "Illegal value for hesstype" << endl;
00538         ad_exit(1);
00539       }
00540       initial_params::straight_through_flag=1;
00541       block_diagonal_flag=3;
00542       local_dtemp.initialize();
00543     
00544       // *****  Note for quadratic prior code: I don't think that this 
00545       // part gets added to the Hessian here.
00546       sparse_count=0;
00547       sparse_count_adjoint=0;
00548       pfmin->user_function();
00549     
00550       // *** Hessian calculated just above did not have quadratic prior
00551       // in it so can save this part for quadratci prioer adjoint calculations
00552       if (quadratic_prior::get_num_quadratic_prior()>0)
00553       {
00554         if (pHess_non_quadprior_part)
00555         {
00556           if (pHess_non_quadprior_part->indexmax() != Hess.indexmax())
00557           {
00558             delete pHess_non_quadprior_part;
00559             pHess_non_quadprior_part=0;
00560           }
00561         }
00562         if (!pHess_non_quadprior_part)
00563         {
00564           pHess_non_quadprior_part=new dmatrix(1,usize,1,usize);
00565           if (!pHess_non_quadprior_part)
00566           {
00567             cerr << "Error allocating memory for Hesssian part" << endl;
00568             ad_exit(1);
00569           }
00570         }
00571         (*pHess_non_quadprior_part)=Hess;
00572       } 
00573            
00574       block_diagonal_flag=0;
00575       initial_params::straight_through_flag=1;
00576     
00577       //dmatrix tHess=dmatrix(*bHess);
00578       initial_params::straight_through_flag=0;
00579       funnel_init_var::lapprox=0;
00580       //cout << initial_df1b2params::cobjfun << endl;
00581       //f=initial_df1b2params::cobjfun;
00582       block_diagonal_flag=0;
00583       dvector scale1(1,nvar);   // need to get scale from somewhere
00584       initial_params::set_inactive_only_random_effects(); 
00585       /*int check=*/initial_params::stddev_scale(scale1,x);
00586       //for (i=1;i<=xadjoint.indexmax();i++)
00587       //  xadjoint(i)*=scale1(i);
00588       laplace_approximation_calculator::where_are_we_flag=0; 
00589   
00590       if (df1b2quadratic_prior::get_num_quadratic_prior()>0)
00591       {  
00592        // !!!! need to fix this!!!!!!!!!!!!!!!!!!!!!!!
00593         laplace_approximation_calculator::where_are_we_flag=3; 
00594         quadratic_prior::in_qp_calculations=1;
00595         funnel_init_var::lapprox=this;
00596         df1b2_gradlist::set_no_derivatives(); 
00597         df1b2quadratic_prior::get_Lxu_contribution(Dux);
00598         quadratic_prior::in_qp_calculations=0;
00599         funnel_init_var::lapprox=0;
00600         laplace_approximation_calculator::where_are_we_flag=0; 
00601       }   
00602       if (initial_df1b2params::separable_flag)
00603       {
00604         dvector scale(1,nvar);   // need to get scale from somewhere
00605         /*int check=*/initial_params::stddev_scale(scale,x);
00606         dvector sscale=scale(1,Dux(1).indexmax());
00607         for (i=1;i<=usize;i++)
00608         {
00609           Dux(i)=elem_prod(Dux(i),sscale);
00610         }
00611         local_dtemp=elem_prod(local_dtemp,sscale);
00612       }
00613       //cout << trans(Dux)(1) << endl;
00614       //cout << trans(Dux)(3) << endl;  
00615       if (quadratic_prior::get_num_quadratic_prior()>0)
00616       {
00617         dvector tmp=evaluate_function_with_quadprior(x,usize,pfmin);
00618         local_dtemp+=tmp;
00619       }
00620           
00621       for (i=1;i<=xsize;i++)
00622       {
00623         xadjoint(i)+=local_dtemp(i);
00624       }    
00625       if (df1b2quadratic_prior::get_num_quadratic_prior()>0)
00626       {
00627        // !!!! need to fix this!!!!!!!!!!!!!!!!!!!!!!!
00628         quadratic_prior::get_cHessian_contribution_from_vHessian(Hess,xsize);
00629       }
00630    
00631       if (hesstype==3)
00632       {
00633         //xadjoint -= uadjoint*solve(*bHess,Dux);
00634         if (bHess_pd_flag==0)
00635         {
00636           xadjoint -= solve(*bHess,uadjoint)*Dux;
00637         }
00638       }
00639       else if (hesstype==4)
00640       {
00641         if (sparse_hessian_flag)
00642         {
00643           //dvector tmp=solve(*sparse_triplet,Hess,uadjoint,*sparse_symbolic)*Dux;
00644           dvector tmp=solve(*sparse_triplet2,uadjoint,*sparse_symbolic2)*Dux;
00645           xadjoint -= tmp;
00646         }
00647         else
00648         {
00649           xadjoint -= solve(Hess,uadjoint)*Dux;
00650         }
00651       }
00652     }
00653     if (bHess_pd_flag==0) break;
00654   }
00655   while(1);
00656 
00657   return xadjoint;
00658 }
00659   //int check_pool_flag1=0;
00660 
00665 void laplace_approximation_calculator::
00666   do_separable_stuff_newton_raphson_banded(df1b2variable& ff)
00667 {
00668   //***********************************************************
00669   //***********************************************************
00670   set_dependent_variable(ff);
00671   df1b2_gradlist::set_no_derivatives();
00672   df1b2variable::passnumber=1;
00673   //if (check_pool_flag1)
00674    // check_pool_depths();
00675   df1b2_gradcalc1();
00676   //if (check_pool_flag1)
00677   //  check_pool_depths();
00678    
00679   init_df1b2vector & locy= *funnel_init_var::py;
00680   imatrix& list=*funnel_init_var::plist;
00681   int num_local_re=0;
00682   int num_fixed_effects=0;
00683 
00684   int i;
00685   ivector lre_index(1,funnel_init_var::num_active_parameters);
00686   ivector lfe_index(1,funnel_init_var::num_active_parameters);
00687 
00688   for (i=1;i<=funnel_init_var::num_active_parameters;i++)
00689   {
00690     if (list(i,1)>xsize) 
00691     {
00692       lre_index(++num_local_re)=i;
00693     }
00694     else if (list(i,1)>0) 
00695     {
00696       lfe_index(++num_fixed_effects)=i;
00697     }
00698   }
00699   
00700   if (num_local_re > 0)
00701   {
00702     int j;
00703     switch(hesstype)
00704     {
00705     case 3:
00706       for (i=1;i<=num_local_re;i++)
00707       {
00708         int lrei=lre_index(i);
00709         for (j=1;j<=num_local_re;j++)
00710         {
00711           int lrej=lre_index(j);
00712           int i1=list(lrei,1)-xsize;
00713           int i2=list(lrei,2);
00714           int j1=list(lrej,1)-xsize;
00715           int j2=list(lrej,2);
00716           if (i1>=j1) (*bHess)(i1,j1)+=locy(i2).u_bar[j2-1];
00717         }
00718       }
00719       break;
00720     case 4:
00721       if (sparse_hessian_flag==0)
00722       {
00723         for (i=1;i<=num_local_re;i++)
00724         {
00725           int lrei=lre_index(i);
00726           for (j=1;j<=num_local_re;j++)
00727           {
00728             int lrej=lre_index(j);
00729             int i1=list(lrei,1)-xsize;
00730             int i2=list(lrei,2);
00731             int j1=list(lrej,1)-xsize;
00732             int j2=list(lrej,2);
00733             Hess(i1,j1)+=locy(i2).u_bar[j2-1];
00734           }
00735         }
00736       }
00737       else
00738       {
00739         for (i=1;i<=num_local_re;i++)
00740         {
00741           int lrei=lre_index(i);
00742           for (j=1;j<=num_local_re;j++)
00743           {
00744             int lrej=lre_index(j);
00745             int i1=list(lrei,1)-xsize;
00746             int i2=list(lrei,2);
00747             int j1=list(lrej,1)-xsize;
00748             int j2=list(lrej,2);
00749             
00750             if (i1<=j1)
00751             {
00752               sparse_count++;
00753               (*sparse_triplet2)((*sparse_iterator)(sparse_count))
00754                 +=locy(i2).u_bar[j2-1];
00755             }
00756           }
00757         }
00758       }
00759       break;
00760     default:
00761       cerr << "illegal value for hesstype" << endl;
00762       ad_exit(1);
00763     }
00764   
00765     for (i=1;i<=num_local_re;i++)
00766     {
00767       int lrei=lre_index(i);
00768       int i1=list(lrei,1);
00769       int i2=list(lrei,2);
00770       //grad(i1-xsize)= re_objective_function_value::pobjfun->u_dot[i2-1];
00771       grad(i1-xsize)+=ff.u_dot[i2-1];
00772     }
00773   } 
00774 
00775   f1b2gradlist->reset();
00776   f1b2gradlist->list.initialize();
00777   f1b2gradlist->list2.initialize();
00778   f1b2gradlist->list3.initialize();
00779   f1b2gradlist->nlist.initialize();
00780   f1b2gradlist->nlist2.initialize();
00781   f1b2gradlist->nlist3.initialize();
00782   funnel_init_var::num_vars=0; 
00783   funnel_init_var::num_active_parameters=0; 
00784   funnel_init_var::num_inactive_vars=0; 
00785 }
00786 int tmp_testcount=0;
00787 df1b2variable * tmp_pen=00;
00788 
00793 dvector laplace_approximation_calculator::
00794   get_newton_raphson_info_banded (function_minimizer * pfmin)
00795 {
00796   //int i,j,ip; 
00797   int ip; 
00798   
00799   int nv=initial_df1b2params::set_index();
00800   if (allocated(used_flags))
00801   {
00802     if (used_flags.indexmax() != nv)
00803     {
00804       used_flags.safe_deallocate();
00805     }
00806   }
00807   if (!allocated(used_flags))
00808   {
00809     used_flags.safe_allocate(1,nv);
00810   }
00811   
00812   for (ip=1;ip<=num_der_blocks;ip++)
00813   {
00814     if (ip>1)   // change to combine sparse matrix stuff with num der blocks
00815     {           // df  3-4-09
00816       sparse_count=0;
00817     }
00818     used_flags.initialize();
00819     // do we need to reallocate memory for df1b2variables?
00820     check_for_need_to_reallocate(ip);
00821     df1b2_gradlist::set_no_derivatives();
00822     //cout << re_objective_function_value::pobjfun << endl;
00823     //cout << re_objective_function_value::pobjfun->ptr << endl;
00824     (*re_objective_function_value::pobjfun)=0;
00825     df1b2variable pen=0.0;
00826     tmp_pen=&pen;
00827   tmp_testcount++;
00828     df1b2variable zz=0.0;
00829 
00830 
00831     initial_df1b2params::reset(y,pen);
00832 
00833 
00834     // call function to do block diagonal newton-raphson
00835     // the step vector from the newton-raphson is in the vector step
00836     df1b2_gradlist::set_yes_derivatives();
00837     
00838     funnel_init_var::lapprox=this;
00839     //cout << funnel_init_var::lapprox << endl;
00840     block_diagonal_flag=1;
00841    /*
00842     if (tmp_testcount>=9)
00843     {
00844       pen.deallocate();
00845     }
00846     */
00847 
00848     if (ip==1)
00849     {
00850       if (sparse_triplet2)
00851         sparse_triplet2->initialize();
00852     }
00853 
00854     pfmin->user_function();
00855     /*
00856     if (tmp_testcount>=9)
00857     {
00858       pen.deallocate();
00859     }
00860      */
00861     funnel_init_var::lapprox=0;
00862     block_diagonal_flag=0;
00863     pen.deallocate();
00864   }
00865 
00866   return step;
00867 }
00868 
00873 void laplace_approximation_calculator::
00874   do_separable_stuff_laplace_approximation_banded(df1b2variable& ff)
00875 {
00876   set_dependent_variable(ff);
00877   //df1b2_gradlist::set_no_derivatives();
00878   df1b2variable::passnumber=1;
00879   df1b2_gradcalc1();
00880    
00881   init_df1b2vector & locy= *funnel_init_var::py;
00882   imatrix& list=*funnel_init_var::plist;
00883 
00884   int i; int j; int us=0; int xs=0;
00885   ivector lre_index(1,funnel_init_var::num_active_parameters);
00886   ivector lfe_index(1,funnel_init_var::num_active_parameters);
00887 
00888   for (i=1;i<=funnel_init_var::num_active_parameters;i++)
00889   {
00890     if (list(i,1)>xsize) 
00891     {
00892       lre_index(++us)=i;
00893     }
00894     else if (list(i,1)>0) 
00895     {
00896       lfe_index(++xs)=i;
00897     }
00898   }
00899   
00900   for (j=1;j<=xs;j++)
00901   {
00902     int j1=list(lfe_index(j),1);
00903     int j2=list(lfe_index(j),2);
00904     xadjoint(j1)+=ff.u_dot[j2-1];
00905   }
00906   
00907   if (us>0)
00908   {
00909     if (hesstype==3)
00910     {
00911       for (i=1;i<=us;i++)
00912       {
00913         for (j=1;j<=us;j++)
00914         {
00915           int i1=list(lre_index(i),1)-xsize;
00916           int i2=list(lre_index(i),2);
00917           int j1=list(lre_index(j),1)-xsize;
00918           int j2=list(lre_index(j),2);
00919           if (i1>=j1) (*bHess)(i1,j1)+=locy(i2).u_bar[j2-1];
00920         }
00921       }
00922     }
00923     else if (hesstype==4)
00924     {
00925       if (sparse_hessian_flag==0)
00926       {
00927         for (i=1;i<=us;i++)
00928         {
00929           for (j=1;j<=us;j++)
00930           {
00931             int i1=list(lre_index(i),1)-xsize;
00932             int i2=list(lre_index(i),2);
00933             int j1=list(lre_index(j),1)-xsize;
00934             int j2=list(lre_index(j),2);
00935             Hess(i1,j1)+=locy(i2).u_bar[j2-1];
00936           }
00937         }
00938       }
00939       else
00940       {
00941         for (i=1;i<=us;i++)
00942         {
00943           for (j=1;j<=us;j++)
00944           {
00945             int i1=list(lre_index(i),1)-xsize;
00946             int i2=list(lre_index(i),2);
00947             int j1=list(lre_index(j),1)-xsize;
00948             int j2=list(lre_index(j),2);
00949             
00950             if (i1<=j1)
00951             {
00952               sparse_count++;
00953               (*sparse_triplet2)((*sparse_iterator)(sparse_count))
00954                 +=locy(i2).u_bar[j2-1];
00955             }
00956           }
00957         }
00958       }
00959     }
00960 
00961     for (i=1;i<=us;i++)
00962     {
00963       int i1=list(lre_index(i),1)-xsize;
00964       int i2=list(lre_index(i),2);
00965       uadjoint(i1)+=ff.u_dot[i2-1];
00966     }
00967   
00968     for (i=1;i<=us;i++)
00969     {
00970       for (j=1;j<=xs;j++)
00971       {
00972         int i1=list(lre_index(i),1)-xsize;
00973         int i2=list(lre_index(i),2);
00974         int j1=list(lfe_index(j),1);
00975         int j2=list(lfe_index(j),2);
00976         Dux(i1,j1)+=locy(i2).u_bar[j2-1];
00977       }
00978     }
00979   } 
00980   f1b2gradlist->reset();
00981   f1b2gradlist->list.initialize();
00982   f1b2gradlist->list2.initialize();
00983   f1b2gradlist->list3.initialize();
00984   f1b2gradlist->nlist.initialize();
00985   f1b2gradlist->nlist2.initialize();
00986   f1b2gradlist->nlist3.initialize();
00987   funnel_init_var::num_vars=0;
00988   funnel_init_var::num_active_parameters=0; 
00989   funnel_init_var::num_inactive_vars=0;
00990 }
00991 
00996 double calculate_laplace_approximation(const dvector& x,const dvector& u0,
00997   const banded_symmetric_dmatrix& bHess,const dvector& _xadjoint,
00998   const dvector& _uadjoint,
00999   const banded_symmetric_dmatrix& _bHessadjoint,function_minimizer * pmin)
01000 {
01001   ADUNCONST(dvector,xadjoint)
01002   ADUNCONST(dvector,uadjoint)
01003   ADUNCONST(banded_symmetric_dmatrix,bHessadjoint)
01004   int bw=bHess.bandwidth();
01005   const int xs=x.size();
01006   const int us=u0.size();
01007   gradient_structure::set_YES_DERIVATIVES();
01008   int nvar=x.size()+u0.size()+((bw+1)*(2*u0.size()-bw))/2;
01009   independent_variables y(1,nvar);
01010   
01011   // need to set random effects active together with whatever
01012   // init parameters should be active in this phase
01013   initial_params::set_inactive_only_random_effects(); 
01014   initial_params::set_active_random_effects(); 
01015   /*int onvar=*/initial_params::nvarcalc(); 
01016   initial_params::xinit(y);    // get the initial values into the
01017   y(1,xs)=x;
01018 
01019   int i,j;
01020 
01021  // Here need hooks for sparse matrix structures
01022   int ii=xs+us+1;
01023   for (i=1;i<=us;i++)
01024   {
01025     int jmin=admax(1,i-bw+1);
01026     for (j=jmin;j<=i;j++)
01027     y(ii++)=bHess(i,j);
01028   }
01029 
01030   dvar_vector vy=dvar_vector(y); 
01031   banded_symmetric_dvar_matrix vHess(1,us,bw);
01032   initial_params::reset(vy);    // get the initial values into the
01033   ii=xs+us+1;
01034   for (i=1;i<=us;i++)
01035   {
01036     int jmin=admax(1,i-bw+1);
01037     for (j=jmin;j<=i;j++)
01038       vHess(i,j)=vy(ii++);
01039   }
01040 
01041    dvariable vf=0.0;
01042 
01043    
01044    *objective_function_value::pobjfun=0.0;
01045    pmin->AD_uf_outer();
01046    vf+=*objective_function_value::pobjfun;
01047 
01048    int sgn=0;
01049    dvariable ld;
01050   
01051    int eigswitch=0;
01052    if (eigswitch)
01053    {
01054      ofstream ofs("ee");
01055      dvector ev(bHess.indexmin(),bHess.indexmax());
01056      dmatrix evecs=eigenvectors(dmatrix(bHess),ev);
01057      ofs << setprecision(3) << setw(12) << setscientific() << dmatrix(bHess) << endl << endl;
01058      ofs << ev << endl << endl << evecs << endl;
01059    }
01060    ld=0.5*ln_det_choleski(vHess,sgn);
01061    if (sgn==1)
01062    {
01063      cout << "Choleski failed" << endl;
01064      pmin->lapprox->bHess_pd_flag=1;
01065    }
01066 
01067    vf+=ld;
01068    const double ltp=0.5*log(2.0*PI);
01069    vf-=us*ltp;
01070    double f=value(vf);
01071    dvector g(1,nvar);
01072    gradcalc(nvar,g);
01073   
01074   ii=1;
01075   for (i=1;i<=xs;i++)
01076     xadjoint(i)=g(ii++);
01077   for (i=1;i<=us;i++)
01078     uadjoint(i)=g(ii++);
01079   for (i=1;i<=us;i++)
01080   {
01081     int jmin=admax(1,i-bw+1);
01082     for (j=jmin;j<=i;j++)
01083       bHessadjoint(i,j)=g(ii++);
01084   }
01085   return f;
01086 }
01087 
01092 dvector laplace_approximation_calculator::
01093   banded_calculations_trust_region_approach(const dvector& _uhat,
01094   function_minimizer * pfmin)
01095 {
01096   dvector& uhat=(dvector&) _uhat;
01097   dvector uhat_old(uhat.indexmin(),uhat.indexmax());
01098   dvector uhat_new(uhat.indexmin(),uhat.indexmax());
01099   dvector uhat_best(uhat.indexmin(),uhat.indexmax());
01100  
01101   double wght=0.0;
01102   double delta=5.e-5;
01103   //do
01104   dvector values(1,300);
01105   double oldfbest=pmin->lapprox->fmc1.fbest; 
01106   double newfbest;
01107   int have_value=0;
01108   //for (int jj=1;jj<=300;jj++)
01109   int jj=1;
01110   double lastval=oldfbest;
01111   do 
01112   {
01113     jj++;
01114     wght+=delta;
01115     //double wght=0.0;
01116     double newval=0.0;
01117     //cout << "Enter weight size " << endl;
01118     //cin >> wght;
01119     if (wght<0.0) 
01120       break;
01121     int mmin=bHess->indexmin();
01122     int mmax=bHess->indexmax();
01123     banded_symmetric_dmatrix tmp(mmin,mmax,bHess->bandwidth());
01124     tmp=*bHess;
01125     uhat_old=uhat;
01126     int ierr=0;
01127     for (int i=mmin;i<=mmax;i++)
01128     {
01129       tmp(i,i)+=wght;
01130     }
01131     banded_lower_triangular_dmatrix bltd=choleski_decomp(tmp,ierr);
01132     if (!ierr)
01133     {
01134       dvector v=solve(bltd,grad);
01135       step=-solve_trans(bltd,v);
01136   
01137       uhat_old=uhat;
01138       uhat+=step;
01139       //cout << "norm(uhat_old) = " << norm(uhat_old) 
01140        //    << "   norm(uhat) = " << norm(uhat)  << endl;
01141        
01142       /*double maxg=*/fabs(evaluate_function(newval,uhat,pfmin));
01143       if (have_value && newval>newfbest) 
01144       {
01145         break;
01146       }
01147       if (jj>1)
01148       { 
01149         if (newval<lastval)  // we are doing better so increasse step size
01150         {
01151           delta*=2;
01152         }
01153         if (newval>lastval && !have_value)  // we have gone to far go back
01154         {
01155           wght-=delta;
01156           delta/=16;
01157         }
01158       }
01159       lastval=newval;
01160         
01161       if (newval<newfbest)
01162       {
01163         newfbest=newval;
01164         uhat_best=uhat;     
01165         have_value=jj;
01166       }
01167       uhat_new=uhat;
01168       uhat=uhat_old;
01169     }
01170     else
01171     {
01172       delta*=2;
01173     }
01174   }
01175   while(jj<10);
01176   if (!have_value)
01177   {
01178     cerr << "can't improve function value in trust region calculations"
01179          << endl;
01180     //ad_exit(1);
01181   }
01182   return uhat_best;
01183   initial_params::set_active_only_random_effects(); 
01184   //int lmn_flag=0;
01185   double maxg;
01186   //double maxg_save;
01187   double f_from_1=0.0;
01188   if (!inner_lmnflag)
01189   {
01190     if (!ADqd_flag)
01191     {
01192       uhat=get_uhat_quasi_newton(uhat_new,pfmin);
01193       maxg=fabs(fmc1.gmax);
01194       f_from_1=fmc1.fbest;
01195     }
01196     else
01197     {
01198       uhat=get_uhat_quasi_newton_qd(uhat_new,pfmin);
01199     }
01200   }
01201   else
01202   {
01203     uhat=get_uhat_lm_newton(uhat_new,pfmin);
01204   }
01205   return uhat;
01206 }
01207 #endif  //#if defined(USE_LAPLACE)