ADMB Documentation  11.5.3197
 All Classes Files Functions Variables Typedefs Friends Defines
df1b2im5.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 calculate_importance_sample_block_diagonal_option_antithetical
00020   (const dvector& x,const dvector& u0,const dmatrix& Hess,
00021   const dvector& _xadjoint,const dvector& _uadjoint,
00022   const dmatrix& _Hessadjoint,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       for (int i=1;i<=lus;i++)
00088       {
00089         for (int j=1;j<=lus;j++)
00090         {
00091           block_diagonal_vhessian(ic,i,j)=vy(ii++);
00092         }
00093       }
00094       block_diagonal_ch(ic)=
00095         choleski_decomp(inv(block_diagonal_vhessian(ic)));
00096     }
00097   }
00098 
00099    int nsamp=pmin->lapprox->num_importance_samples;
00100 
00101    dvariable vf=0.0;
00102 
00103    dvar_vector sample_value(1,nsamp);
00104    sample_value.initialize();
00105 
00106    dvar_vector tau(1,us);;
00107    int is;
00108    for (is=1;is<=nsamp;is++)
00109    {
00110      int offset=0;
00111      pmin->lapprox->importance_sampling_counter=is;
00112      for (int ic=1;ic<=nsc;ic++)
00113      {
00114        int lus=lrea(ic);
00115        tau(offset+1,offset+lus).shift(1)=block_diagonal_ch(ic)*
00116          (*pmin->lapprox->antiepsilon)(is);
00117        offset+=lus;
00118      }
00119 
00120      // have to reorder the terms to match the block diagonal hessian
00121      imatrix & ls=*(pmin->lapprox->block_diagonal_re_list);
00122      int mmin=ls.indexmin();
00123      int mmax=ls.indexmax();
00124 
00125      ii=1;
00126      for (int i=mmin;i<=mmax;i++)
00127      {
00128        int cmin=ls(i).indexmin();
00129        int cmax=ls(i).indexmax();
00130        for (int j=cmin;j<=cmax;j++)
00131        {
00132          vy(ls(i,j))+=tau(ii++);
00133        }
00134      }
00135      if (ii-1 != us)
00136      {
00137        cerr << "error in interface" << endl;
00138        ad_exit(1);
00139      }
00140      initial_params::reset(vy);    // get the values into the model
00141      ii=1;
00142      for (int 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 
00152      *objective_function_value::pobjfun=0.0;
00153      pmin->AD_uf_outer();
00154 
00155      if (pmin->lapprox->use_outliers==0)
00156      {
00157        // assumes that all separable calls have the same number
00158        // of random effects
00159        double neps=0.5*nsc*norm2((*pmin->lapprox->antiepsilon)(is));
00160 
00161        (*pmin->lapprox->importance_sampling_values)(is)=
00162            neps-value(*objective_function_value::pobjfun);
00163 
00164        (*pmin->lapprox->importance_sampling_weights)(is)=neps;
00165 
00166         sample_value(is)=*objective_function_value::pobjfun
00167           -neps;
00168      }
00169      else
00170      {
00171        dvector& e=pmin->lapprox->epsilon(is);
00172        double neps=-sum(log(.95*exp(-0.5*square(e))+
00173           0.0166666667*exp(-square(e)/18.0)));
00174 
00175        (*pmin->lapprox->importance_sampling_values)(is)=
00176          neps-value(*objective_function_value::pobjfun);
00177 
00178        sample_value(is)=*objective_function_value::pobjfun
00179          -neps;
00180      }
00181    }
00182 
00183    nsc=pmin->lapprox->num_separable_calls;
00184    dmatrix weights(1,nsc,1,nsamp);
00185    for (is=1;is<=nsamp;is++)
00186    {
00187      int offset=0;
00188      for (int ic=1;ic<=nsc;ic++)
00189      {
00190        int lus=lrea(ic);
00191        // assumes that all spearable calls have the same number of
00192        // random effects
00193        dvector e= (*pmin->lapprox->antiepsilon)(is);
00194        offset+=lus;
00195        if (pmin->lapprox->use_outliers==0)
00196        {
00197          weights(ic,is)=0.5*norm2(e);
00198        }
00199        else
00200        {
00201          weights(ic,is)=-sum(log(.95*exp(-0.5*square(e))+
00202           0.0166666667*exp(-square(e)/18.0)));
00203        }
00204      }
00205    }
00206    dvar_vector lcomp(1,nsc);
00207    for (int isc=1;isc<=nsc;isc++)
00208    {
00209      dvar_vector & comps=
00210        (*pmin->lapprox->importance_sampling_components)(isc);
00211 
00212      dvar_vector diff=comps-weights(isc);
00213      double dmin=min(value(diff));
00214      lcomp(isc)=dmin-log(mean(exp(dmin-diff)));
00215    }
00216 
00217    //double ns=lcomp.indexmax()-lcomp.indexmin()+1;
00218    //double min_vf=min(value(lcomp));
00219    vf= sum(lcomp);
00220    vf-=us*0.91893853320467241;
00221 
00222    int sgn=0;
00223    dvariable ld=0.0;
00224    if (ad_comm::no_ln_det_choleski_flag)
00225    {
00226      for (int ic=1;ic<=nsc;ic++)
00227      {
00228        ld+=ln_det(block_diagonal_vhessian(ic),sgn);
00229      }
00230      ld*=0.5;
00231    }
00232    else
00233    {
00234      for (int ic=1;ic<=nsc;ic++)
00235      {
00236        ld+=ln_det_choleski(block_diagonal_vhessian(ic));
00237      }
00238      ld*=0.5;
00239    }
00240 
00241    vf+=ld;
00242 
00243    double f=value(vf);
00244    gradcalc(nvar,g);
00245 
00246    // put uhat back into the model
00247    gradient_structure::set_NO_DERIVATIVES();
00248    vy(xs+1,xs+us).shift(1)=u0;
00249    initial_params::reset(vy);    // get the values into the model
00250    gradient_structure::set_YES_DERIVATIVES();
00251 
00252   ii=1;
00253   for (int i=1;i<=xs;i++)
00254     xadjoint(i)=g(ii++);
00255   for (int i=1;i<=us;i++)
00256     uadjoint(i)=g(ii++);
00257   for (int ic=1;ic<=nsc;ic++)
00258   {
00259     int lus=lrea(ic);
00260     for (int i=1;i<=lus;i++)
00261     {
00262       for (int j=1;j<=lus;j++)
00263       {
00264         (*pmin->lapprox->block_diagonal_vhessianadjoint)(ic)(i,j)=g(ii++);
00265       }
00266     }
00267   }
00268   return f;
00269 }