ADMB Documentation  11.1.1016
 All Classes Files Functions Variables Typedefs Friends Defines
df1b2lp2.cpp
Go to the documentation of this file.
00001 /*
00002  * $Id: df1b2lp2.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 #if defined(USE_LAPLACE)
00012 #  include <admodel.h>
00013 #  include <df1b2fun.h>
00014 #  include <adrndeff.h>
00015 void get_second_ders(int xs,int us,const init_df1b2vector y,dmatrix& Hess,
00016   dmatrix& Dux, df1b2_gradlist * f1b2gradlist,function_minimizer * pfmin,
00017   laplace_approximation_calculator* lap);
00018 double calculate_laplace_approximation(const dvector& x,const dvector& u0,
00019   const dmatrix& Hess,const dvector& _xadjoint,const dvector& _uadjoint,
00020   const dmatrix& _Hessadjoint,function_minimizer * pmin);
00021 
00022 static void xxx(ivector re_list,ivector fe_list){}
00023 
00028 dvector laplace_approximation_calculator::block_diagonal_calculations
00029   (const dvector& _x,const double& _f,function_minimizer * pfmin)
00030 {
00031   // for use when there is no separability
00032   ADUNCONST(dvector,x)
00033   ADUNCONST(double,f)
00034   //int i,j;
00035   int i;
00036 
00037   initial_params::set_inactive_only_random_effects(); 
00038   gradient_structure::set_NO_DERIVATIVES();
00039   initial_params::reset(x);    // get current x values into the model
00040   gradient_structure::set_YES_DERIVATIVES();
00041 
00042   initial_params::set_active_only_random_effects(); 
00043   //int lmn_flag=0;
00044   if (!inner_lmnflag)
00045   {
00046     if (!ADqd_flag)
00047       uhat=get_uhat_quasi_newton_block_diagonal(x,pfmin);
00048     else
00049       uhat=get_uhat_quasi_newton_qd(x,pfmin);
00050   }
00051   else
00052   {
00053     uhat=get_uhat_lm_newton(x,pfmin);
00054   }
00055   if (!allocated(scale))
00056   {
00057     scale.allocate(1,uhat.indexmax());
00058   }
00059   else
00060   {
00061     if (scale.indexmax() != uhat.indexmax())
00062     {
00063       scale.deallocate();
00064       scale.allocate(1,uhat.indexmax());
00065     }
00066   }
00067 
00068   if (!allocated(curv))
00069   {
00070     curv.allocate(1,uhat.indexmax());
00071   }
00072   else
00073   {
00074     if (curv.indexmax() != uhat.indexmax())
00075     {
00076       curv.deallocate();
00077       curv.allocate(1,uhat.indexmax());
00078     }
00079   }
00080 
00081   if (sparse_hessian_flag==0)
00082   {
00083     for (i=1;i<=xsize;i++)
00084     {  
00085       y(i)=x(i);
00086     }
00087     for (int i=1;i<=usize;i++)
00088     {
00089       y(i+xsize)=uhat(i);
00090     }
00091   }
00092   else  
00093   {
00094     for (i=1;i<=xsize;i++)
00095     {  
00096       value(y(i))=x(i);
00097     }
00098     for (int i=1;i<=usize;i++)
00099     {
00100       value(y(i+xsize))=uhat(i);
00101     }
00102   }
00103   //cout << y << endl;
00104         
00105   for(int ii=1;ii<=num_nr_iters;ii++)
00106   {  
00107     {   
00108       // test newton raphson
00109       //Hess.initialize();
00110       int check=initial_params::stddev_scale(scale,uhat);
00111       check=initial_params::stddev_curvscale(curv,uhat);
00112       max_separable_g=0.0;
00113       pmin->inner_opt_flag=1;
00114       step=get_newton_raphson_info_block_diagonal(pfmin);
00115       cout << "max separable g " << max_separable_g << endl; 
00116       cout << "Newton raphson " << ii << endl;
00117       uhat+=step;
00118     
00119       evaluate_function(uhat,pfmin);
00120       pmin->inner_opt_flag=0;
00121     }
00122   
00123     if (sparse_hessian_flag==0)
00124     {
00125       for (int i=1;i<=usize;i++)
00126       {
00127         y(i+xsize)=uhat(i);
00128       }
00129     }
00130     else  
00131     {
00132       for (int i=1;i<=usize;i++)
00133       {
00134         value(y(i+xsize))=uhat(i);
00135       }
00136     }
00137   }
00138  
00139   cout << initial_df1b2params::cobjfun << endl;
00140   xadjoint.initialize();
00141   uadjoint.initialize();
00142   block_diagonal_flag=2;
00143   used_flags.initialize();
00144   funnel_init_var::lapprox=this;
00145   if (use_gauss_hermite>0)
00146   {
00147     df1b2variable pen=0.0;
00148     initial_df1b2params::reset(y,pen);
00149     initial_params::straight_through_flag=0;
00150     block_diagonal_flag=6;
00151     num_separable_calls=0;
00152     // get the block diagonal hessians to use in importance sampling
00153     pfmin->user_function();
00154     //cout << (*pfmin->lapprox->block_diagonal_hessian) << endl;
00155     block_diagonal_flag=2;
00156     initial_params::straight_through_flag=0;
00157  
00158     // do importance sampling and get ders bakc to Hessian adjoint
00159     // new stuff for more than one random effect in each separable call
00160     //  Apr 17 07
00161     if (multi_random_effects==0)
00162     {
00163       f=do_gauss_hermite_block_diagonal(x,uhat,Hess,xadjoint,
00164         uadjoint,Hessadjoint,pfmin);
00165     }
00166     else
00167     {
00168       f=do_gauss_hermite_block_diagonal_multi(x,uhat,Hess,xadjoint,
00169         uadjoint,Hessadjoint,pfmin);
00170     }
00171     int xmax=xadjoint.indexmax();
00172     dvector x_con(1,xmax);
00173     x_con.initialize();
00174     dvector dscale(1,nvar);   // need to get scale from somewhere
00175     dscale=1.0;
00176     /*int check=*/initial_params::stddev_scale(dscale,x);
00177     dvector sscale=dscale(1,xsize);
00178     // *******************************************************
00179     // *******************************************************
00180     // *******************************************************
00181     // derivatives from hessian adjoint back
00182     {
00183       x_con.initialize();
00184       
00185       for (int i=1;i<=num_separable_calls;i++)
00186       {
00187         ivector& re_list=(*block_diagonal_re_list)(i);
00188         ivector& fe_list=(*block_diagonal_fe_list)(i);
00189         dmatrix& Dux=(*block_diagonal_Dux)(i);
00190         dmatrix& H=(*block_diagonal_hessian)(i);
00191         xxx(re_list,fe_list);
00192         int mmax=re_list.indexmax();
00193         dvector tmp(1,mmax);
00194         
00195         int j;
00196         for (j=1;j<=re_list.indexmax();j++)
00197         {
00198           tmp(j)=uadjoint(re_list(j)-xmax);
00199         } 
00200 
00201         if (allocated(fe_list))
00202         {
00203           if (allocated(H))
00204           {
00205             dvector tmp1=solve(H,tmp);
00206             dvector xtmp=tmp1*Dux;
00207             for (j=1;j<=fe_list.indexmax();j++)
00208             {
00209               x_con(fe_list(j))+=xtmp(j);
00210             }
00211           }
00212         }
00213       }
00214       if (initial_df1b2params::separable_flag)
00215       {
00216         x_con=elem_prod(x_con,sscale);
00217       }
00218     }
00219     xadjoint-=x_con;
00220     // *******************************************************
00221     // *******************************************************
00222     // *******************************************************
00223     
00224     block_diagonal_flag=3;
00225     //pfmin->lapprox->xadjoint.initialize();
00226     //pfmin->lapprox->uadjoint.initialize();
00227     pfmin->lapprox->num_separable_calls=0;
00228     pfmin->lapprox->check_local_xadjoint.initialize();
00229     pfmin->lapprox->check_local_xadjoint2.initialize();
00230     pfmin->lapprox->check_local_uadjoint.initialize();
00231     pfmin->lapprox->check_local_uadjoint2.initialize();
00232     //df1b2_gradlist::set_yes_derivatives();
00233     //initial_df1b2params::reset(y,pen);
00234     pfmin->user_function();
00235     dvector lcx=elem_prod(check_local_xadjoint2,sscale);
00236     xadjoint+=lcx;
00237     //df1b2_gradlist::set_no_derivatives();
00238     funnel_init_var::lapprox=0;
00239     block_diagonal_flag=0;
00240     initial_params::set_inactive_only_random_effects(); 
00241   }
00242   else if (num_importance_samples>0)
00243   {
00244     df1b2variable pen=0.0;
00245     initial_df1b2params::reset(y,pen);
00246     initial_params::straight_through_flag=0;
00247     block_diagonal_flag=6;
00248     num_separable_calls=0;
00249     // get the block diagonal hessians to use in importance sampling
00250     pfmin->user_function();
00251     //cout << (*pfmin->lapprox->block_diagonal_hessian) << endl;
00252     block_diagonal_flag=2;
00253     initial_params::straight_through_flag=0;
00254     // do importance sampling and get ders bakc to Hessian adjoint
00255     if (isfunnel_flag==0)
00256     {
00257       if (antiflag==0)
00258       {
00259         f=calculate_importance_sample_block_diagonal_option2(x,uhat,Hess,
00260           xadjoint,uadjoint,Hessadjoint,pfmin);
00261       }
00262       else
00263       {
00264         f=calculate_importance_sample_block_diagonal_option_antithetical(x,
00265           uhat,Hess,xadjoint,uadjoint,Hessadjoint,pfmin);
00266       }
00267     }
00268     else
00269     {
00270       f=calculate_importance_sample_block_diagonal_funnel(x,uhat,Hess,xadjoint,
00271         uadjoint,Hessadjoint,pfmin);
00272     }
00273 
00274     int xmax=xadjoint.indexmax();
00275     dvector x_con(1,xmax);
00276     x_con.initialize();
00277     dvector dscale(1,nvar);   // need to get scale from somewhere
00278     dscale=1.0;
00279     /*int check=*/initial_params::stddev_scale(dscale,x);
00280     dvector sscale=dscale(1,xsize);
00281     // *******************************************************
00282     // *******************************************************
00283     // *******************************************************
00284     // derivatives from hessian adjoint back
00285     {
00286       x_con.initialize();
00287       
00288       for (int i=1;i<=num_separable_calls;i++)
00289       {
00290         dmatrix& H=(*block_diagonal_hessian)(i);
00291         if (allocated(H))
00292         {
00293           ivector& re_list=(*block_diagonal_re_list)(i);
00294           ivector& fe_list=(*block_diagonal_fe_list)(i);
00295           dmatrix& Dux=(*block_diagonal_Dux)(i);
00296           xxx(re_list,fe_list);
00297           int mmax=re_list.indexmax();
00298           dvector tmp(1,mmax);
00299           
00300           int j;
00301           for (j=1;j<=re_list.indexmax();j++)
00302           {
00303             tmp(j)=uadjoint(re_list(j)-xmax);
00304           } 
00305   
00306           if (allocated(fe_list))
00307           {
00308             dvector tmp1=solve(H,tmp);
00309             dvector xtmp=tmp1*Dux;
00310             for (j=1;j<=fe_list.indexmax();j++)
00311             {
00312               x_con(fe_list(j))+=xtmp(j);
00313             }
00314           }
00315         }
00316 
00317       }
00318       if (initial_df1b2params::separable_flag)
00319       {
00320         //for (i=1;i<=usize;i++)
00321         //{
00322         //  Dux(i)=elem_prod(Dux(i),sscale);
00323         //}
00324         x_con=elem_prod(x_con,sscale);
00325       }
00326     }
00327     xadjoint-=x_con;
00328     // *******************************************************
00329     // *******************************************************
00330     // *******************************************************
00331     
00332     block_diagonal_flag=3;
00333     //pfmin->lapprox->xadjoint.initialize();
00334     //pfmin->lapprox->uadjoint.initialize();
00335     pfmin->lapprox->num_separable_calls=0;
00336     pfmin->lapprox->check_local_xadjoint.initialize();
00337     pfmin->lapprox->check_local_xadjoint2.initialize();
00338     pfmin->lapprox->check_local_uadjoint.initialize();
00339     pfmin->lapprox->check_local_uadjoint2.initialize();
00340     //df1b2_gradlist::set_yes_derivatives();
00341     //initial_df1b2params::reset(y,pen);
00342     pfmin->user_function();
00343     dvector lcx=elem_prod(check_local_xadjoint2,sscale);
00344     xadjoint+=lcx;
00345     //df1b2_gradlist::set_no_derivatives();
00346     funnel_init_var::lapprox=0;
00347     block_diagonal_flag=0;
00348     initial_params::set_inactive_only_random_effects(); 
00349   }
00350   else
00351   {
00352     if (function_minimizer::first_hessian_flag)
00353     {
00354       // need hessin of random effects for stddeV report
00355       df1b2variable pen=0.0;
00356       initial_df1b2params::reset(y,pen);
00357       initial_params::straight_through_flag=0;
00358       block_diagonal_flag=6;
00359       allocate_block_diagonal_stuff();
00360       num_separable_calls=0;
00361       // get the block diagonal hessians to use in importance sampling
00362       pfmin->user_function();
00363       //cout << (*pfmin->lapprox->block_diagonal_hessian) << endl;
00364       block_diagonal_flag=2;
00365     }
00366     initial_params::straight_through_flag=1;
00367     df1b2variable pen=0.0;
00368     if (saddlepointflag==2)
00369     {
00370       pmin->inner_opt_flag=0;
00371       f=get_fx_fu(pfmin);
00372     }
00373     initial_df1b2params::reset(y,pen);
00374     pmin->inner_opt_flag=1;
00375     pfmin->pre_user_function();
00376     pmin->inner_opt_flag=0;
00377     initial_params::straight_through_flag=0;
00378     if (saddlepointflag!=2)
00379     {
00380       f=initial_df1b2params::cobjfun;
00381     }
00382     else
00383     {
00384       xadjoint=(*grad_x_u)(1,xsize)-(*grad_x);
00385     }
00386 
00387     if (saddlepointflag!=2 && pmin->multinomial_weights==0)
00388     {
00389       f-=usize*.91893853320467241; 
00390     }
00391 
00392     funnel_init_var::lapprox=0;
00393     block_diagonal_flag=0;
00394     dvector scale1(1,xsize);   // need to get scale from somewhere
00395     initial_params::set_inactive_only_random_effects(); 
00396     /*int check=*/initial_params::stddev_scale(scale1,x);
00397     for (i=1;i<=xadjoint.indexmax();i++)
00398       xadjoint(i)*=scale1(i);
00399   }
00400   //cout << initial_df1b2params::cobjfun << endl;
00401   //f=initial_df1b2params::cobjfun;
00402   return xadjoint;
00403 }
00404 
00409 dvector laplace_approximation_calculator::get_newton_raphson_info_block_diagonal
00410   (function_minimizer * pfmin)
00411 {
00412   //int i,j,ip; 
00413   int ip; 
00414   
00415   int nv=initial_df1b2params::set_index();
00416   if (allocated(used_flags))
00417   {
00418     if (used_flags.indexmax() != nv)
00419     {
00420       used_flags.safe_deallocate();
00421     }
00422   }
00423   if (!allocated(used_flags))
00424   {
00425     used_flags.safe_allocate(1,nv);
00426   }
00427   
00428   //for (ip=1;ip<=num_der_blocks;ip++)
00429   {
00430     used_flags.initialize();
00431     // do we need to reallocate memory for df1b2variables?
00432     check_for_need_to_reallocate(ip);
00433     df1b2_gradlist::set_no_derivatives();
00434     //cout << re_objective_function_value::pobjfun << endl;
00435     //cout << re_objective_function_value::pobjfun->ptr << endl;
00436     (*re_objective_function_value::pobjfun)=0;
00437     df1b2variable pen=0.0;
00438     df1b2variable zz=0.0;
00439 
00440     initial_df1b2params::reset(y,pen);
00441 
00442     // call function to do block diagonal newton-raphson
00443     // the step vector from the newton-raphson is in the vector step
00444     df1b2_gradlist::set_yes_derivatives();
00445     
00446     funnel_init_var::lapprox=this;
00447     //cout << funnel_init_var::lapprox << endl;
00448     block_diagonal_flag=1;
00449     pfmin->pre_user_function();
00450     //pfmin->user_function();
00451     funnel_init_var::lapprox=0;
00452     block_diagonal_flag=0;
00453   }
00454   return step;
00455 }
00456 
00457 #endif //#if defined(USE_LAPLACE)