ADMB Documentation  11.5.3197
 All Classes Files Functions Variables Typedefs Friends Defines
df1b2gh.cpp
Go to the documentation of this file.
00001 /*
00002  * $Id$
00003  *
00004  * Author: David Fournier
00005  * Copyright (c) 2008-2012 Regents of the University of California
00006  */
00011 #  include <admodel.h>
00012 #  include <df1b2fun.h>
00013 #  include <adrndeff.h>
00014 
00019 double do_gauss_hermite_block_diagonal(const dvector& x,
00020   const dvector& u0,const dmatrix& Hess,const dvector& _xadjoint,
00021   const dvector& _uadjoint,const dmatrix& _Hessadjoint,
00022   function_minimizer * pmin)
00023 {
00024   ADUNCONST(dvector,xadjoint)
00025   ADUNCONST(dvector,uadjoint)
00026   //ADUNCONST(dmatrix,Hessadjoint)
00027   const int xs=x.size();
00028   const int us=u0.size();
00029   gradient_structure::set_NO_DERIVATIVES();
00030   int nsc=pmin->lapprox->num_separable_calls;
00031   const ivector lrea = (*pmin->lapprox->num_local_re_array)(1,nsc);
00032   int hroom =  sum(square(lrea));
00033   int nvar=x.size()+u0.size()+hroom;
00034   independent_variables y(1,nvar);
00035 
00036   // need to set random effects active together with whatever
00037   // init parameters should be active in this phase
00038   initial_params::set_inactive_only_random_effects();
00039   initial_params::set_active_random_effects();
00040   /*int onvar=*/initial_params::nvarcalc();
00041   initial_params::xinit(y);    // get the initial values into the
00042   // do we need this next line?
00043   y(1,xs)=x;
00044 
00045   // contribution for quadratic prior
00046   if (quadratic_prior::get_num_quadratic_prior()>0)
00047   {
00048     //Hess+=quadratic_prior::get_cHessian_contribution();
00049     int & vxs = (int&)(xs);
00050     quadratic_prior::get_cHessian_contribution(Hess,vxs);
00051   }
00052  // Here need hooks for sparse matrix structures
00053 
00054   dvar3_array & block_diagonal_vhessian=
00055     *pmin->lapprox->block_diagonal_vhessian;
00056   block_diagonal_vhessian.initialize();
00057   dvar3_array& block_diagonal_ch=
00058     *pmin->lapprox->block_diagonal_vch;
00059     //dvar3_array(*pmin->lapprox->block_diagonal_ch);
00060   int ii=xs+us+1;
00061   d3_array& bdH=(*pmin->lapprox->block_diagonal_hessian);
00062   for (int ic=1;ic<=nsc;ic++)
00063   {
00064     int lus=lrea(ic);
00065     for (int i=1;i<=lus;i++)
00066       for (int j=1;j<=lus;j++)
00067         y(ii++)=bdH(ic)(i,j);
00068   }
00069 
00070   dvector g(1,nvar);
00071   gradcalc(0,g);
00072   gradient_structure::set_YES_DERIVATIVES();
00073   dvar_vector vy=dvar_vector(y);
00074   //initial_params::stddev_vscale(d,vy);
00075   ii=xs+us+1;
00076   if (initial_df1b2params::have_bounded_random_effects)
00077   {
00078     cerr << "can't do importance sampling with bounded random effects"
00079      " at present" << endl;
00080     ad_exit(1);
00081   }
00082   else
00083   {
00084     for (int ic=1;ic<=nsc;ic++)
00085     {
00086       int lus=lrea(ic);
00087       if (lus>0)
00088       {
00089         for (int i=1;i<=lus;i++)
00090         {
00091           for (int j=1;j<=lus;j++)
00092           {
00093             block_diagonal_vhessian(ic,i,j)=vy(ii++);
00094           }
00095         }
00096         block_diagonal_ch(ic)=
00097           choleski_decomp(inv(block_diagonal_vhessian(ic)));
00098       }
00099     }
00100   }
00101 
00102    int nsamp=pmin->lapprox->use_gauss_hermite;
00103    pmin->lapprox->in_gauss_hermite_phase=1;
00104    dvar_vector sample_value(1,nsamp);
00105    sample_value.initialize();
00106 
00107    dvar_vector tau(1,us);;
00108    // !!! This only works for one random efect in each separable call
00109    // at present.
00110    for (int is=1;is<=nsamp;is++)
00111    {
00112      int offset=0;
00113      pmin->lapprox->num_separable_calls=0;
00114      pmin->lapprox->gh->is=is;
00115      for (int ic=1;ic<=nsc;ic++)
00116      {
00117        int lus=lrea(ic);
00118        // will need vector stuff here when more than one random effect
00119        if (lus>1)
00120        {
00121          cerr << "error not implemented" << endl;
00122          ad_exit(1);
00123        }
00124        if (lus>0)
00125        {
00126          tau(offset+1,offset+lus).shift(1)=block_diagonal_ch(ic)(1,1)*
00127            pmin->lapprox->gh->x(is);
00128          offset+=lus;
00129        }
00130      }
00131 
00132      // have to reorder the terms to match the block diagonal hessian
00133      imatrix & ls=*(pmin->lapprox->block_diagonal_re_list);
00134      int mmin=ls.indexmin();
00135      int mmax=ls.indexmax();
00136 
00137      ii=1;
00138      for (int i=mmin;i<=mmax;i++)
00139      {
00140        int cmin=ls(i).indexmin();
00141        int cmax=ls(i).indexmax();
00142        for (int j=cmin;j<=cmax;j++)
00143        {
00144          vy(ls(i,j))+=tau(ii++);
00145        }
00146      }
00147      if (ii-1 != us)
00148      {
00149        cerr << "error in interface" << endl;
00150        ad_exit(1);
00151      }
00152      initial_params::reset(vy);    // get the values into the model
00153      ii=1;
00154      for (int i=mmin;i<=mmax;i++)
00155      {
00156        int cmin=ls(i).indexmin();
00157        int cmax=ls(i).indexmax();
00158        for (int j=cmin;j<=cmax;j++)
00159        {
00160          vy(ls(i,j))-=tau(ii++);
00161        }
00162      }
00163 
00164      *objective_function_value::pobjfun=0.0;
00165      pmin->AD_uf_outer();
00166    }
00167 
00168    nsc=pmin->lapprox->num_separable_calls;
00169 
00170    dvariable vf=pmin->do_gauss_hermite_integration();
00171 
00172    int sgn=0;
00173    dvariable ld=0.0;
00174    if (ad_comm::no_ln_det_choleski_flag)
00175    {
00176      for (int ic=1;ic<=nsc;ic++)
00177      {
00178        if (allocated(block_diagonal_vhessian(ic)))
00179        {
00180          ld+=ln_det(block_diagonal_vhessian(ic),sgn);
00181        }
00182      }
00183      ld*=0.5;
00184    }
00185    else
00186    {
00187      for (int ic=1;ic<=nsc;ic++)
00188      {
00189        if (allocated(block_diagonal_vhessian(ic)))
00190        {
00191          ld+=ln_det_choleski(block_diagonal_vhessian(ic));
00192        }
00193      }
00194      ld*=0.5;
00195    }
00196 
00197    vf+=ld;
00198    //vf+=us*0.91893853320467241;
00199 
00200    double f=value(vf);
00201    gradcalc(nvar,g);
00202 
00203    // put uhat back into the model
00204    gradient_structure::set_NO_DERIVATIVES();
00205    vy(xs+1,xs+us).shift(1)=u0;
00206    initial_params::reset(vy);    // get the values into the model
00207    gradient_structure::set_YES_DERIVATIVES();
00208 
00209    pmin->lapprox->in_gauss_hermite_phase=0;
00210 
00211   ii=1;
00212   for (int i=1;i<=xs;i++)
00213     xadjoint(i)=g(ii++);
00214   for (int i=1;i<=us;i++)
00215     uadjoint(i)=g(ii++);
00216   for (int ic=1;ic<=nsc;ic++)
00217   {
00218     int lus=lrea(ic);
00219     for (int i=1;i<=lus;i++)
00220     {
00221       for (int j=1;j<=lus;j++)
00222       {
00223         (*pmin->lapprox->block_diagonal_vhessianadjoint)(ic)(i,j)=g(ii++);
00224       }
00225     }
00226   }
00227   return f;
00228 }