/** * Description not yet available. * \param */ double ghk(const dvector& lower,const dvector& upper,const dmatrix& Sigma, const dmatrix& eps) { int m=eps.indexmax(); int n=lower.indexmax(); double ssum=0.0; dmatrix ch=choleski_decomp(Sigma); dvector l(1,n); dvector u(1,n); for (int k=1;k<=m;k++) { double weight=1.0; l=lower; u=upper; for (int j=1;j<=n;j++) { l(j)/=ch(j,j); u(j)/=ch(j,j); double Phiu=cumd_norm(u(j)); double Phil=cumd_norm(l(j)); weight*=Phiu-Phil; double eta=inv_cumd_norm((Phiu-Phil)*eps(k,j)+Phil); for (int i=j+1;i<=n;i++) { double tmp=ch(i,j)*eta; l(i)-=tmp; u(i)-=tmp; } } ssum+=weight; } return ssum/m; }
/** * Description not yet available. * \param */ double ghk(const dvector& lower,const dvector& upper,const dmatrix& Sigma, const dmatrix& eps,int _i) { int n=lower.indexmax(); dmatrix ch=choleski_decomp(Sigma); dvector l(1,n); dvector u(1,n); ghk_test(eps,_i); // test for valid i range double weight=1.0; int k=_i; { l=lower; u=upper; for (int j=1;j<=n;j++) { l(j)/=ch(j,j); u(j)/=ch(j,j); double Phiu=cumd_norm(u(j)); double Phil=cumd_norm(l(j)); weight*=Phiu-Phil; double eta=inv_cumd_norm((Phiu-Phil)*eps(k,j)+Phil); for (int i=j+1;i<=n;i++) { double tmp=ch(i,j)*eta; l(i)-=tmp; u(i)-=tmp; } } } return weight; }
/** * Description not yet available. * \param */ dvariable ghk_m(const dvar_vector& upper,const dvar_matrix& Sigma, const dmatrix& eps) { RETURN_ARRAYS_INCREMENT(); int n=upper.indexmax(); int m=eps.indexmax(); dvariable ssum=0.0; dvar_vector u(1,n); dvar_matrix ch=choleski_decomp(Sigma); for (int k=1;k<=m;k++) { dvariable weight=1.0; u=upper; for (int j=1;j<=n;j++) { u(j)/=ch(j,j); dvariable Phiu=cumd_norm(u(j)); weight*=Phiu; dvariable eta=inv_cumd_norm(1.e-30+Phiu*eps(k,j)); for (int i=j+1;i<=n;i++) { dvariable tmp=ch(i,j)*eta; u(i)-=tmp; } } ssum+=weight; } RETURN_ARRAYS_DECREMENT(); return ssum/m; }
/** * Description not yet available. * \param */ dvariable ghk(const dvar_vector& lower,const dvar_vector& upper, const dvar_matrix& Sigma, const dmatrix& eps,int _i) { RETURN_ARRAYS_INCREMENT(); int n=lower.indexmax(); dvar_matrix ch=choleski_decomp(Sigma); dvar_vector l(1,n); dvar_vector u(1,n); ghk_test(eps,_i); dvariable weight=1.0; int k=_i; { l=lower; u=upper; for (int j=1;j<=n;j++) { l(j)/=ch(j,j); u(j)/=ch(j,j); dvariable Phiu=cumd_norm(u(j)); dvariable Phil=cumd_norm(l(j)); weight*=Phiu-Phil; dvariable eta=inv_cumd_norm((Phiu-Phil)*eps(k,j)+Phil+1.e-30); for (int i=j+1;i<=n;i++) { dvariable tmp=ch(i,j)*eta; l(i)-=tmp; u(i)-=tmp; } } } RETURN_ARRAYS_DECREMENT(); return weight; }
/** * Description not yet available. * \param */ void laplace_approximation_calculator:: do_newton_raphson_banded(function_minimizer * pfmin,double f_from_1, int& no_converge_flag) { //quadratic_prior * tmpptr=quadratic_prior::ptr[0]; //cout << tmpptr << endl; laplace_approximation_calculator::where_are_we_flag=2; double maxg=fabs(evaluate_function(uhat,pfmin)); laplace_approximation_calculator::where_are_we_flag=0; dvector uhat_old(1,usize); for(int ii=1;ii<=num_nr_iters;ii++) { // test newton raphson switch(hesstype) { case 3: bHess->initialize(); break; case 4: Hess.initialize(); break; default: cerr << "Illegal value for hesstype here" << endl; ad_exit(1); } grad.initialize(); //int check=initial_params::stddev_scale(scale,uhat); //check=initial_params::stddev_curvscale(curv,uhat); //max_separable_g=0.0; sparse_count = 0; step=get_newton_raphson_info_banded(pfmin); //if (bHess) // cout << "norm(*bHess) = " << norm(*bHess) << endl; //cout << "norm(Hess) = " << norm(Hess) << endl; //cout << grad << endl; //check_pool_depths(); if (!initial_params::mc_phase) cout << "Newton raphson " << ii << " "; if (quadratic_prior::get_num_quadratic_prior()>0) { quadratic_prior::get_cHessian_contribution(Hess,xsize); quadratic_prior::get_cgradient_contribution(grad,xsize); } int ierr=0; if (hesstype==3) { if (use_dd_nr==0) { banded_lower_triangular_dmatrix bltd=choleski_decomp(*bHess,ierr); if (ierr && no_converge_flag ==0) { no_converge_flag=1; //break; } if (ierr) { double oldval; evaluate_function(oldval,uhat,pfmin); uhat=banded_calculations_trust_region_approach(uhat,pfmin); } else { if (dd_nr_flag==0) { dvector v=solve(bltd,grad); step=-solve_trans(bltd,v); //uhat_old=uhat; uhat+=step; } else { #if defined(USE_DD_STUFF) int n=grad.indexmax(); maxg=fabs(evaluate_function(uhat,pfmin)); uhat=dd_newton_raphson2(grad,*bHess,uhat); #else cerr << "high precision Newton Raphson not implemented" << endl; ad_exit(1); #endif } maxg=fabs(evaluate_function(uhat,pfmin)); if (f_from_1< pfmin->lapprox->fmc1.fbest) { uhat=banded_calculations_trust_region_approach(uhat,pfmin); maxg=fabs(evaluate_function(uhat,pfmin)); } } } else { cout << "error not used" << endl; ad_exit(1); /* banded_symmetric_ddmatrix bHessdd=banded_symmetric_ddmatrix(*bHess); ddvector gradd=ddvector(grad); //banded_lower_triangular_ddmatrix bltdd=choleski_decomp(bHessdd,ierr); if (ierr && no_converge_flag ==0) { no_converge_flag=1; break; } if (ierr) { double oldval; evaluate_function(oldval,uhat,pfmin); uhat=banded_calculations_trust_region_approach(uhat,pfmin); maxg=fabs(evaluate_function(uhat,pfmin)); } else { ddvector v=solve(bHessdd,gradd); step=-make_dvector(v); //uhat_old=uhat; uhat=make_dvector(ddvector(uhat)+step); maxg=fabs(evaluate_function(uhat,pfmin)); if (f_from_1< pfmin->lapprox->fmc1.fbest) { uhat=banded_calculations_trust_region_approach(uhat,pfmin); maxg=fabs(evaluate_function(uhat,pfmin)); } } */ } if (maxg < 1.e-13) { break; } } else if (hesstype==4) { dvector step; # if defined(USE_ATLAS) if (!ad_comm::no_atlas_flag) { step=-atlas_solve_spd(Hess,grad,ierr); } else { dmatrix A=choleski_decomp_positive(Hess,ierr); if (!ierr) { step=-solve(Hess,grad); //step=-solve(A*trans(A),grad); } } if (!ierr) break; # else if (sparse_hessian_flag) { //step=-solve(*sparse_triplet,Hess,grad,*sparse_symbolic); dvector temp=solve(*sparse_triplet2,grad,*sparse_symbolic2,ierr); if (ierr) { step=-temp; } else { cerr << "matrix not pos definite in sparse choleski" << endl; pfmin->bad_step_flag=1; int on; int nopt; if ((on=option_match(ad_comm::argc,ad_comm::argv,"-ieigvec",nopt)) >-1) { dmatrix M=make_dmatrix(*sparse_triplet2); ofstream ofs3("inner-eigvectors"); ofs3 << "eigenvalues and eigenvectors " << endl; dvector v=eigenvalues(M); dmatrix ev=trans(eigenvectors(M)); ofs3 << "eigenvectors" << endl; int i; for (i=1;i<=ev.indexmax();i++) { ofs3 << setw(4) << i << " " << setshowpoint() << setw(14) << setprecision(10) << v(i) << " " << setshowpoint() << setw(14) << setprecision(10) << ev(i) << endl; } } } //cout << norm2(step-tmpstep) << endl; //dvector step1=-solve(Hess,grad); //cout << norm2(step-step1) << endl; } else { step=-solve(Hess,grad); } # endif if (pmin->bad_step_flag) break; uhat_old=uhat; uhat+=step; double maxg_old=maxg; maxg=fabs(evaluate_function(uhat,pfmin)); if (maxg>maxg_old) { uhat=uhat_old; evaluate_function(uhat,pfmin); break; } if (maxg < 1.e-13) { break; } } if (sparse_hessian_flag==0) { for (int i=1;i<=usize;i++) { y(i+xsize)=uhat(i); } } else { for (int i=1;i<=usize;i++) { value(y(i+xsize))=uhat(i); } } } }
/** Symmetrize and invert the hessian */ void function_minimizer::hess_inv(void) { initial_params::set_inactive_only_random_effects(); int nvar=initial_params::nvarcalc(); // get the number of active parameters independent_variables x(1,nvar); initial_params::xinit(x); // get the initial values into the x vector //double f; dmatrix hess(1,nvar,1,nvar); uistream ifs("admodel.hes"); int file_nvar = 0; ifs >> file_nvar; if (nvar != file_nvar) { cerr << "Number of active variables in file mod_hess.rpt is wrong" << endl; } for (int i = 1;i <= nvar; i++) { ifs >> hess(i); if (!ifs) { cerr << "Error reading line " << i << " of the hessian" << " in routine hess_inv()" << endl; exit(1); } } int hybflag = 0; ifs >> hybflag; dvector sscale(1,nvar); ifs >> sscale; if (!ifs) { cerr << "Error reading sscale" << " in routine hess_inv()" << endl; } double maxerr=0.0; for (int i = 1;i <= nvar; i++) { for (int j=1;j<i;j++) { double tmp=(hess(i,j)+hess(j,i))/2.; double tmp1=fabs(hess(i,j)-hess(j,i)); tmp1/=(1.e-4+fabs(hess(i,j))+fabs(hess(j,i))); if (tmp1>maxerr) maxerr=tmp1; hess(i,j)=tmp; hess(j,i)=tmp; } } /* if (maxerr>1.e-2) { cerr << "warning -- hessian aprroximation is poor" << endl; } */ for (int i = 1;i <= nvar; i++) { int zero_switch=0; for (int j=1;j<=nvar;j++) { if (hess(i,j)!=0.0) { zero_switch=1; } } if (!zero_switch) { cerr << " Hessian is 0 in row " << i << endl; cerr << " This means that the derivative if probably identically 0 " " for this parameter" << endl; } } int ssggnn; ln_det(hess,ssggnn); int on1=0; { ofstream ofs3((char*)(ad_comm::adprogram_name + adstring(".eva"))); { dvector se=eigenvalues(hess); ofs3 << setshowpoint() << setw(14) << setprecision(10) << "unsorted:\t" << se << endl; se=sort(se); ofs3 << setshowpoint() << setw(14) << setprecision(10) << "sorted:\t" << se << endl; if (se(se.indexmin())<=0.0) { negative_eigenvalue_flag=1; cout << "Warning -- Hessian does not appear to be" " positive definite" << endl; } } ivector negflags(0,hess.indexmax()); int num_negflags=0; { int on = option_match(ad_comm::argc,ad_comm::argv,"-eigvec"); on1=option_match(ad_comm::argc,ad_comm::argv,"-spmin"); if (on > -1 || on1 >-1 ) { ofs3 << setshowpoint() << setw(14) << setprecision(10) << eigenvalues(hess) << endl; dmatrix ev=trans(eigenvectors(hess)); ofs3 << setshowpoint() << setw(14) << setprecision(10) << ev << endl; for (int i=1;i<=ev.indexmax();i++) { double lam=ev(i)*hess*ev(i); ofs3 << setshowpoint() << setw(14) << setprecision(10) << lam << " " << ev(i)*ev(i) << endl; if (lam<0.0) { num_negflags++; negflags(num_negflags)=i; } } if ( (on1>-1) && (num_negflags>0)) // we will try to get away from { // saddle point negative_eigenvalue_flag=0; spminflag=1; if(negdirections) { delete negdirections; } negdirections = new dmatrix(1,num_negflags); for (int i=1;i<=num_negflags;i++) { (*negdirections)(i)=ev(negflags(i)); } } int on2 = option_match(ad_comm::argc,ad_comm::argv,"-cross"); if (on2>-1) { // saddle point dmatrix cross(1,ev.indexmax(),1,ev.indexmax()); for (int i = 1;i <= ev.indexmax(); i++) { for (int j=1;j<=ev.indexmax();j++) { cross(i,j)=ev(i)*ev(j); } } ofs3 << endl << " e(i)*e(j) "; ofs3 << endl << cross << endl; } } } if (spminflag==0) { if (num_negflags==0) { hess=inv(hess); int on=0; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-eigvec"))>-1) { int i; ofs3 << "choleski decomp of correlation" << endl; dmatrix ch=choleski_decomp(hess); for (i=1;i<=ch.indexmax();i++) ofs3 << ch(i)/norm(ch(i)) << endl; ofs3 << "parameterization of choleski decomnp of correlation" << endl; for (i=1;i<=ch.indexmax();i++) { dvector tmp=ch(i)/norm(ch(i)); ofs3 << tmp(1,i)/tmp(i) << endl; } } } } } if (spminflag==0) { if (on1<0) { for (int i = 1;i <= nvar; i++) { if (hess(i,i) <= 0.0) { hess_errorreport(); ad_exit(1); } } } { adstring tmpstring="admodel.cov"; if (ad_comm::wd_flag) tmpstring = ad_comm::adprogram_name + ".cov"; uostream ofs((char*)tmpstring); ofs << nvar << hess; ofs << gradient_structure::Hybrid_bounded_flag; ofs << sscale; } } }
void function_minimizer::shmc_mcmc_routine(int nmcmc,int iseed0,double dscale, int restart_flag) { if (nmcmc<=0) { cerr << endl << "Error: Negative iterations for MCMC not meaningful" << endl; ad_exit(1); } uostream * pofs_psave=NULL; if (mcmc2_flag==1) { initial_params::restore_start_phase(); } initial_params::set_inactive_random_effects(); initial_params::set_active_random_effects(); int nvar_re=initial_params::nvarcalc(); int nvar=initial_params::nvarcalc(); // get the number of active parameters if (mcmc2_flag==0) { initial_params::set_inactive_random_effects(); nvar=initial_params::nvarcalc(); // get the number of active parameters } initial_params::restore_start_phase(); independent_variables parsave(1,nvar_re); // dvector x(1,nvar); // initial_params::xinit(x); // dvector pen_vector(1,nvar); // { // initial_params::reset(dvar_vector(x),pen_vector); // } initial_params::mc_phase=1; int old_Hybrid_bounded_flag=-1; int on,nopt = 0; //// ------------------------------ Parse input options // Step size. If not specified, will be adapted. If specified must be >0 // and will not be adapted. double eps=0.1; double _eps=-1.0; int useDA=1; // whether to adapt step size if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-hyeps",nopt))>-1) { if (!nopt) // not specified means to adapt, using function below to find reasonable one { cerr << "Warning: No step size given after -hyeps, ignoring" << endl; useDA=1; } else // read in specified value and do not adapt { istringstream ist(ad_comm::argv[on+1]); ist >> _eps; if (_eps<=0) { cerr << "Error: step size (-hyeps argument) needs positive number"; ad_exit(1); } else { eps=_eps; useDA=0; } } } // Chain number -- for console display purposes only int chain=1; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-chain",nopt))>-1) { if (nopt) { int iii=atoi(ad_comm::argv[on+1]); if (iii <1) { cerr << "Error: chain must be >= 1" << endl; ad_exit(1); } else { chain=iii; } } } // Number of leapfrog steps. Defaults to 10. int L=10; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-hynstep",nopt))>-1) { if (nopt) { int _L=atoi(ad_comm::argv[on+1]); if (_L < 1 ) { cerr << "Error: hynstep argument must be integer > 0 " << endl; ad_exit(1); } else { L=_L; } } } // Number of warmup samples if using adaptation of step size. Defaults to // half of iterations. int nwarmup= (int)nmcmc/2; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-nwarmup",nopt))>-1) { if (nopt) { int iii=atoi(ad_comm::argv[on+1]); if (iii <=0 || iii > nmcmc) { cerr << "Error: nwarmup must be 0 < nwarmup < nmcmc" << endl; ad_exit(1); } else { nwarmup=iii; } } } // Target acceptance rate for step size adaptation. Must be // 0<adapt_delta<1. Defaults to 0.8. double adapt_delta=0.8; // target acceptance rate specified by the user if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-adapt_delta",nopt))>-1) { if (nopt) { istringstream ist(ad_comm::argv[on+1]); double _adapt_delta; ist >> _adapt_delta; if (_adapt_delta < 0 || _adapt_delta > 1 ) { cerr << "Error: adapt_delta must be between 0 and 1" " using default of 0.8" << endl; } else { adapt_delta=_adapt_delta; } } } // Use diagnoal covariance (identity mass matrix) int diag_option=0; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcdiag"))>-1) { diag_option=1; cout << " Setting covariance matrix to diagonal with entries " << dscale << endl; } // Restart chain from previous run? int mcrestart_flag=option_match(ad_comm::argc,ad_comm::argv,"-mcr"); if(mcrestart_flag > -1){ cerr << endl << "Error: -mcr option not implemented for HMC" << endl; ad_exit(1); } if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcec"))>-1) { cerr << endl << "Error: -mcec option not yet implemented with HMC" << endl; ad_exit(1); // use_empirical_flag=1; // read_empirical_covariance_matrix(nvar,S,ad_comm::adprogram_name); } // Prepare the mass matrix for use. Depends on many factors below. dmatrix S(1,nvar,1,nvar); dvector old_scale(1,nvar); int old_nvar; // Need to grab old_scale values still, since it is scaled below read_covariance_matrix(S,nvar,old_Hybrid_bounded_flag,old_scale); if (diag_option) // set covariance to be diagonal { S.initialize(); for (int i=1;i<=nvar;i++) { S(i,i)=dscale; } } // How much to thin, for now fixed at 1. if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcsave"))>-1) { cerr << "Option -mcsave does not currently work with HMC -- every iteration is saved" << endl; ad_exit(1); } //// ------------------------------ End of input processing //// Setup more inputs and outputs pofs_psave= new uostream((char*)(ad_comm::adprogram_name + adstring(".psv"))); if (!pofs_psave|| !(*pofs_psave)) { cerr << "Error trying to open file" << ad_comm::adprogram_name + adstring(".psv") << endl; ad_exit(1); } if (mcrestart_flag == -1 ) { (*pofs_psave) << nvar; } // need to rescale the hessian // get the current scale dvector x0(1,nvar); dvector current_scale(1,nvar); initial_params::xinit(x0); int mctmp=initial_params::mc_phase; initial_params::mc_phase=0; initial_params::stddev_scale(current_scale,x0); initial_params::mc_phase=mctmp; // cout << "old scale=" << old_scale << endl; // cout << "current scale=" << current_scale << endl; // cout << "S before=" << S << endl; // I think this is only needed if mcmc2 is used?? // for (int i=1;i<=nvar;i++) // { // for (int j=1;j<=nvar;j++) // { // S(i,j)*=old_scale(i)*old_scale(j); // } // } if(diag_option){ for (int i=1;i<=nvar;i++) { for (int j=1;j<=nvar;j++) { S(i,j)*=current_scale(i)*current_scale(j); } } } // cout << "S after=" << S << endl; gradient_structure::set_NO_DERIVATIVES(); if (mcmc2_flag==0) { initial_params::set_inactive_random_effects(); } // Setup random number generator, based on seed passed int iseed=2197; if (iseed0) iseed=iseed0; random_number_generator rng(iseed); gradient_structure::set_YES_DERIVATIVES(); initial_params::xinit(x0); // Dual averaging components dvector epsvec(1,nmcmc+1), epsbar(1,nmcmc+1), Hbar(1,nmcmc+1); epsvec.initialize(); epsbar.initialize(); Hbar.initialize(); double time_warmup=0; double time_total=0; std::clock_t start = clock(); time_t now = time(0); tm* localtm = localtime(&now); cout << endl << "Starting static HMC for model '" << ad_comm::adprogram_name << "' at " << asctime(localtm); // write sampler parameters ofstream adaptation("adaptation.csv", ios::trunc); adaptation << "accept_stat__,stepsize__,int_time__,energy__,lp__" << endl; // Declare and initialize the variables needed for the algorithm dmatrix chd = choleski_decomp(S); // cholesky decomp of mass matrix dvector y(1,nvar); // unbounded parameters y.initialize(); // transformed params independent_variables z(1,nvar); z=chd*y; dvector gr(1,nvar); // gradients in unbounded space // Need to run this to fill gr with current gradients and initial NLL. double nllbegin=get_hybrid_monte_carlo_value(nvar,z,gr); if(std::isnan(nllbegin)){ cerr << "Starting MCMC trajectory at NaN -- something is wrong!" << endl; ad_exit(1); } // initial rotated gradient dvector gr2(1,nvar); gr2=gr*chd; dvector p(1,nvar); // momentum vector p.fill_randn(rng); // Copy initial value to parsave in case first trajectory rejected initial_params::copy_all_values(parsave,1.0); double iaccept=0.0; // The gradient and params at beginning of trajectory, in case rejected. dvector gr2begin(1,nvar); gr2begin=gr2; dvector ybegin(1,nvar); ybegin=y; double nll=nllbegin; // if(useDA){ // eps=find_reasonable_stepsize(nvar,y,p,chd); // epsvec(1)=eps; epsbar(1)=eps; Hbar(1)=0; // } double mu=log(10*eps); // Start of MCMC chain for (int is=1;is<=nmcmc;is++) { // Random momentum for next iteration, only affects Ham values p.fill_randn(rng); double H0=nll+0.5*norm2(p); // Generate trajectory int divergence=0; for (int i=1;i<=L;i++) { // leapfrog updates gr, p, y, and gr2 by reference nll=leapfrog(nvar, gr, chd, eps, p, y, gr2); // Break trajectory early if a divergence occurs to save computation if(std::isnan(nll)){ divergence=1; break; } } // end of trajectory // Test whether to accept the proposed state double Ham=nll+0.5*norm2(p); // Update Hamiltonian for proposed set double alpha=min(1.0, exp(H0-Ham)); // acceptance ratio double rr=randu(rng); // Runif(1) if (rr<alpha && !divergence){ // accept iaccept++; // Update for next iteration: params, Hamiltonian and gr2 ybegin=y; gr2begin=gr2; nllbegin=nll; initial_params::copy_all_values(parsave,1.0); } else { // Reject and don't update anything to reuse initials for next trajectory y=ybegin; gr2=gr2begin; nll=nllbegin; } // Save parameters to psv file, duplicated if rejected (*pofs_psave) << parsave; // Adaptation of step size (eps). if(useDA && is <= nwarmup){ eps=adapt_eps(is, eps, alpha, adapt_delta, mu, epsvec, epsbar, Hbar); } adaptation << alpha << "," << eps << "," << eps*L << "," << H0 << "," << -nll << endl; if(is ==nwarmup) time_warmup = ( std::clock()-start)/(double) CLOCKS_PER_SEC; print_mcmc_progress(is, nmcmc, nwarmup, chain); } // end of MCMC chain // This final ratio should closely match adapt_delta if(useDA){ cout << "Final acceptance ratio=" << iaccept/nmcmc << " and target is " << adapt_delta<<endl; cout << "Final step size=" << eps << "; after " << nwarmup << " warmup iterations"<< endl; } else { cout << "Final acceptance ratio=" << iaccept/nmcmc << endl; } time_total = ( std::clock() - start ) / (double) CLOCKS_PER_SEC; print_mcmc_timing(time_warmup, time_total); // I assume this closes the connection to the file?? if (pofs_psave) { delete pofs_psave; pofs_psave=NULL; } } // end of HMC function
/** * Description not yet available. * \param */ double do_gauss_hermite_block_diagonal_multi(const dvector& x, const dvector& u0,const dmatrix& Hess,const dvector& _xadjoint, const dvector& _uadjoint,const dmatrix& _Hessadjoint, function_minimizer * pmin) { ADUNCONST(dvector,xadjoint) ADUNCONST(dvector,uadjoint) //ADUNCONST(dmatrix,Hessadjoint) dvector & w= *(pmin->multinomial_weights); const int xs=x.size(); const int us=u0.size(); gradient_structure::set_NO_DERIVATIVES(); int nsc=pmin->lapprox->num_separable_calls; const ivector lrea = (*pmin->lapprox->num_local_re_array)(1,nsc); int hroom = sum(square(lrea)); int nvar=x.size()+u0.size()+hroom; independent_variables y(1,nvar); // need to set random effects active together with whatever // init parameters should be active in this phase initial_params::set_inactive_only_random_effects(); initial_params::set_active_random_effects(); /*int onvar=*/initial_params::nvarcalc(); initial_params::xinit(y); // get the initial values into the // do we need this next line? y(1,xs)=x; int i,j; // contribution for quadratic prior if (quadratic_prior::get_num_quadratic_prior()>0) { //Hess+=quadratic_prior::get_cHessian_contribution(); int & vxs = (int&)(xs); quadratic_prior::get_cHessian_contribution(Hess,vxs); } // Here need hooks for sparse matrix structures dvar3_array & block_diagonal_vhessian= *pmin->lapprox->block_diagonal_vhessian; block_diagonal_vhessian.initialize(); dvar3_array& block_diagonal_ch= *pmin->lapprox->block_diagonal_vch; //dvar3_array(*pmin->lapprox->block_diagonal_ch); int ii=xs+us+1; d3_array& bdH=(*pmin->lapprox->block_diagonal_hessian); int ic; for (ic=1;ic<=nsc;ic++) { int lus=lrea(ic); for (i=1;i<=lus;i++) for (j=1;j<=lus;j++) y(ii++)=bdH(ic)(i,j); } dvector g(1,nvar); gradcalc(0,g); gradient_structure::set_YES_DERIVATIVES(); dvar_vector vy=dvar_vector(y); //initial_params::stddev_vscale(d,vy); ii=xs+us+1; if (initial_df1b2params::have_bounded_random_effects) { cerr << "can't do importance sampling with bounded random effects" " at present" << endl; ad_exit(1); } else { for (int ic=1;ic<=nsc;ic++) { int lus=lrea(ic); if (lus>0) { for (i=1;i<=lus;i++) { for (j=1;j<=lus;j++) { block_diagonal_vhessian(ic,i,j)=vy(ii++); } } block_diagonal_ch(ic)= choleski_decomp(inv(block_diagonal_vhessian(ic))); } } } int nsamp=pmin->lapprox->use_gauss_hermite; pmin->lapprox->in_gauss_hermite_phase=1; dvar_vector sample_value(1,nsamp); sample_value.initialize(); dvar_vector tau(1,us);; // !!! This only works for one random efect in each separable call // at present. if (pmin->lapprox->gh->mi) { delete pmin->lapprox->gh->mi; pmin->lapprox->gh->mi=0; } pmin->lapprox->gh->mi=new multi_index(1,nsamp, pmin->lapprox->multi_random_effects); multi_index & mi = *(pmin->lapprox->gh->mi); //for (int is=1;is<=nsamp;is++) dvector& xx=pmin->lapprox->gh->x; do { int offset=0; pmin->lapprox->num_separable_calls=0; //pmin->lapprox->gh->is=is; for (ic=1;ic<=nsc;ic++) { int lus=lrea(ic); // will need vector stuff here when more than one random effect if (lus>0) { //tau(offset+1,offset+lus).shift(1)=block_diagonal_ch(ic)(1,1)* // pmin->lapprox->gh->x(is); dvector xv(1,lus); for (int iu=1;iu<=lus;iu++) { xv(iu)= xx(mi()(iu)); } tau(offset+1,offset+lus).shift(1)=block_diagonal_ch(ic)*xv; offset+=lus; } } // have to reorder the terms to match the block diagonal hessian imatrix & ls=*(pmin->lapprox->block_diagonal_re_list); int mmin=ls.indexmin(); int mmax=ls.indexmax(); int ii=1; int i; for (i=mmin;i<=mmax;i++) { int cmin=ls(i).indexmin(); int cmax=ls(i).indexmax(); for (int j=cmin;j<=cmax;j++) { vy(ls(i,j))+=tau(ii++); } } if (ii-1 != us) { cerr << "error in interface" << endl; ad_exit(1); } initial_params::reset(vy); // get the values into the model ii=1; for (i=mmin;i<=mmax;i++) { int cmin=ls(i).indexmin(); int cmax=ls(i).indexmax(); for (int j=cmin;j<=cmax;j++) { vy(ls(i,j))-=tau(ii++); } } *objective_function_value::pobjfun=0.0; pmin->AD_uf_outer(); ++mi; } while(mi.get_depth()<=pmin->lapprox->multi_random_effects); nsc=pmin->lapprox->num_separable_calls; dvariable vf=pmin->do_gauss_hermite_integration(); int sgn=0; dvariable ld=0.0; if (ad_comm::no_ln_det_choleski_flag) { for (int ic=1;ic<=nsc;ic++) { if (allocated(block_diagonal_vhessian(ic))) { ld+=w(2*ic)*ln_det(block_diagonal_vhessian(ic),sgn); } } ld*=0.5; } else { for (int ic=1;ic<=nsc;ic++) { if (allocated(block_diagonal_vhessian(ic))) { ld+=w(2*ic)*ln_det_choleski(block_diagonal_vhessian(ic)); } } ld*=0.5; } vf+=ld; //vf+=us*0.91893853320467241; double f=value(vf); gradcalc(nvar,g); // put uhat back into the model gradient_structure::set_NO_DERIVATIVES(); vy(xs+1,xs+us).shift(1)=u0; initial_params::reset(vy); // get the values into the model gradient_structure::set_YES_DERIVATIVES(); pmin->lapprox->in_gauss_hermite_phase=0; ii=1; for (i=1;i<=xs;i++) xadjoint(i)=g(ii++); for (i=1;i<=us;i++) uadjoint(i)=g(ii++); for (ic=1;ic<=nsc;ic++) { int lus=lrea(ic); for (i=1;i<=lus;i++) { for (j=1;j<=lus;j++) { (*pmin->lapprox->block_diagonal_vhessianadjoint)(ic)(i,j)=g(ii++); } } } return f; }
void function_minimizer::pvm_master_mcmc_routine(int nmcmc,int iseed0,double dscale, int restart_flag) { uostream * pofs_psave=NULL; dmatrix mcmc_display_matrix; int mcmc_save_index=1; int mcmc_wrap_flag=0; int mcmc_gui_length=10000; int no_sd_mcmc=0; int on2=-1; if ( (on2=option_match(ad_comm::argc,ad_comm::argv,"-nosdmcmc"))>-1) no_sd_mcmc=1; if (stddev_params::num_stddev_params==0) { cerr << " You must declare at least one object of type sdreport " << endl << " to do the mcmc calculations" << endl; return; } { //ofstream of_bf("testbf"); //if (adjm_ptr) set_labels_for_mcmc(); ivector number_offsets; dvector lkvector; //double current_bf=0; double lcurrent_bf=0; double size_scale=1.0; double total_spread=200; //double total_spread=2500; uostream * pofs_sd = NULL; int nvar=initial_params::nvarcalc(); // get the number of active parameters int scov_option=0; dmatrix s_covar; dvector s_mean; int on=-1; int ncsim=25000; int nslots=800; //int nslots=3600; int initial_nsim=4800; int ntmp=0; int ncor=0; double bfsum=0; int ibfcount=0; double llbest; double lbmax; //if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcscov",ntmp))>-1) //{ scov_option=1; s_covar.allocate(1,nvar,1,nvar); s_mean.allocate(1,nvar); s_mean.initialize(); s_covar.initialize(); int ndvar=stddev_params::num_stddev_calc(); int numdvar=stddev_params::num_stddev_number_calc(); /* if (adjm_ptr) { mcmc_display_matrix.allocate(1,numdvar,1,mcmc_gui_length); number_offsets.allocate(1,numdvar); number_offsets=stddev_params::copy_all_number_offsets(); } */ dvector x(1,nvar); dvector scale(1,nvar); dmatrix values; int have_hist_flag=0; initial_params::xinit(x); dvector pen_vector(1,nvar); { initial_params::reset(dvar_vector(x),pen_vector); cout << pen_vector << endl << endl; } initial_params::mc_phase=0; initial_params::stddev_scale(scale,x); initial_params::mc_phase=1; dvector bmn(1,nvar); dvector mean_mcmc_values(1,ndvar); dvector s(1,ndvar); dvector h(1,ndvar); //dvector h; dvector square_mcmc_values(1,ndvar); square_mcmc_values.initialize(); mean_mcmc_values.initialize(); bmn.initialize(); int use_empirical_flag=0; int diag_option=0; int topt=0; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcdiag"))>-1) { diag_option=1; cout << " Setting covariance matrix to diagonal with entries " << dscale << endl; } dmatrix S(1,nvar,1,nvar); dvector sscale(1,nvar); if (!diag_option) { int on,nopt; int rescale_bounded_flag=0; double rescale_bounded_power=0.5; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcrb",nopt))>-1) { if (nopt) { int iii=atoi(ad_comm::argv[on+1]); if (iii < 1 || iii > 9) { cerr << " -mcrb argument must be integer between 1 and 9 --" " using default of 5" << endl; rescale_bounded_power=0.5; } else rescale_bounded_power=iii/10.0; } else { rescale_bounded_power=0.5; } rescale_bounded_flag=1; } if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcec"))>-1) { use_empirical_flag=1; } if (use_empirical_flag) { read_empirical_covariance_matrix(nvar,S,ad_comm::adprogram_name); } else if (!rescale_bounded_flag) { int tmp; read_covariance_matrix(S,nvar,tmp,sscale); } else { read_hessian_matrix_and_scale1(nvar,S,rescale_bounded_power, mcmc2_flag); //read_hessian_matrix_and_scale(nvar,S,pen_vector); } { // scale covariance matrix for model space dmatrix tmp(1,nvar,1,nvar); for (int i=1;i<=nvar;i++) { tmp(i,i)=S(i,i)*(scale(i)*scale(i)); for (int j=1;j<i;j++) { tmp(i,j)=S(i,j)*(scale(i)*scale(j)); tmp(j,i)=tmp(i,j); } } S=tmp; } } else { S.initialize(); for (int i=1;i<=nvar;i++) { S(i,i)=dscale; } } cout << sort(eigenvalues(S)) << endl; dmatrix chd = choleski_decomp( (dscale*2.4/sqrt(double(nvar))) * S); dmatrix chdinv=inv(chd); int sgn; dmatrix symbds(1,2,1,nvar); initial_params::set_all_simulation_bounds(symbds); ofstream ofs_sd1((char*)(ad_comm::adprogram_name + adstring(".mc2"))); { long int iseed=0; int number_sims; if (nmcmc<=0) { number_sims= 100000; } else { number_sims= nmcmc; } //cin >> iseed; if (iseed0<=0) { iseed=-36519; } else { iseed=-iseed0; } if (iseed>0) { iseed=-iseed; } cout << "Initial seed value " << iseed << endl; random_number_generator rng(iseed); rng.better_rand(); //better_rand(iseed); double lprob=0.0; double lpinv=0.0; double lprob3=0.0; // get lower and upper bounds independent_variables y(1,nvar); independent_variables parsave(1,nvar); // read in the mcmc values to date int ii=1; dmatrix hist; if (restart_flag) { int tmp=0; if (!no_sd_mcmc) { hist.allocate(1,ndvar,-nslots,nslots); tmp=read_hist_data(hist,h,mean_mcmc_values,s,parsave,iseed, size_scale); values.allocate(1,ndvar,-nslots,nslots); for (int i=1;i<=ndvar;i++) { values(i).fill_seqadd(mean_mcmc_values(i)-0.5*total_spread*s(i) +.5*h(i),h(i)); } } if (iseed>0) { iseed=-iseed; } double br=rng.better_rand(); if (tmp) have_hist_flag=1; chd=size_scale*chd; chdinv=chdinv/size_scale; } else { int on=-1; int nopt=0; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcpin",nopt))>-1) { if (nopt) { cifstream cif((char *)ad_comm::argv[on+1]); if (!cif) { cerr << "Error trying to open mcmc par input file " << ad_comm::argv[on+1] << endl; exit(1); } cif >> parsave; if (!cif) { cerr << "Error reading from mcmc par input file " << ad_comm::argv[on+1] << endl; exit(1); } } else { cerr << "Illegal option with -mcpin" << endl; } } else { ii=1; initial_params::copy_all_values(parsave,ii); } } ii=1; initial_params::restore_all_values(parsave,ii); gradient_structure::set_NO_DERIVATIVES(); ofstream ogs("sims"); ogs << nvar << " " << number_sims << endl; initial_params::xinit(y); send_int_to_slaves(1); double llc=-pvm_master_get_monte_carlo_value(nvar,y); send_int_to_slaves(1); llbest=-pvm_master_get_monte_carlo_value(nvar,y); lbmax=llbest; // store current mcmc variable values in param_values //dmatrix store_mcmc_values(1,number_sims,1,ndvar); #if defined(USE_BAYES_FACTORS) lkvector.allocate(1,number_sims); #endif dvector mcmc_values(1,ndvar); dvector mcmc_number_values; //if (adjm_ptr) mcmc_number_values.allocate(1,numdvar); int offs=1; stddev_params::copy_all_values(mcmc_values,offs); /* if (adjm_ptr) { offs=1; stddev_params::copy_all_number_values(mcmc_number_values,offs); } */ int change_ball=2500; int nopt; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcscale",nopt))>-1) { if (nopt) { int iii=atoi(ad_comm::argv[on+1]); if (iii <=0) { cerr << " Invalid option following command line option -mcball -- " << endl << " ignored" << endl; } else change_ball=iii; } } int iac=0; int liac=0; int isim=0; int itmp=0; double logr; int u_option=0; double ll; int s_option=1; int psvflag=0; if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcu"))>-1) { u_option=1; } if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcnoscale"))>-1) { s_option=0; } //cout << llc << " " << llc << endl; int iac_old=0; int i_old=0; { if (!restart_flag) { pofs_sd = new uostream((char*)(ad_comm::adprogram_name + adstring(".mcm"))); } int mcsave_flag=0; int mcrestart_flag=option_match(ad_comm::argc,ad_comm::argv,"-mcr"); if ( (on=option_match(ad_comm::argc,ad_comm::argv,"-mcsave"))>-1) { int jj=(int)atof(ad_comm::argv[on+1]); if (jj <=0) { cerr << " Invalid option following command line option -mcsave -- " << endl; } else { mcsave_flag=jj; if ( mcrestart_flag>-1) { // check that nvar is correct { uistream uis((char*)(ad_comm::adprogram_name + adstring(".psv"))); if (!uis) { cerr << "Error trying to open file" << ad_comm::adprogram_name + adstring(".psv") << " for mcrestart" << endl; cerr << " I am starting a new file " << endl; psvflag=1; } else { int nv1; uis >> nv1; if (nv1 !=nvar) { cerr << "wrong number of independent variables in" << ad_comm::adprogram_name + adstring(".psv") << cerr << " I am starting a new file " << endl; psvflag=1; } } } if (!psvflag) { pofs_psave= new uostream( (char*)(ad_comm::adprogram_name + adstring(".psv")),ios::app); } else { pofs_psave= new uostream((char*)(ad_comm::adprogram_name + adstring(".psv"))); } } else {