/*! */ void Postprocess_method::read(vector <string> words, unsigned int & pos, Program_options & options) { if(!readvalue(words, pos=0, configfile, "READCONFIG")) error("Need READCONFIG in Postprocess"); if(!readvalue(words, pos=0, nskip, "NSKIP")) nskip=0; evaluate_energy=true; if(haskeyword(words,pos=0,"NOENERGY")) evaluate_energy=false; vector <vector < string> > dens_words; vector<string> tmp_dens; pos=0; while(readsection(words, pos, tmp_dens, "DENSITY")) { dens_words.push_back(tmp_dens); } vector <vector < string> > avg_words; pos=0; while(readsection(words, pos, tmp_dens, "AVERAGE")) { avg_words.push_back(tmp_dens); } sys=NULL; allocate(options.systemtext[0], sys); sys->generatePseudo(options.pseudotext, pseudo); debug_write(cout, "wfdata allocate\n"); wfdata=NULL; if(options.twftext.size() < 1) error("Need TRIALFUNC section for POSTPROCESS"); allocate(options.twftext[0], sys, wfdata); average_var.Resize(avg_words.size()); average_var=NULL; for(int i=0; i< average_var.GetDim(0); i++) { allocate(avg_words[i], sys, wfdata, average_var(i)); } densplt.Resize(dens_words.size()); for(int i=0; i< densplt.GetDim(0); i++) { allocate(dens_words[i], sys, options.runid,densplt(i)); } }
//---------------------------------------------------------------------- int Postprocess_method::gen_point(Wavefunction * wf, Sample_point * sample, Config_save_point & configpos, doublevar weight, Properties_point & pt) { Primary guide; Properties_gather gather; configpos.restorePos(sample); pseudo->randomize(); if(evaluate_energy) gather.gatherData(pt,pseudo,sys,wfdata,wf,sample,&guide); pt.avgrets.Resize(1,average_var.GetDim(0)); //cout << mpi_info.node << " generating a point with " << average_var.GetDim(0) << " avgrets " << endl; for(int i=0; i< average_var.GetDim(0); i++) { average_var(i)->randomize(wfdata,wf,sys,sample); average_var(i)->evaluate(wfdata, wf, sys, pseudo,sample,pt, pt.avgrets(0,i)); } for(int i=0; i< densplt.GetDim(0); i++) densplt(i)->accumulate(sample,weight); pt.weight(0)=weight; }
int Dmc_method::generateVariables(Program_options & options) { if(!have_read_options) error("need to call Dmc_method::read before generateVariables"); if(have_allocated_variables) error("already allocated variables in Dmc_method"); have_allocated_variables=1; allocate(options.systemtext[0], mysys); mysys->generatePseudo(options.pseudotext, mypseudo); allocate(options.twftext[0], mysys, mywfdata); densplt.Resize(dens_words.size()); for(int i=0; i< densplt.GetDim(0); i++) { allocate(dens_words[i], mysys, options.runid,densplt(i)); } nldensplt.Resize(nldens_words.size()); for(int i=0; i< nldensplt.GetDim(0); i++) { allocate(nldens_words[i], mysys, options.runid,nldensplt(i)); } return 1; }
int Reptation_method::generateVariables(Program_options & options) { if(!have_read_options) { error("Need to call Reptation_method::read() before generateVariables()"); } have_generated_variables=1; debug_write(cout, "properties\n"); allocate(options.systemtext[0], sys ); allocate(options.twftext[0], sys, mywfdata); debug_write(cout, "Pseudopotential\n"); sys->generatePseudo(options.pseudotext, pseudo); densplt.Resize(dens_words.size()); for(int i=0; i< densplt.GetDim(0); i++) { allocate(dens_words[i], sys, options.runid,densplt(i)); } return 1; }
/*! */ void Dmc_method::runWithVariables(Properties_manager & prop, System * sys, Wavefunction_data * wfdata, Pseudopotential * pseudo, ostream & output) { allocateIntermediateVariables(sys, wfdata); if(!wfdata->supports(laplacian_update)) error("DMC doesn't support all-electron moves..please" " change your wave function to use the new Jastrow"); cout.precision(15); output.precision(10); prop.setSize(wf->nfunc(), nblock, nstep, nconfig, sys, wfdata); restorecheckpoint(readconfig, sys, wfdata, pseudo); prop.initializeLog(average_var); //MB: new properties manager for forward walking (one per each length) Array1 <Properties_manager> prop_fw; prop_fw.Resize(fw_length.GetSize()); if(max_fw_length){ for(int s=0;s<fw_length.GetSize();s++){ string logfile, label_temp; prop.getLog(logfile, label_temp); char strbuff[40]; sprintf(strbuff, "%d", fw_length(s)); label_temp+="_fw"; label_temp+=strbuff; prop_fw(s).setLog(logfile, label_temp); prop_fw(s).setSize(wf->nfunc(), nblock, nstep, nconfig, sys, wfdata); prop_fw(s).initializeLog(average_var); } } if(nhist==-1) nhist=1; doublevar teff=timestep; for(int block=0; block < nblock; block++) { int totkilled=0; int totbranch=0; int totpoints=0; for(int step=0; step < nstep; ) { int npsteps=min(feedback_interval, nstep-step); Dynamics_info dinfo; doublevar acsum=0; doublevar deltar2=0; Array1 <doublevar> epos(3); doublevar avg_acceptance=0; for(int walker=0; walker < nconfig; walker++) { pts(walker).config_pos.restorePos(sample); wf->updateLap(wfdata, sample); //------Do several steps without branching for(int p=0; p < npsteps; p++) { pseudo->randomize(); for(int e=0; e< nelectrons; e++) { int acc; acc=dyngen->sample(e, sample, wf, wfdata, guidingwf, dinfo, timestep); if(dinfo.accepted) deltar2+=dinfo.diffusion_rate/(nconfig*nelectrons*npsteps); if(dinfo.accepted) { pts(walker).age(e)=0; } else { pts(walker).age(e)++; } avg_acceptance+=dinfo.acceptance/(nconfig*nelectrons*npsteps); if(acc>0) acsum++; } totpoints++; Properties_point pt; if(tmoves or tmoves_sizeconsistent) { //------------------T-moves doTmove(pt,pseudo,sys,wfdata,wf,sample,guidingwf); } ///---------------------------------done with the T-moves else { mygather.gatherData(pt, pseudo, sys, wfdata, wf, sample, guidingwf); } Dmc_history new_hist; new_hist.main_en=pts(walker).prop.energy(0); pts(walker).past_energies.push_front(new_hist); deque<Dmc_history> & past(pts(walker).past_energies); if(past.size() > nhist) past.erase(past.begin()+nhist, past.end()); pts(walker).prop=pt; if(!pure_dmc) { pts(walker).weight*=getWeight(pts(walker),teff,etrial); //Introduce potentially a small bias to avoid instability. if(pts(walker).weight>max_poss_weight) pts(walker).weight=max_poss_weight; } else pts(walker).weight=getWeightPURE_DMC(pts(walker),teff,etrial); if(pts(walker).ignore_walker) { pts(walker).ignore_walker=0; pts(walker).weight=1; pts(walker).prop.count=0; } pts(walker).prop.weight=pts(walker).weight; //This is somewhat inaccurate..will need to change it later //For the moment, the autocorrelation will be slightly //underestimated pts(walker).prop.parent=walker; pts(walker).prop.nchildren=1; pts(walker).prop.children(0)=walker; pts(walker).prop.avgrets.Resize(1,average_var.GetDim(0)); for(int i=0; i< average_var.GetDim(0); i++) { average_var(i)->randomize(wfdata,wf,sys,sample); average_var(i)->evaluate(wfdata, wf, sys, pseudo, sample, pts(walker).prop.avgrets(0,i)); } prop.insertPoint(step+p, walker, pts(walker).prop); for(int i=0; i< densplt.GetDim(0); i++) densplt(i)->accumulate(sample,pts(walker).prop.weight(0)); for(int i=0; i< nldensplt.GetDim(0); i++) nldensplt(i)->accumulate(sample,pts(walker).prop.weight(0), wfdata,wf); //MB: making the history of prop.avgrets for forward walking if(max_fw_length){ forwardWalking(walker, step+p,prop_fw); }//if FW } pts(walker).config_pos.savePos(sample); } //---Finished moving all walkers doublevar accept_ratio=acsum/(nconfig*nelectrons*npsteps); teff=timestep*accept_ratio; //deltar2/rf_diffusion; updateEtrial(feedback); step+=npsteps; int nkilled; if(!pure_dmc) nkilled=calcBranch(); else nkilled=0; totkilled+=nkilled; totbranch+=nkilled; } ///----Finished block if(!low_io || block==nblock-1) { savecheckpoint(storeconfig,sample); for(int i=0; i< densplt.GetDim(0); i++) densplt(i)->write(); for(int i=0; i< nldensplt.GetDim(0); i++) nldensplt(i)->write(log_label); } if(!pure_dmc){ prop.endBlock(); } else prop.endBlockSHDMC(); if(max_fw_length){ for(int s=0;s<fw_length.GetSize();s++){ //prop_fw(s).endBlock(); prop_fw(s).endBlock_per_step(); } } totbranch=parallel_sum(totbranch); totkilled=parallel_sum(totkilled); totpoints=parallel_sum(totpoints); Properties_final_average finavg; prop.getFinal(finavg); eref=finavg.avg(Properties_types::total_energy,0); updateEtrial(feedback); doublevar maxage=0; doublevar avgage=0; for(int w=0;w < nconfig; w++) { for(int e=0; e< nelectrons; e++) { if(maxage<pts(w).age(e)) maxage=pts(w).age(e); avgage+=pts(w).age(e); } } avgage/=(nconfig*nelectrons); if(output) { //cout << "Block " << block // << " nconfig " << totconfig // << " etrial " << etrial << endl; if(global_options::rappture ) { ofstream rapout("rappture.log"); rapout << "=RAPPTURE-PROGRESS=>" << int(100.0*doublevar(block+1)/doublevar(nblock)) << " Diffusion Monte Carlo" << endl; cout << "=RAPPTURE-PROGRESS=>" << int(100.0*doublevar(block+1)/doublevar(nblock)) << " Diffusion Monte Carlo" << endl; rapout.close(); } output << "***" << endl; output << "Block " << block << " etrial " << etrial << endl; output << "maximum age " << maxage << " average age " << avgage << endl; dyngen->showStats(output); prop.printBlockSummary(output); output << "Branched " << totbranch << " times. So a branch every " << doublevar(totpoints)/doublevar(totbranch) << " steps " << endl; } dyngen->resetStats(); } if(output) { output << "\n ----------Finished DMC------------\n\n"; prop.printSummary(output,average_var); //MB: final printout for FW if(max_fw_length){ for(int s=0;s<fw_length.GetSize();s++) prop_fw(s).printSummary(output,average_var); } } wfdata->clearObserver(); deallocateIntermediateVariables(); }
/*! */ void Reptation_method::runWithVariables(Properties_manager & prop, System * sys, Wavefunction_data * wfdata, Pseudopotential * psp, ostream & output) { allocateIntermediateVariables(sys, wfdata); prop.setSize(wf->nfunc(), nblock, nstep, 1, sys, wfdata); prop.initializeLog(average_var); Properties_manager prop_center; string logfile, label_temp; prop.getLog(logfile, label_temp); label_temp+="_cen"; prop_center.setLog(logfile, label_temp); prop_center.setSize(wf->nfunc(), nblock, nstep, 1, sys, wfdata); prop_center.initializeLog(average_var); cout.precision(10); output.precision(10); Sample_point * center_samp(NULL); sys->generateSample(center_samp); Reptile_point pt; Array1 <Reptile> reptiles; int nreptile=1; if(!readcheck(readconfig,reptiles)) { Array1 <Config_save_point> configs; generate_sample(sample,wf,wfdata,guidewf,nreptile,configs); reptiles.Resize(nreptile); for(int r=0; r< nreptile; r++) { reptiles[r].direction=1; configs(r).restorePos(sample); wf->notify(all_electrons_move,0); wf->updateLap(wfdata,sample); for(int i=0; i< reptile_length; i++) { doublevar main_diffusion; slither(1,reptiles[r].reptile, mygather,pt,main_diffusion); reptiles[r].reptile.push_back(pt); } } } nreptile=reptiles.GetDim(0); //assert(reptile.size()==reptile_length); //Branch limiting variables //we start off with no limiting, and establish the parameters after the //first block. This seems to be reasonably stable, since it's mostly //to keep the reptile from getting stuck. eref=0; energy_cutoff=1e16; //--------begin averaging.. Array3 <doublevar> derivatives_block(nblock, sys->nIons(), 3); for(int block=0; block< nblock; block++) { //clock_t block_start_time=clock(); doublevar avg_age=0; doublevar max_age=0; doublevar main_diff=0; double ntry=0, naccept=0; double nbounce=0; for(int r=0; r< nreptile; r++) { Reptile & curr_reptile=reptiles[r]; //Control variable that will be set to one when //we change direction, which signals to recalculate //the wave function int recalc=1; for(int step=0; step< nstep; step++) { psp->randomize(); if(recalc) { if(curr_reptile.direction==1) curr_reptile.reptile[reptile_length-1].restorePos(sample); else curr_reptile.reptile[0].restorePos(sample); } doublevar main_diffusion; doublevar accept=slither(curr_reptile.direction, curr_reptile.reptile,mygather, pt, main_diffusion); ntry++; if(accept+rng.ulec() > 1.0) { recalc=0; naccept++; main_diff+=main_diffusion; if(curr_reptile.direction==1) { curr_reptile.reptile.pop_front(); curr_reptile.reptile.push_back(pt); } else { curr_reptile.reptile.pop_back(); curr_reptile.reptile[0].branching=pt.branching; curr_reptile.reptile.push_front(pt); } } else { recalc=1; curr_reptile.direction*=-1; nbounce++; } for(deque<Reptile_point>::iterator i=curr_reptile.reptile.begin(); i!=curr_reptile.reptile.end(); i++) { i->age++; avg_age+=i->age/reptile_length; if(i->age > max_age) max_age=i->age; } Properties_point avgpt; get_avg(curr_reptile.reptile, avgpt); avgpt.parent=0; avgpt.nchildren=1; //just one walker avgpt.children(0)=0; prop.insertPoint(step, 0, avgpt); int cpt=reptile_length/2+1; Properties_point centpt; get_center_avg(curr_reptile.reptile, centpt); centpt.parent=0; centpt.nchildren=1; centpt.children(0)=0; prop_center.insertPoint(step, 0, centpt); curr_reptile.reptile[cpt].restorePos(center_samp); for(int i=0; i< densplt.GetDim(0); i++) densplt(i)->accumulate(center_samp,1.0); if(center_trace != "" && (block*nstep+step)%trace_wait==0) { ofstream checkfile(center_trace.c_str(), ios::app); if(!checkfile)error("Couldn't open ", center_trace); checkfile << "SAMPLE_POINT { \n"; write_config(checkfile, sample); checkfile << "}\n\n"; } } //step } //reptile prop.endBlock(); prop_center.endBlock(); double ntot=parallel_sum(nstep); Properties_block lastblock; prop.getLastBlock(lastblock); eref=lastblock.avg(Properties_types::total_energy,0); energy_cutoff=10*sqrt(lastblock.var(Properties_types::total_energy,0)); nbounce=parallel_sum(nbounce); naccept=parallel_sum(naccept); ntry=parallel_sum(ntry); avg_age=parallel_sum(avg_age); for(int i=0; i< densplt.GetDim(0); i++) densplt(i)->write(); storecheck(reptiles, storeconfig); main_diff=parallel_sum(main_diff); if(output) { output << "****Block " << block << " acceptance " << naccept/ntry << " average steps before bounce " << ntot/nbounce << endl; output << "average age " << avg_age/ntot << " max age " << max_age << endl; output << "eref " << eref << " cutoff " << energy_cutoff << endl; output << "Green's function sampler:" << endl; sampler->showStats(output); prop.printBlockSummary(output); output << "Center averaging: " << endl; prop_center.printBlockSummary(output); } sampler->resetStats(); //clock_t block_end_time=clock(); //cout << mpi_info.node << ":CPU block time " //// << double(block_end_time-block_start_time)/double(CLOCKS_PER_SEC) // << endl; } //block if(output) { output << "############## Reptation Done ################\n"; output << "End averages " << endl; prop.printSummary(output,average_var); output << "Center averages " << endl; prop_center.printSummary(output,average_var); //Print out a PDB file with one of the reptiles, for visualization purposes if(print_pdb) { ofstream pdbout("rmc.pdb"); pdbout.precision(3); pdbout << "REMARK 4 Mode COMPLIES WITH FORMAT V. 2.0\n"; int nelectrons=sample->electronSize(); int counter=1; string name="H"; for(int e=0; e<nelectrons; e++) { for(deque<Reptile_point>::iterator i=reptiles[0].reptile.begin(); i!=reptiles[0].reptile.end(); i++) { pdbout<<"ATOM"<<setw(7)<< counter <<" " <<name<<" UNK 1" <<setw(12)<< i->electronpos[e][0] <<setw(8)<< i->electronpos[e][1] <<setw(8)<< i->electronpos[e][2] << " 1.00 0.00\n"; counter++; } } int nions=sys->nIons(); Array1 <doublevar> ionpos(3); vector <string> atomnames; sys->getAtomicLabels(atomnames); for(int i=0; i< nions; i++) { sys->getIonPos(i,ionpos); pdbout<<"ATOM"<<setw(7)<< counter <<" " <<atomnames[i]<<" UNK 1" <<setw(12)<< ionpos[0] <<setw(8)<< ionpos[1] <<setw(8)<< ionpos[2] << " 1.00 0.00\n"; } counter=1; for(int e=0; e<nelectrons; e++) { for(deque<Reptile_point>::iterator i=reptiles[0].reptile.begin(); i!=reptiles[0].reptile.end(); i++) { if(i != reptiles[0].reptile.begin()) { pdbout << "CONECT" << setw(5) << counter << setw(5) << counter-1 << endl; } counter++; } } } //------------Done PDB file } delete center_samp; wfdata->clearObserver(); deallocateIntermediateVariables(); }
void Postprocess_method::run(Program_options & options, ostream & output) { Sample_point * sample=NULL; Wavefunction * wf=NULL; sys->generateSample(sample); wfdata->generateWavefunction(wf); sample->attachObserver(wf); Properties_gather gather; Primary guide; int nelec=sample->electronSize(); int ndim=3; int npoints_tot=0; FILE * f; if(mpi_info.node==0) { f=fopen(configfile.c_str(),"r"); if(ferror(f)) error("Could not open",configfile); fseek(f,0,SEEK_END); long int lSize=ftell(f); rewind(f); npoints_tot=lSize/(sizeof(doublevar)*(nelec*3+1+4)); output << "Estimated number of samples in this file: " << npoints_tot << endl; output << "We are skipping the first " << nskip << " of these " << endl; Config_save_point tmpconfig; for(int i=0; i< nskip; i++) { doublevar weight; tmpconfig.readBinary(f,nelec,ndim,weight); // doublevar weight; // if(!fread(&weight,sizeof(doublevar),1,f)) error("Misformatting in binary file",configfile, " perhaps nskip is too large?"); } } #ifdef USE_MPI if(mpi_info.nprocs<2) error("POSTPROCESS must be run with at least 2 processes if it is run in parallel."); if(mpi_info.node==0) { master(wf,sample,f,output); } else { worker(wf,sample); } #else Config_save_point tmpconfig; Properties_point pt; pt.setSize(1); int npoints=0; Postprocess_average postavg(average_var.GetDim(0)); doublevar weight; while(tmpconfig.readBinary(f,nelec,ndim,weight)) { tmpconfig.restorePos(sample); gen_point(wf,sample,tmpconfig,weight,pt); postavg.update_average(pt); npoints++; doublevar progress=doublevar(npoints)/doublevar(npoints_tot); if(fabs(progress*10-int(progress*10)) < 0.5/npoints_tot) { cout << "progress: " << progress*100 << "% done" << endl; } } postavg.print(average_var,output); #endif //USE_MPI for(int i=0; i< densplt.GetDim(0); i++) densplt(i)->write(); if(mpi_info.node==0) fclose(f); delete sample; delete wf; }