ADMB Documentation  11.1.2428
 All Classes Files Functions Variables Typedefs Friends Defines
df1b2gh.cpp
Go to the documentation of this file.
00001 /*
00002  * $Id: df1b2gh.cpp 1935 2014-04-26 02:02:58Z johnoel $
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   int i,j;
00046 
00047   // contribution for quadratic prior
00048   if (quadratic_prior::get_num_quadratic_prior()>0)
00049   {
00050     //Hess+=quadratic_prior::get_cHessian_contribution();
00051     int & vxs = (int&)(xs);
00052     quadratic_prior::get_cHessian_contribution(Hess,vxs);
00053   }
00054  // Here need hooks for sparse matrix structures
00055 
00056   dvar3_array & block_diagonal_vhessian=
00057     *pmin->lapprox->block_diagonal_vhessian;
00058   block_diagonal_vhessian.initialize();
00059   dvar3_array& block_diagonal_ch=
00060     *pmin->lapprox->block_diagonal_vch;
00061     //dvar3_array(*pmin->lapprox->block_diagonal_ch);
00062   int ii=xs+us+1;
00063   d3_array& bdH=(*pmin->lapprox->block_diagonal_hessian);
00064   int ic;
00065   for (ic=1;ic<=nsc;ic++)
00066   {
00067     int lus=lrea(ic);
00068     for (i=1;i<=lus;i++)
00069       for (j=1;j<=lus;j++)
00070         y(ii++)=bdH(ic)(i,j);
00071   }
00072 
00073   dvector g(1,nvar);
00074   gradcalc(0,g);
00075   gradient_structure::set_YES_DERIVATIVES();
00076   dvar_vector vy=dvar_vector(y);
00077   //initial_params::stddev_vscale(d,vy);
00078   ii=xs+us+1;
00079   if (initial_df1b2params::have_bounded_random_effects)
00080   {
00081     cerr << "can't do importance sampling with bounded random effects"
00082      " at present" << endl;
00083     ad_exit(1);
00084   }
00085   else
00086   {
00087     for (int ic=1;ic<=nsc;ic++)
00088     {
00089       int lus=lrea(ic);
00090       if (lus>0)
00091       {
00092         for (i=1;i<=lus;i++)
00093         {
00094           for (j=1;j<=lus;j++)
00095           {
00096             block_diagonal_vhessian(ic,i,j)=vy(ii++);
00097           }
00098         }
00099         block_diagonal_ch(ic)=
00100           choleski_decomp(inv(block_diagonal_vhessian(ic)));
00101       }
00102     }
00103   }
00104 
00105    int nsamp=pmin->lapprox->use_gauss_hermite;
00106    pmin->lapprox->in_gauss_hermite_phase=1;
00107    dvar_vector sample_value(1,nsamp);
00108    sample_value.initialize();
00109 
00110    dvar_vector tau(1,us);;
00111    // !!! This only works for one random efect in each separable call
00112    // at present.
00113    for (int is=1;is<=nsamp;is++)
00114    {
00115      int offset=0;
00116      pmin->lapprox->num_separable_calls=0;
00117      pmin->lapprox->gh->is=is;
00118      for (ic=1;ic<=nsc;ic++)
00119      {
00120        int lus=lrea(ic);
00121        // will need vector stuff here when more than one random effect
00122        if (lus>1)
00123        {
00124          cerr << "error not implemented" << endl;
00125          ad_exit(1);
00126        }
00127        if (lus>0)
00128        {
00129          tau(offset+1,offset+lus).shift(1)=block_diagonal_ch(ic)(1,1)*
00130            pmin->lapprox->gh->x(is);
00131          offset+=lus;
00132        }
00133      }
00134 
00135      // have to reorder the terms to match the block diagonal hessian
00136      imatrix & ls=*(pmin->lapprox->block_diagonal_re_list);
00137      int mmin=ls.indexmin();
00138      int mmax=ls.indexmax();
00139 
00140      int ii=1;
00141      int i;
00142      for (i=mmin;i<=mmax;i++)
00143      {
00144        int cmin=ls(i).indexmin();
00145        int cmax=ls(i).indexmax();
00146        for (int j=cmin;j<=cmax;j++)
00147        {
00148          vy(ls(i,j))+=tau(ii++);
00149        }
00150      }
00151      if (ii-1 != us)
00152      {
00153        cerr << "error in interface" << endl;
00154        ad_exit(1);
00155      }
00156      initial_params::reset(vy);    // get the values into the model
00157      ii=1;
00158      for (i=mmin;i<=mmax;i++)
00159      {
00160        int cmin=ls(i).indexmin();
00161        int cmax=ls(i).indexmax();
00162        for (int j=cmin;j<=cmax;j++)
00163        {
00164          vy(ls(i,j))-=tau(ii++);
00165        }
00166      }
00167 
00168      *objective_function_value::pobjfun=0.0;
00169      pmin->AD_uf_outer();
00170    }
00171 
00172    nsc=pmin->lapprox->num_separable_calls;
00173 
00174    dvariable vf=pmin->do_gauss_hermite_integration();
00175 
00176    int sgn=0;
00177    dvariable ld=0.0;
00178    if (ad_comm::no_ln_det_choleski_flag)
00179    {
00180      for (int ic=1;ic<=nsc;ic++)
00181      {
00182        if (allocated(block_diagonal_vhessian(ic)))
00183        {
00184          ld+=ln_det(block_diagonal_vhessian(ic),sgn);
00185        }
00186      }
00187      ld*=0.5;
00188    }
00189    else
00190    {
00191      for (int ic=1;ic<=nsc;ic++)
00192      {
00193        if (allocated(block_diagonal_vhessian(ic)))
00194        {
00195          ld+=ln_det_choleski(block_diagonal_vhessian(ic));
00196        }
00197      }
00198      ld*=0.5;
00199    }
00200 
00201    vf+=ld;
00202    //vf+=us*0.91893853320467241;
00203 
00204    double f=value(vf);
00205    gradcalc(nvar,g);
00206 
00207    // put uhat back into the model
00208    gradient_structure::set_NO_DERIVATIVES();
00209    vy(xs+1,xs+us).shift(1)=u0;
00210    initial_params::reset(vy);    // get the values into the model
00211    gradient_structure::set_YES_DERIVATIVES();
00212 
00213    pmin->lapprox->in_gauss_hermite_phase=0;
00214 
00215   ii=1;
00216   for (i=1;i<=xs;i++)
00217     xadjoint(i)=g(ii++);
00218   for (i=1;i<=us;i++)
00219     uadjoint(i)=g(ii++);
00220   for (ic=1;ic<=nsc;ic++)
00221   {
00222     int lus=lrea(ic);
00223     for (i=1;i<=lus;i++)
00224     {
00225       for (j=1;j<=lus;j++)
00226       {
00227         (*pmin->lapprox->block_diagonal_vhessianadjoint)(ic)(i,j)=g(ii++);
00228       }
00229     }
00230   }
00231   return f;
00232 }