ADMB Documentation  11.5.3197
 All Classes Files Functions Variables Typedefs Friends Defines
mod_hess.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  */
00007 #include <df1b2fun.h>
00008 #include <admodel.h>
00009 //#include <parallel.h>
00010 
00011 void hess_calcreport(int i,int nvar);
00012 void hess_errorreport(void);
00013 void set_labels_for_hess(int);
00014 
00015 // estimate the matrix of second derivatives
00016 void ad_update_hess_stats_report(int i,int nvar);
00017 
00018 void function_minimizer::hess_routine(void)
00019 {
00020   if (random_effects_flag && lapprox != 0)
00021   {
00022     if (laplace_approximation_calculator::alternative_user_function_flag == 1)
00023     {
00024       laplace_approximation_calculator::alternative_user_function_flag = 2;
00025     }
00026     hess_routine_random_effects();
00027     if (laplace_approximation_calculator::alternative_user_function_flag == 2)
00028     {
00029       laplace_approximation_calculator::alternative_user_function_flag = 1;
00030     }
00031   }
00032 #if defined(USE_ADPVM)
00033   else if (ad_comm::pvm_manager)
00034   {
00035     switch (ad_comm::pvm_manager->mode)
00036     {
00037     case 1: // master
00038       hess_routine_master();
00039       break;
00040     case 2: // slave
00041       hess_routine_slave();
00042       break;
00043     default:
00044       cerr << "Error: Illegal value for pvm_manager->mode." << endl;
00045       ad_exit(1);
00046     }
00047     cout << "finished hess routine" << endl;
00048   }
00049 #endif
00050   else
00051   {
00052     hess_routine_noparallel();
00053   }
00054 }
00055 void function_minimizer::hess_routine_noparallel(void)
00056 {
00057   int nvar=initial_params::nvarcalc(); // get the number of active parameters
00058   //if (adjm_ptr) set_labels_for_hess(nvar);
00059   independent_variables x(1,nvar);
00060   initial_params::xinit(x);        // get the initial values into the x vector
00061   double delta=1.e-5;
00062   dvector g1(1,nvar);
00063   dvector g2(1,nvar);
00064   dvector gbest(1,nvar);
00065   dvector hess(1,nvar);
00066   dvector hess1(1,nvar);
00067   dvector hess2(1,nvar);
00068   double eps=.1;
00069   double eps2=eps*eps;
00070   gradient_structure::set_YES_DERIVATIVES();
00071   gbest.fill_seqadd(1.e+50,0.);
00072 
00073   adstring tmpstring="admodel.hes";
00074   if (ad_comm::wd_flag)
00075      tmpstring = ad_comm::adprogram_name + ".hes";
00076   uostream ofs((char*)tmpstring);
00077 
00078   ofs << nvar;
00079   {
00080     {
00081       dvariable vf=0.0;
00082       vf=initial_params::reset(dvar_vector(x));
00083       *objective_function_value::pobjfun=0.0;
00084       pre_userfunction();
00085       vf+=*objective_function_value::pobjfun;
00086       gradcalc(nvar, g1, vf);
00087     }
00088     double sdelta1;
00089     double sdelta2;
00090     for (int i=1;i<=nvar;i++)
00091     {
00092       hess_calcreport(i,nvar);
00093 
00094       double xsave=x(i);
00095       sdelta1=x(i)+delta;
00096       sdelta1-=x(i);
00097       x(i)=xsave+sdelta1;
00098       dvariable vf=0.0;
00099       vf=initial_params::reset(dvar_vector(x));
00100       *objective_function_value::pobjfun=0.0;
00101       pre_userfunction();
00102       vf+=*objective_function_value::pobjfun;
00103       gradcalc(nvar, g1, vf);
00104 
00105       sdelta2=x(i)-delta;
00106       sdelta2-=x(i);
00107       x(i)=xsave+sdelta2;
00108       vf=0.0;
00109       vf=initial_params::reset(dvar_vector(x));
00110       *objective_function_value::pobjfun=0.0;
00111       pre_userfunction();
00112       vf+=*objective_function_value::pobjfun;
00113       gradcalc(nvar, g2, vf);
00114       x(i)=xsave;
00115       hess1=(g1-g2)/(sdelta1-sdelta2);
00116 
00117       sdelta1=x(i)+eps*delta;
00118       sdelta1-=x(i);
00119       x(i)=xsave+sdelta1;
00120       vf=0.0;
00121       vf=initial_params::reset(dvar_vector(x));
00122       *objective_function_value::pobjfun=0.0;
00123       pre_userfunction();
00124       vf+=*objective_function_value::pobjfun;
00125       gradcalc(nvar, g1, vf);
00126 
00127       x(i)=xsave-eps*delta;
00128       sdelta2=x(i)-eps*delta;
00129       sdelta2-=x(i);
00130       x(i)=xsave+sdelta2;
00131       vf=0.0;
00132       vf=initial_params::reset(dvar_vector(x));
00133       *objective_function_value::pobjfun=0.0;
00134       pre_userfunction();
00135       vf+=*objective_function_value::pobjfun;
00136       gradcalc(nvar, g2, vf);
00137       x(i)=xsave;
00138 
00139       hess2=(g1-g2)/(sdelta1-sdelta2);
00140 
00141       hess=(eps2*hess1-hess2) /(eps2-1.);
00142 
00143       ofs << hess;
00144       //if (adjm_ptr) ad_update_hess_stats_report(nvar,i);
00145     }
00146   }
00147   initial_params::reset(dvar_vector(x));
00148   ofs << gradient_structure::Hybrid_bounded_flag;
00149   dvector tscale(1,nvar);   // need to get scale from somewhere
00150   /*int check=*/initial_params::stddev_scale(tscale,x);
00151   ofs << tscale;
00152 }
00153 
00154 void function_minimizer::hess_routine_and_constraint(int iprof,
00155   const dvector& g, dvector& fg)
00156 {
00157   int nvar=initial_params::nvarcalc(); // get the number of active parameters
00158   independent_variables x(1,nvar);
00159   initial_params::xinit(x);        // get the initial values into the x vector
00160   double delta=1.e-6;
00161   dvector g1(1,nvar);
00162   dvector g2(1,nvar);
00163   dvector gbest(1,nvar);
00164   dvector hess(1,nvar);
00165   dvector hess1(1,nvar);
00166   dvector hess2(1,nvar);
00167   //double eps=.1;
00168   gradient_structure::set_YES_DERIVATIVES();
00169   gbest.fill_seqadd(1.e+50,0.);
00170   uostream ofs("admodel.hes");
00171   //ofstream ofs5("tmphess");
00172   double lambda=fg*g/norm2(g);
00173   cout << fg-lambda*g << endl;
00174   cout << norm(fg-lambda*g) << " " << fg*g/(norm(g)*norm(fg)) << endl;
00175   ofs << nvar;
00176   {
00177     {
00178       dvariable vf=0.0;
00179       vf=initial_params::reset(dvar_vector(x));
00180       *objective_function_value::pobjfun=0.0;
00181       pre_userfunction();
00182       vf+=*objective_function_value::pobjfun;
00183       vf-=lambda*likeprof_params::likeprofptr[iprof]->variable();
00184       gradcalc(nvar, g1, vf);
00185     }
00186     double sdelta1;
00187     double sdelta2;
00188 
00189     for (int i=1;i<=nvar;i++)
00190     {
00191       hess_calcreport(i,nvar);
00192 
00193       double xsave=x(i);
00194       sdelta1=x(i)+delta;
00195       sdelta1-=x(i);
00196       x(i)=xsave+sdelta1;
00197       dvariable vf=0.0;
00198       vf=initial_params::reset(dvar_vector(x));
00199       *objective_function_value::pobjfun=0.0;
00200       pre_userfunction();
00201       vf+=*objective_function_value::pobjfun;
00202       vf-=lambda*likeprof_params::likeprofptr[iprof]->variable();
00203       gradcalc(nvar, g1, vf);
00204 
00205       sdelta2=x(i)-delta;
00206       sdelta2-=x(i);
00207       x(i)=xsave+sdelta2;
00208       vf=0.0;
00209       vf=initial_params::reset(dvar_vector(x));
00210       *objective_function_value::pobjfun=0.0;
00211       pre_userfunction();
00212       vf+=*objective_function_value::pobjfun;
00213       vf-=lambda*likeprof_params::likeprofptr[iprof]->variable();
00214       gradcalc(nvar, g2, vf);
00215       x(i)=xsave;
00216       hess1=(g1-g2)/(sdelta1-sdelta2);
00217   /*
00218       sdelta1=x(i)+eps*delta;
00219       sdelta1-=x(i);
00220       x(i)=xsave+sdelta1;
00221       vf=0.0;
00222       vf=initial_params::reset(dvar_vector(x));
00223       *objective_function_value::pobjfun=0.0;
00224       pre_userfunction();
00225       vf+=*objective_function_value::pobjfun;
00226       vf-=lambda*likeprof_params::likeprofptr[iprof]->variable();
00227       f=value(vf);
00228       gradcalc(nvar,g1);
00229 
00230       x(i)=xsave-eps*delta;
00231       sdelta2=x(i)-eps*delta;
00232       sdelta2-=x(i);
00233       x(i)=xsave+sdelta2;
00234       vf=0.0;
00235       vf=0.0;
00236       vf=initial_params::reset(dvar_vector(x));
00237       *objective_function_value::pobjfun=0.0;
00238       pre_userfunction();
00239       vf+=*objective_function_value::pobjfun;
00240       vf-=lambda*likeprof_params::likeprofptr[iprof]->variable();
00241       f=value(vf);
00242       gradcalc(nvar,g2);
00243       x(i)=xsave;
00244 
00245       double eps2=eps*eps;
00246       hess2=(g1-g2)/(sdelta1-sdelta2);
00247       hess=(eps2*hess1-hess2) /(eps2-1.);
00248     */
00249       hess=hess1;
00250       ofs << hess;
00251     }
00252   }
00253   gradient_structure::set_NO_DERIVATIVES();
00254 }
00255 /*
00256 void function_minimizer::hess_routine_and_constraint(int iprof)
00257 {
00258   int nvar=initial_params::nvarcalc(); // get the number of active parameters
00259   independent_variables x(1,nvar);
00260   initial_params::xinit(x);        // get the initial values into the x vector
00261   double f;
00262   double delta=1.e-6;
00263   dvector g1(1,nvar);
00264   dvector g2(1,nvar);
00265   dvector gbest(1,nvar);
00266   dvector hess(1,nvar);
00267   dvector hess1(1,nvar);
00268   dvector hess2(1,nvar);
00269   double eps=.1;
00270   gradient_structure::set_YES_DERIVATIVES();
00271   gbest.fill_seqadd(1.e+50,0.);
00272   uostream ofs("admodel.hes");
00273   //ofstream ofs5("tmphess");
00274   ofs << nvar;
00275   {
00276     {
00277       dvariable vf=0.0;
00278       vf=initial_params::reset(dvar_vector(x));
00279       *objective_function_value::pobjfun=0.0;
00280       pre_userfunction();
00281       vf+=*objective_function_value::pobjfun;
00282       f=value(vf);
00283       gradcalc(nvar,g1);
00284     }
00285     for (int i=1;i<=nvar;i++)
00286     {
00287       cout << "Estimating row " << i << " out of " << nvar
00288            << " for hessian" << endl;
00289 
00290       double f=0.0;
00291       double xsave=x(i);
00292       double sdelta1=x(i)+delta;
00293       sdelta1-=x(i);
00294       x(i)=xsave+sdelta1;
00295       dvariable vf=0.0;
00296       vf=initial_params::reset(dvar_vector(x));
00297       *objective_function_value::pobjfun=0.0;
00298       pre_userfunction();
00299       vf+=*objective_function_value::pobjfun;
00300       f=value(vf);
00301       gradcalc(nvar,g1);
00302 
00303       double sdelta2=x(i)-delta;
00304       sdelta2-=x(i);
00305       x(i)=xsave+sdelta2;
00306       vf=0.0;
00307       vf=initial_params::reset(dvar_vector(x));
00308       *objective_function_value::pobjfun=0.0;
00309       pre_userfunction();
00310       vf+=*objective_function_value::pobjfun;
00311       f=value(vf);
00312       gradcalc(nvar,g2);
00313       x(i)=xsave;
00314       hess1=(g1-g2)/(sdelta1-sdelta2);
00315 
00316       sdelta1=x(i)+eps*delta;
00317       sdelta1-=x(i);
00318       x(i)=xsave+sdelta1;
00319       vf=0.0;
00320       vf=initial_params::reset(dvar_vector(x));
00321       *objective_function_value::pobjfun=0.0;
00322       pre_userfunction();
00323       vf+=*objective_function_value::pobjfun;
00324       f=value(vf);
00325       gradcalc(nvar,g1);
00326 
00327       x(i)=xsave-eps*delta;
00328       sdelta2=x(i)-eps*delta;
00329       sdelta2-=x(i);
00330       x(i)=xsave+sdelta2;
00331       vf=0.0;
00332       vf=0.0;
00333       vf=initial_params::reset(dvar_vector(x));
00334       *objective_function_value::pobjfun=0.0;
00335       pre_userfunction();
00336       vf+=*objective_function_value::pobjfun;
00337       f=value(vf);
00338       gradcalc(nvar,g2);
00339       x(i)=xsave;
00340 
00341       double eps2=eps*eps;
00342       hess2=(g1-g2)/(sdelta1-sdelta2);
00343       hess=(eps2*hess1-hess2) /(eps2-1.);
00344       ofs << hess;
00345     }
00346     // do the hessian for the constraint function
00347     uostream cofs("constrnt.hes");
00348     cofs << nvar;
00349     for (i=1;i<=nvar;i++)
00350     {
00351       cout << "Estimating row " << i << " out of " << nvar
00352            << " for hessian" << endl;
00353 
00354       double f=0.0;
00355       double xsave=x(i);
00356       sdelta1=x(i)+delta;
00357       sdelta1-=x(i);
00358       x(i)=xsave+sdelta1;
00359       dvariable vf=0.0;
00360       vf=initial_params::reset(dvar_vector(x));
00361       *objective_function_value::pobjfun=0.0;
00362       pre_userfunction();
00363       vf=likeprof_params::likeprofptr[iprof]->variable();
00364       f=value(vf);
00365       gradcalc(nvar,g1);
00366 
00367       sdelta2=x(i)-delta;
00368       sdelta2-=x(i);
00369       x(i)=xsave+sdelta2;
00370       vf=0.0;
00371       vf=initial_params::reset(dvar_vector(x));
00372       *objective_function_value::pobjfun=0.0;
00373       pre_userfunction();
00374       vf=likeprof_params::likeprofptr[iprof]->variable();
00375       f=value(vf);
00376       gradcalc(nvar,g2);
00377       x(i)=xsave;
00378       hess1=(g1-g2)/(sdelta1-sdelta2);
00379 
00380       sdelta1=x(i)+eps*delta;
00381       sdelta1-=x(i);
00382       x(i)=xsave+sdelta1;
00383       vf=0.0;
00384       vf=initial_params::reset(dvar_vector(x));
00385       *objective_function_value::pobjfun=0.0;
00386       pre_userfunction();
00387       vf=likeprof_params::likeprofptr[iprof]->variable();
00388       f=value(vf);
00389       gradcalc(nvar,g1);
00390 
00391       x(i)=xsave-eps*delta;
00392       sdelta2=x(i)-eps*delta;
00393       sdelta2-=x(i);
00394       x(i)=xsave+sdelta2;
00395       vf=0.0;
00396       vf=initial_params::reset(dvar_vector(x));
00397       *objective_function_value::pobjfun=0.0;
00398       pre_userfunction();
00399       vf=likeprof_params::likeprofptr[iprof]->variable();
00400       f=value(vf);
00401       gradcalc(nvar,g2);
00402       x(i)=xsave;
00403 
00404       double eps2=eps*eps;
00405       hess2=(g1-g2)/(sdelta1-sdelta2);
00406       hess=(eps2*hess1-hess2) /(eps2-1.);
00407       cofs << hess;
00408     }
00409   }
00410   gradient_structure::set_NO_DERIVATIVES();
00411 }
00412 */
00413 
00418 void function_minimizer::depvars_routine(void)
00419 {
00420   reset_gradient_stack();
00421   dvector ggg(1,1);
00422   gradcalc(0,ggg);
00423   gradient_structure::set_YES_DERIVATIVES();
00424   initial_params::restore_start_phase();
00425   if (lapprox && lapprox->no_re_ders_flag)
00426   {
00427     initial_params::set_inactive_random_effects();
00428   }
00429   int nvar=initial_params::nvarcalc(); // get the number of active parameters
00430   int ndvar=stddev_params::num_stddev_calc();
00431   independent_variables x(1,nvar);
00432   initial_params::xinit(x);        // get the initial values into the x vector
00433   //double f;
00434   //double delta=1.e-7;
00435   adstring tmpstring="admodel.dep";
00436   if (ad_comm::wd_flag)
00437      tmpstring = ad_comm::adprogram_name + ".dep";
00438   ofstream ofs((char*)tmpstring);
00439   if (lapprox)
00440   {
00441     lapprox->no_function_component_flag=1;
00442   }
00443 
00444   dvariable vf;
00445   vf=initial_params::reset(dvar_vector(x));
00446   *objective_function_value::pobjfun=0.0;
00447   pre_userfunction();
00448   vf+=*objective_function_value::pobjfun;
00449 
00450   ofs << nvar << "  "  << ndvar << endl;
00451   int i;
00452   for (i=0;i< stddev_params::num_stddev_params;i++)
00453   {
00454      stddev_params::stddevptr[i]->set_dependent_variables();
00455   }
00456   gradient_structure::jacobcalc(nvar,ofs);
00457   for (i=0;i< stddev_params::num_stddev_params;i++)
00458   {
00459      ofs << stddev_params::stddevptr[i]->label() << "  ";
00460      ofs << stddev_params::stddevptr[i]->size_count() << endl;
00461   }
00462   if (lapprox)
00463   {
00464     lapprox->no_function_component_flag=0;
00465   }
00466   gradient_structure::set_NO_DERIVATIVES();
00467 }
00471 void function_minimizer::hess_inv(void)
00472 {
00473   initial_params::set_inactive_only_random_effects();
00474   int nvar=initial_params::nvarcalc(); // get the number of active parameters
00475   independent_variables x(1,nvar);
00476 
00477   initial_params::xinit(x);        // get the initial values into the x vector
00478   //double f;
00479   dmatrix hess(1,nvar,1,nvar);
00480   uistream ifs("admodel.hes");
00481   int file_nvar = 0;
00482   ifs >> file_nvar;
00483   if (nvar != file_nvar)
00484   {
00485     cerr << "Number of active variables in file mod_hess.rpt is wrong"
00486          << endl;
00487   }
00488 
00489   for (int i = 1;i <= nvar; i++)
00490   {
00491     ifs >> hess(i);
00492     if (!ifs)
00493     {
00494       cerr << "Error reading line " << i  << " of the hessian"
00495            << " in routine hess_inv()" << endl;
00496       exit(1);
00497     }
00498   }
00499   int hybflag = 0;
00500   ifs >> hybflag;
00501   dvector sscale(1,nvar);
00502   ifs >> sscale;
00503   if (!ifs)
00504   {
00505     cerr << "Error reading sscale"
00506          << " in routine hess_inv()" << endl;
00507   }
00508 
00509   double maxerr=0.0;
00510   for (int i = 1;i <= nvar; i++)
00511   {
00512     for (int j=1;j<i;j++)
00513     {
00514       double tmp=(hess(i,j)+hess(j,i))/2.;
00515       double tmp1=fabs(hess(i,j)-hess(j,i));
00516       tmp1/=(1.e-4+fabs(hess(i,j))+fabs(hess(j,i)));
00517       if (tmp1>maxerr) maxerr=tmp1;
00518       hess(i,j)=tmp;
00519       hess(j,i)=tmp;
00520     }
00521   }
00522   /*
00523   if (maxerr>1.e-2)
00524   {
00525     cerr << "warning -- hessian aprroximation is poor" << endl;
00526   }
00527  */
00528 
00529   for (int i = 1;i <= nvar; i++)
00530   {
00531     int zero_switch=0;
00532     for (int j=1;j<=nvar;j++)
00533     {
00534       if (hess(i,j)!=0.0)
00535       {
00536         zero_switch=1;
00537       }
00538     }
00539     if (!zero_switch)
00540     {
00541       cerr << " Hessian is 0 in row " << i << endl;
00542       cerr << " This means that the derivative if probably identically 0 "
00543               " for this parameter" << endl;
00544     }
00545   }
00546 
00547   int ssggnn;
00548   ln_det(hess,ssggnn);
00549   int on1=0;
00550   {
00551     ofstream ofs3((char*)(ad_comm::adprogram_name + adstring(".eva")));
00552     {
00553       dvector se=eigenvalues(hess);
00554       ofs3 << setshowpoint() << setw(14) << setprecision(10)
00555            << "unsorted:\t" << se << endl;
00556      se=sort(se);
00557      ofs3 << setshowpoint() << setw(14) << setprecision(10)
00558      << "sorted:\t" << se << endl;
00559      if (se(se.indexmin())<=0.0)
00560       {
00561         negative_eigenvalue_flag=1;
00562         cout << "Warning -- Hessian does not appear to be"
00563          " positive definite" << endl;
00564       }
00565     }
00566     ivector negflags(0,hess.indexmax());
00567     int num_negflags=0;
00568     {
00569       int on = option_match(ad_comm::argc,ad_comm::argv,"-eigvec");
00570       on1=option_match(ad_comm::argc,ad_comm::argv,"-spmin");
00571       if (on > -1 || on1 >-1 )
00572       {
00573         ofs3 << setshowpoint() << setw(14) << setprecision(10)
00574           << eigenvalues(hess) << endl;
00575         dmatrix ev=trans(eigenvectors(hess));
00576         ofs3 << setshowpoint() << setw(14) << setprecision(10)
00577           << ev << endl;
00578         for (int i=1;i<=ev.indexmax();i++)
00579         {
00580           double lam=ev(i)*hess*ev(i);
00581           ofs3 << setshowpoint() << setw(14) << setprecision(10)
00582             << lam << "  "  << ev(i)*ev(i) << endl;
00583           if (lam<0.0)
00584           {
00585             num_negflags++;
00586             negflags(num_negflags)=i;
00587           }
00588         }
00589         if ( (on1>-1) && (num_negflags>0))   // we will try to get away from
00590         {                                     // saddle point
00591           negative_eigenvalue_flag=0;
00592           spminflag=1;
00593           if(negdirections)
00594           {
00595             delete negdirections;
00596           }
00597           negdirections = new dmatrix(1,num_negflags);
00598           for (int i=1;i<=num_negflags;i++)
00599           {
00600             (*negdirections)(i)=ev(negflags(i));
00601           }
00602         }
00603         int on2 = option_match(ad_comm::argc,ad_comm::argv,"-cross");
00604         if (on2>-1)
00605         {                                     // saddle point
00606           dmatrix cross(1,ev.indexmax(),1,ev.indexmax());
00607           for (int i = 1;i <= ev.indexmax(); i++)
00608           {
00609             for (int j=1;j<=ev.indexmax();j++)
00610             {
00611               cross(i,j)=ev(i)*ev(j);
00612             }
00613           }
00614           ofs3 <<  endl << "  e(i)*e(j) ";
00615           ofs3 << endl << cross << endl;
00616         }
00617       }
00618     }
00619 
00620     if (spminflag==0)
00621     {
00622       if (num_negflags==0)
00623       {
00624         hess=inv(hess);
00625         int on=0;
00626         if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-eigvec"))>-1)
00627         {
00628           int i;
00629           ofs3 << "choleski decomp of correlation" << endl;
00630           dmatrix ch=choleski_decomp(hess);
00631           for (i=1;i<=ch.indexmax();i++)
00632             ofs3 << ch(i)/norm(ch(i)) << endl;
00633           ofs3 << "parameterization of choleski decomnp of correlation" << endl;
00634           for (i=1;i<=ch.indexmax();i++)
00635           {
00636             dvector tmp=ch(i)/norm(ch(i));
00637             ofs3 << tmp(1,i)/tmp(i) << endl;
00638           }
00639         }
00640       }
00641     }
00642   }
00643   if (spminflag==0)
00644   {
00645     if (on1<0)
00646     {
00647       for (int i = 1;i <= nvar; i++)
00648       {
00649         if (hess(i,i) <= 0.0)
00650         {
00651           hess_errorreport();
00652           ad_exit(1);
00653         }
00654       }
00655     }
00656     {
00657       adstring tmpstring="admodel.cov";
00658       if (ad_comm::wd_flag)
00659         tmpstring = ad_comm::adprogram_name + ".cov";
00660       uostream ofs((char*)tmpstring);
00661       ofs << nvar << hess;
00662       ofs << gradient_structure::Hybrid_bounded_flag;
00663       ofs << sscale;
00664     }
00665   }
00666 }
00667 void hess_calcreport(int i,int nvar)
00668 {
00669   if (ad_printf)
00670     (*ad_printf)("Estimating row %d out of %d for hessian\n",i,nvar);
00671   else
00672     cout << "Estimating row " << i << " out of " << nvar << " for hessian"
00673          << endl;
00674 }
00675 void hess_errorreport(void)
00676 {
00677   if (ad_printf)
00678     (*ad_printf)("Hessian does not appear to be positive definite\n");
00679   else
00680     cout << "Hessian does not appear to be positive definite\n" << endl;
00681 }