ADMB Documentation  11.1.1017
 All Classes Files Functions Variables Typedefs Friends Defines
mod_rhes.cpp
Go to the documentation of this file.
00001 /*
00002  * $Id: mod_rhes.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 #include <sstream>
00012 using std::istringstream;
00013 
00014 #if defined(USE_LAPLACE)
00015 #include <admodel.h>
00016 #include <df1b2fun.h>
00017 #include <adrndeff.h>
00018 #ifdef __GNUDOS__
00019   #include <gccmanip.h>
00020 #endif
00021 
00022 void get_inverse_sparse_hessian(dcompressed_triplet & st, hs_symbolic& S,
00023   uostream& ofs1,ofstream& ofs,int usize,int xsize,dvector& u);
00024 
00029 banded_lower_triangular_dmatrix quiet_choleski_decomp(
00030   const banded_symmetric_dmatrix& _M,const int& _ierr)
00031 {
00032   int & ierr = (int &) _ierr;
00033   ADUNCONST(banded_symmetric_dmatrix,M)
00034   int minsave=M.indexmin();
00035   M.shift(1);
00036   int n=M.indexmax();
00037   
00038   int bw=M.bandwidth();
00039   banded_lower_triangular_dmatrix L(1,n,bw);
00040 #ifndef SAFE_INITIALIZE
00041     L.initialize();
00042 #endif
00043 
00044   int i,j,k;
00045   double tmp;
00046     if (M(1,1)<=0)
00047     {
00048       ierr=1;
00049       return L;
00050     }
00051   L(1,1)=sqrt(M(1,1));
00052   for (i=2;i<=bw;i++)
00053   {
00054     L(i,1)=M(i,1)/L(1,1);
00055   }
00056 
00057   for (i=2;i<=n;i++)
00058   {
00059     for (j=i-bw+1;j<=i-1;j++)
00060     {
00061       if (j>1)
00062       { 
00063         tmp=M(i,j);
00064         for (k=i-bw+1;k<=j-1;k++)
00065         {
00066     if (k>0 && k>j-bw)
00067             tmp-=L(i,k)*L(j,k);
00068         }
00069         L(i,j)=tmp/L(j,j);
00070       }
00071     }
00072     tmp=M(i,i);
00073     for (k=i-bw+1;k<=i-1;k++)
00074     {
00075       if (k>0)  
00076         tmp-=L(i,k)*L(i,k);
00077     }
00078     if (tmp<=0)
00079     {
00080       ierr=1;
00081       return L;
00082     }
00083     L(i,i)=sqrt(tmp);
00084   }
00085   M.shift(minsave);
00086   L.shift(minsave);
00087 
00088   return L;
00089 }
00090 
00091 void useless(const double& sdelta2);
00092 
00097 void function_minimizer::hess_routine_random_effects(void)
00098 {
00099 #if defined(USE_ADPVM)
00100   if (ad_comm::pvm_manager)
00101   {
00102     switch (ad_comm::pvm_manager->mode)
00103     {
00104     case 1: //master
00105       hess_routine_noparallel_random_effects();
00106       break;
00107     case 2: //slave
00108       hess_routine_slave_random_effects();
00109       break;
00110     default:
00111       cerr << "Illegal value for mode" << endl;
00112       ad_exit(1);
00113     }
00114   }
00115   else
00116 #endif
00117   {
00118       hess_routine_noparallel_random_effects();
00119   }
00120 }
00121 dvector get_solution_vector(int npts);
00122 
00127 void function_minimizer::hess_routine_noparallel_random_effects(void)
00128 {
00129 
00130   int nvar=initial_params::nvarcalc(); // get the number of active parameters
00131   //if (adjm_ptr) set_labels_for_hess(nvar);
00132   independent_variables x(1,nvar);
00133   initial_params::xinit(x);        // get the initial values into the x vector
00134   double f;
00135   double delta=1.e-4;
00136   dvector g1(1,nvar);
00137   dvector g0(1,nvar);
00138   dvector g2(1,nvar);
00139   dvector gbest(1,nvar);
00140   dvector hess(1,nvar);
00141   dvector hess1(1,nvar);
00142   dvector hess2(1,nvar);
00143   //double eps=.1;
00144   gradient_structure::set_YES_DERIVATIVES();
00145   gbest.fill_seqadd(1.e+50,0.);
00146 
00147     
00148   dvector ddd(1,nvar);
00149   gradcalc(0,ddd);
00150   adstring tmpstring;
00151 
00152   {
00153     first_hessian_flag=1;
00154     {
00155       g1=(*lapprox)(x,f,this);
00156       g0=g1;
00157     }
00158     // modify so thaqt we have l_uu and dux for delta method
00159     // DF feb 15 05
00160     //if (lapprox->hesstype==2 || lapprox->hesstype==3) 
00161     if (lapprox->hesstype==2 ) 
00162     {
00163       if (lapprox->block_diagonal_hessian)
00164       {
00165         //if (ad_comm::wd_flag)
00166         tmpstring = ad_comm::adprogram_name + ".rhes";
00167         ofstream ofs((char*)(tmpstring));
00168             ofs << "   value      std.dev" << endl;
00169         int mmin=lapprox->block_diagonal_hessian->indexmin();
00170         int mmax=lapprox->block_diagonal_hessian->indexmax();
00171         int i,j;
00172         int ii=1;
00173         dvector & u= lapprox->uhat;
00174         for (i=mmin;i<=mmax;i++)
00175         {
00176           if (allocated((*(lapprox->block_diagonal_hessian))(i)))
00177           {
00178             dmatrix m= inv((*(lapprox->block_diagonal_hessian))(i));
00179             dvector d=sqrt(diagonal(m));
00180             int jmin=d.indexmin();
00181             int jmax=d.indexmax();
00182             for (j=jmin;j<=jmax;j++)
00183             {
00184               //if (ii<=u.indexmax())
00185               {
00186                 ofs << setprecision(5) << setscientific() 
00187                     << setw(14) << u(ii++) << " " << d(j) << endl;;
00188               }
00189             }
00190           }
00191         }
00192       }
00193       else if (lapprox->bHess)
00194       {
00195         //if (ad_comm::wd_flag)
00196         tmpstring = ad_comm::adprogram_name + ".rhes";
00197         ofstream ofs((char*)(tmpstring));
00198             ofs << "   value      std.dev" << endl;
00199         int mmin=lapprox->bHess->indexmin();
00200         int mmax=lapprox->bHess->indexmax();
00201         //int i,j;
00202         int i;
00203         //int ii=1;
00204         dvector & u= lapprox->uhat;
00205         dvector e(mmin,mmax);
00206         //choleski_decomp(*lapprox->bHess);
00207         int ierr;
00208         
00209         banded_lower_triangular_dmatrix tmp=choleski_decomp(*lapprox->bHess,
00210           ierr);
00211         e.initialize();
00212         for (i=mmin;i<=mmax;i++)
00213         {
00214           e(i)=1.0;
00215           dvector v=solve(tmp,e);
00216           e(i)=0;
00217           
00218           double d=sqrt(v*v);
00219             ofs << setprecision(5) << setscientific() 
00220                 << setw(14) << u(i) << " " << d << endl;;
00221         }
00222       }
00223     }
00224     else
00225     {
00226       int i;
00227       //if (ad_comm::wd_flag)
00228       dmatrix m;
00229       tmpstring = ad_comm::adprogram_name + ".rhes";
00230       ofstream ofs((char*)(tmpstring));
00231           ofs << "   value      std.dev" << endl;
00232       //int ii=1;
00233       tmpstring = ad_comm::adprogram_name + ".luu";
00234       uostream ofs1((char*)(tmpstring));
00235       dvector & u= lapprox->uhat;
00236       if (lapprox->hesstype !=3)
00237       {
00238         if (allocated(lapprox->Hess))
00239         {
00240           m= inv(lapprox->Hess);
00241           int mmin=m.indexmin();
00242           int mmax=m.indexmax();
00243           for (i=mmin;i<=mmax;i++)
00244           {
00245             ofs << setprecision(5) << setscientific() 
00246                 << setw(14) << u(i) << " " << sqrt(m(i,i)) << endl;;
00247           }
00248           // save l_uu and l_xu for covariance calculations
00249           ofs1 << lapprox->usize << lapprox->xsize;
00250           ofs1 << m;
00251         }
00252         else if (lapprox->sparse_triplet2)
00253         {
00254           dcompressed_triplet & st= *(lapprox->sparse_triplet2);
00255           hs_symbolic& S= *(lapprox->sparse_symbolic2);
00256           get_inverse_sparse_hessian(st,S,ofs1,ofs,lapprox->usize,
00257             lapprox->xsize,u);
00258           // save l_uu and l_xu for covariance calculations
00259         }
00260       }
00261       else
00262       {
00263         if (lapprox->bHess)
00264         {
00265           int ierr=0;
00266           int mmin=lapprox->bHess->indexmin();
00267           int mmax=lapprox->bHess->indexmax();
00268           const banded_lower_triangular_dmatrix& C=
00269             quiet_choleski_decomp(*lapprox->bHess,ierr);
00270           ivector e(mmin,mmax);
00271           e.initialize();
00272           if (ierr==0)
00273           {
00274             ofs1 << lapprox->usize << lapprox->xsize;
00275             for (int i=mmin;i<=mmax;i++)
00276             {
00277               if (i>1) e(i-1)=0;
00278               e(i)=1;
00279               dvector w=solve_trans(C,solve(C,e));
00280               ofs << setprecision(5) << setscientific() 
00281                   << setw(14) << u(i) << " " << sqrt(w(i)) << endl;;
00282               ofs1 << w;
00283             }
00284           }
00285           else
00286           {
00287           }
00288         }
00289       }
00290       if (!ofs)
00291       {
00292         cerr << "Error writing to file " << tmpstring << endl;
00293         ad_exit(1);
00294       }
00295       // save l_uu and l_xu for covariance calculations
00296       ofs1 << lapprox->Dux;
00297       if (!ofs1)
00298       {
00299         cerr << "Error writing to file " << tmpstring << endl;
00300         ad_exit(1);
00301       }
00302       ofs1.close();
00303     }
00304 
00305     {
00306       int i,j;
00307       tmpstring = ad_comm::adprogram_name + ".luu";
00308       uistream uis1((char*)(tmpstring));
00309       uis1 >> i >> j;
00310       cout << i << " " << j << endl;
00311     }
00312 
00313     int npts=2;
00314     int on,nopt;
00315     if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-hpts",nopt))>-1)
00316     {
00317       if (nopt !=1)
00318       {
00319         cerr << "Usage -hpts option needs non-negative integer  -- ignored" << endl;
00320       }
00321       else
00322       {   
00323         npts=atoi(ad_comm::argv[on+1]);
00324       }
00325     }
00326 
00327     double _delta=0.0;
00328 
00329     if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-hsize",nopt))>-1)
00330     {
00331       if (!nopt)
00332       {
00333         cerr << "Usage -hsize option needs number  -- ignored" << endl;
00334       }
00335       else
00336       {   
00337     
00338         istringstream ist(ad_comm::argv[on+1]);
00339         ist >> _delta;
00340     
00341         if (_delta<=0)
00342         {
00343           cerr << "Usage -hsize option needs positive number  -- ignored" << endl;
00344           _delta=0.0;
00345         } 
00346       }
00347       if (_delta>0) 
00348       {
00349         delta=_delta;
00350       }
00351     }
00352 
00353     // get a number which is exactly representable
00354     double sdelta=1.0+delta;
00355     useless(sdelta);
00356     sdelta-=1.0;
00357 
00358     
00359     {
00360       //
00361       uostream uos("hessian.bin");
00362       uos << npts;
00363       for (int i=1;i<=nvar;i++)
00364       {
00365         cout << "Estimating row " << i << " out of " << nvar
00366        << " for hessian" << endl;
00367   
00368         for (int j=-npts;j<=npts;j++)
00369         {
00370           if (j !=0)
00371           {
00372             double f=0.0;
00373             double xsave=x(i);
00374             x(i)=xsave+j*sdelta;
00375             g1=(*lapprox)(x,f,this);
00376             x(i)=xsave;
00377             uos << i << j << sdelta << g1;
00378           }
00379           else
00380           {
00381             uos << i << j << sdelta << g0;
00382           }
00383         }
00384       }
00385     }
00386     // check for accuracy
00387     {
00388       double sd;
00389       uistream uis("hessian.bin");
00390       uis >> npts;
00391       dvector v=get_solution_vector(npts);
00392       v.shift(-npts);
00393       dmatrix tmp(-npts,npts,1,nvar);
00394       dmatrix hess(1,nvar,1,nvar);
00395       ivector iind(-npts,npts);
00396       ivector jind(-npts,npts);
00397       int i;
00398       for (i=1;i<=nvar;i++)
00399       {
00400         for (int j=-npts;j<=npts;j++)
00401         {
00402           uis >> iind(j) >> jind(j) >> sd >> tmp(j);
00403         }
00404         hess(i)=(v*tmp).shift(1);
00405         hess(i)/=sd;
00406       }
00407       {
00408         adstring tmpstring="admodel.hes";
00409         if (ad_comm::wd_flag)
00410         {
00411           tmpstring = ad_comm::adprogram_name + ".hes";
00412         }
00413         uostream ofs((char*)tmpstring);
00414         ofs << nvar;
00415         dmatrix shess(1,nvar,1,nvar);
00416         double maxerr=0.0;
00417         for (i=1;i<=nvar;i++)
00418         {
00419           for (int j=1;j<=nvar;j++)
00420           {
00421             shess(i,j)=(hess(i,j)-hess(j,i))/
00422              (1.e-3+sfabs(hess(i,j))+fabs(hess(j,i)));
00423             if (shess(i,j)>maxerr) maxerr=shess(i,j);
00424           }
00425         }
00426         ofstream ofs1("hesscheck");
00427         ofs1 << "maxerr = " << maxerr << endl << endl;
00428         ofs1 << setw(10) << hess << endl << endl;
00429         ofs1 << setw(10) << shess << endl;
00430         ofs << hess;
00431         ofs << gradient_structure::Hybrid_bounded_flag;
00432         initial_params::set_inactive_only_random_effects(); 
00433         dvector tscale(1,nvar);   // need to get scale from somewhere
00434         /*int check=*/initial_params::stddev_scale(tscale,x);
00435         ofs << tscale;
00436       }
00437     }
00438    /*
00439     first_hessian_flag=0;
00440     double sdelta1;
00441     double sdelta2;
00442     lapprox->fmc1.fringe=1.e-9;
00443     for (int i=1;i<=nvar;i++)
00444     {
00445 #if defined (__SPDLL__)
00446       hess_calcreport(i,nvar);
00447 #else
00448       cout << "Estimating row " << i << " out of " << nvar
00449      << " for hessian" << endl;
00450 #endif
00451 
00452       double f=0.0;
00453       double xsave=x(i);
00454       sdelta1=x(i)+delta;
00455       useless(sdelta1);
00456       sdelta1-=x(i);
00457       x(i)=xsave+sdelta1;
00458       g1=(*lapprox)(x,f,this);
00459 
00460       sdelta2=x(i)-delta;
00461       useless(sdelta2);
00462       sdelta2-=x(i);
00463       x(i)=xsave+sdelta2;
00464       g2=(*lapprox)(x,f,this);
00465       x(i)=xsave;
00466       hess1=(g1-g2)/(sdelta1-sdelta2);
00467  
00468       sdelta1=x(i)+eps*delta;
00469       useless(sdelta1);
00470       sdelta1-=x(i);
00471       x(i)=xsave+sdelta1;
00472       g1=(*lapprox)(x,f,this);
00473 
00474       x(i)=xsave-eps*delta;
00475       sdelta2=x(i)-eps*delta;
00476       useless(sdelta2);
00477       sdelta2-=x(i);
00478       x(i)=xsave+sdelta2;
00479       g2=(*lapprox)(x,f,this);
00480       x(i)=xsave;
00481 
00482       initial_params::set_inactive_only_random_effects(); 
00483       initial_params::reset(dvar_vector(x));
00484       double eps2=eps*eps;
00485       hess2=(g1-g2)/(sdelta1-sdelta2);
00486       hess=(eps2*hess1-hess2) /(eps2-1.);
00487    
00488       ofs << hess;
00489     }
00490    */
00491   }
00492 }
00493 
00494 #if defined(USE_ADPVM)
00495 
00500 void function_minimizer::hess_routine_slave_random_effects(void)
00501 {
00502 
00503   int nvar=initial_params::nvarcalc(); // get the number of active parameters
00504   //if (adjm_ptr) set_labels_for_hess(nvar);
00505   independent_variables x(1,nvar);
00506   initial_params::xinit(x);        // get the initial values into the x vector
00507   double f;
00508   double delta=1.e-6;
00509   double eps=.1;
00510   gradient_structure::set_YES_DERIVATIVES();
00511 
00512     
00513   dvector ddd(1,nvar);
00514   gradcalc(0,ddd);
00515   {
00516     {
00517       (*lapprox)(x,f,this);
00518     }
00519     double sdelta1;
00520     double sdelta2;
00521     lapprox->fmc1.fringe=1.e-9;
00522     for (int i=1;i<=nvar;i++)
00523     {
00524       double f=0.0;
00525       double xsave=x(i);
00526       sdelta1=x(i)+delta;
00527       useless(sdelta1);
00528       sdelta1-=x(i);
00529       x(i)=xsave+sdelta1;
00530       (*lapprox)(x,f,this);
00531 
00532       sdelta2=x(i)-delta;
00533       useless(sdelta2);
00534       sdelta2-=x(i);
00535       x(i)=xsave+sdelta2;
00536       (*lapprox)(x,f,this);
00537       x(i)=xsave;
00538  
00539       sdelta1=x(i)+eps*delta;
00540       useless(sdelta1);
00541       sdelta1-=x(i);
00542       x(i)=xsave+sdelta1;
00543       (*lapprox)(x,f,this);
00544 
00545       x(i)=xsave-eps*delta;
00546       sdelta2=x(i)-eps*delta;
00547       useless(sdelta2);
00548       sdelta2-=x(i);
00549       x(i)=xsave+sdelta2;
00550       (*lapprox)(x,f,this);
00551       x(i)=xsave;
00552 
00553       initial_params::set_inactive_only_random_effects(); 
00554       initial_params::reset(dvar_vector(x));
00555    
00556     }
00557   }
00558 }
00559 #endif // #if defined(USE_ADPVM)
00560 
00565 dvector get_solution_vector(int n)
00566 {
00567   int i;
00568   int n1=2*n+1;
00569   dmatrix tmp(1,n1,1,n1);
00570   dvector v(1,n1);
00571   v.initialize();
00572   v(2)=1.0; 
00573   dvector c(1,n1);
00574   c.fill_seqadd(-n,1);
00575   tmp.initialize();
00576   tmp(1)=1;
00577   tmp(2)=c;
00578   for (i=3;i<=n1;i++)
00579   {
00580     tmp(i)=elem_prod(tmp(i-1),c);
00581   }
00582   dmatrix tmp1=inv(tmp);
00583   return tmp1*v;
00584 }
00585 
00586 #endif