/*! */ 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 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; }