///initialize polymers void ReptationMC::initReptile() { //overwrite the number of cuts for Bounce algorithm if(UseBounce) NumCuts = 1; app_log() << "Moving " << NumCuts << " for each reptation step" << endl; //Reptile is NOT allocated. Create one. if(Reptile == 0) { MCWalkerConfiguration::iterator it(W.begin()); Walker_t* cur(*it); cur->Weight=1.0; W.R = cur->R; //DistanceTable::update(W); W.update(); RealType logpsi(Psi.evaluateLog(W)); RealType eloc_cur = H.evaluate(W); cur->resetProperty(logpsi,Psi.getPhase(),eloc_cur); H.saveProperty(cur->getPropertyBase()); cur->Drift = W.G; Reptile = new PolymerChain(cur,PolymerLength,NumCuts); Reptile->resizeArrays(1); } //If ClonePolyer==false, generate configuration using diffusion //Not so useful if(!ClonePolymer) { Walker_t* cur((*Reptile)[0]); RealType g = std::sqrt(Tau); for(int i=0; i<NumCuts-1; i++ ) { //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); W.R = cur->R + g*deltaR + Tau*cur->Drift; //update the distance table associated with W //DistanceTable::update(W); W.update(); //evaluate wave function RealType logpsic(Psi.evaluateLog(W)); RealType e0=H.evaluate(W); cur = (*Reptile)[i+1]; cur->resetProperty(logpsic,Psi.getPhase(),e0); H.saveProperty(cur->getPropertyBase()); cur->R = W.R; cur->Drift = W.G; } } }
bool DummyQMC::run() { m_oneover2tau = 0.5/Tau; m_sqrttau = sqrt(Tau); cout << "Lattice of ParticleSet " << endl; W.Lattice.print(cout); //property container to hold temporary properties, such as a local energy //MCWalkerConfiguration::PropertyContainer_t Properties; MCWalkerConfiguration::Walker_t& thisWalker(**W.begin()); //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); //new poosition W.R = m_sqrttau*deltaR + thisWalker.R + thisWalker.Drift; //apply boundary condition: put everything back to a unit cell W.applyBC(W.R); //update the distance table associated with W W.update(); //evaluate wave function //update the properties: note that we are getting \f$\sum_i \ln(|psi_i|)\f$ and catching the sign separately ValueType logpsi(Psi.evaluateLog(W)); RealType eloc=H.evaluate(W); return true; }
void VMCParticleByParticle::runBlockWithDrift() { IndexType step = 0; MCWalkerConfiguration::iterator it_end(W.end()); do { //advanceWalkerByWalker(); MCWalkerConfiguration::iterator it(W.begin()); while(it != it_end) { //MCWalkerConfiguration::WalkerData_t& w_buffer = *(W.DataSet[iwalker]); Walker_t& thisWalker(**it); Buffer_t& w_buffer(thisWalker.DataSet); W.R = thisWalker.R; w_buffer.rewind(); W.copyFromBuffer(w_buffer); Psi.copyFromBuffer(W,w_buffer); RealType psi_old = thisWalker.Properties(SIGN); RealType psi = psi_old; //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); bool moved = false; for(int iat=0; iat<W.getTotalNum(); iat++) { PosType dr = m_sqrttau*deltaR[iat]+thisWalker.Drift[iat]; PosType newpos = W.makeMove(iat,dr); //RealType ratio = Psi.ratio(W,iat); RealType ratio = Psi.ratio(W,iat,dG,dL); G = W.G+dG; //RealType forwardGF = exp(-0.5*dot(deltaR[iat],deltaR[iat])); //dr = (*it)->R[iat]-newpos-Tau*G[iat]; //RealType backwardGF = exp(-oneover2tau*dot(dr,dr)); RealType logGf = -0.5e0*dot(deltaR[iat],deltaR[iat]); RealType scale=getDriftScale(Tau,G); //COMPLEX WARNING //dr = thisWalker.R[iat]-newpos-scale*G[iat]; dr = thisWalker.R[iat]-newpos-scale*real(G[iat]); RealType logGb = -m_oneover2tau*dot(dr,dr); RealType prob = std::min(1.0e0,ratio*ratio*exp(logGb-logGf)); //alternatively if(Random() < prob) { moved = true; ++nAccept; W.acceptMove(iat); Psi.acceptMove(W,iat); W.G = G; W.L += dL; //thisWalker.Drift = scale*G; assignDrift(scale,G,thisWalker.Drift); } else { ++nReject; W.rejectMove(iat); Psi.rejectMove(iat); } } if(moved) { w_buffer.rewind(); W.copyToBuffer(w_buffer); psi = Psi.evaluate(W,w_buffer); thisWalker.R = W.R; RealType eloc=H.evaluate(W); thisWalker.resetProperty(log(abs(psi)), psi,eloc); H.saveProperty(thisWalker.getPropertyBase()); } else { ++nAllRejected; } ++it; } ++step;++CurrentStep; Estimators->accumulate(W); //if(CurrentStep%100 == 0) updateWalkers(); } while(step<nSteps); }
void VMCParticleByParticle::runBlockWithoutDrift() { IndexType step = 0; MCWalkerConfiguration::iterator it_end(W.end()); do { //advanceWalkerByWalker(); MCWalkerConfiguration::iterator it(W.begin()); while(it != it_end) { //MCWalkerConfiguration::WalkerData_t& w_buffer = *(W.DataSet[iwalker]); Walker_t& thisWalker(**it); Buffer_t& w_buffer(thisWalker.DataSet); W.R = thisWalker.R; w_buffer.rewind(); W.copyFromBuffer(w_buffer); Psi.copyFromBuffer(W,w_buffer); RealType psi_old = thisWalker.Properties(SIGN); RealType psi = psi_old; for(int iter=0; iter<nSubSteps; iter++) { makeGaussRandom(deltaR); bool stucked=true; for(int iat=0; iat<W.getTotalNum(); iat++) { PosType dr = m_sqrttau*deltaR[iat]; PosType newpos = W.makeMove(iat,dr); RealType ratio = Psi.ratio(W,iat); //RealType ratio = Psi.ratio(W,iat,dG,dL); RealType prob = std::min(1.0e0,ratio*ratio); //alternatively if(Random() < prob) { stucked=false; ++nAccept; W.acceptMove(iat); Psi.acceptMove(W,iat); //W.G+=dG; //W.L+=dL; } else { ++nReject; W.rejectMove(iat); Psi.rejectMove(iat); } } if(stucked) { ++nAllRejected; } } thisWalker.R = W.R; w_buffer.rewind(); W.updateBuffer(w_buffer); RealType logpsi = Psi.updateBuffer(W,w_buffer); //W.copyToBuffer(w_buffer); //RealType logpsi = Psi.evaluate(W,w_buffer); RealType eloc=H.evaluate(W); thisWalker.resetProperty(logpsi,Psi.getPhase(),eloc); H.saveProperty(thisWalker.getPropertyBase()); ++it; } ++step;++CurrentStep; Estimators->accumulate(W); } while(step<nSteps); }
void ReptationMC::moveReptile(){ //RealType oneovertau = 1.0/Tau; //RealType oneover2tau = 0.5*oneovertau; RealType tauover2 = 0.5*Tau; RealType g = std::sqrt(Tau); typedef MCWalkerConfiguration::PropertyContainer_t PropertyContainer_t; if(!UseBounce && Random()<0.5) { Reptile->flip(); NumTurns++; } Walker_t* anchor = Reptile->makeEnds(); //save the local energies of the anchor and tails //eloc_xp = the energy of the front //eloc_yp = the energy of the proposed move //eloc_x = the energy of the tail //eloc_y = the energy of the tail-1 RealType eloc_xp = anchor->Properties(LOCALENERGY); RealType eloc_x = Reptile->tails[0]->Properties(LOCALENERGY); RealType eloc_y = Reptile->tails[1]->Properties(LOCALENERGY); NumCuts = Reptile->NumCuts; RealType Wpolymer=0.0; for(int i=0; i<NumCuts; ) { Walker_t* head=Reptile->heads[i]; //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); W.R = anchor->R + g*deltaR + Tau* anchor->Drift; //update the distance table associated with W //DistanceTable::update(W); W.update(); //evaluate wave function RealType logpsi(Psi.evaluateLog(W)); //update the properties of the front chain //RealType eloc_yp = head->Properties(LOCALENERGY) = H.evaluate(W); //H.copy(head->getEnergyBase()); //head->Properties(LOCALPOTENTIAL) = H.getLocalPotential(); RealType eloc_yp = H.evaluate(W); head->resetProperty(logpsi,Psi.getPhase(),eloc_yp); H.saveProperty(head->getPropertyBase()); head->R = W.R; //ValueType vsq = Dot(W.G,W.G); //ValueType scale = ((-1.0+sqrt(1.0+2.0*Tau*vsq))/vsq); //head->Drift = scale*W.G; head->Drift = W.G; //\f${x-y-\tau\nabla \ln \Psi_{T}(y))\f$ //deltaR = anchor->R - W.R - heads[i]->Drift; //Gdrift *= exp(-oneover2tau*Dot(deltaR,deltaR)); /* \f$ X= \{R_0, R_1, ... , R_M\}\f$ \f$ X' = \{R_1, .., R_M, R_{M+1}\}\f$ \f[ G_B(R_{M+1}\leftarrow R_{M}, \tau)/G_B(R_{0}\leftarrow R_{1}, \tau) = exp\(-\tau/2[E_L(R_{M+1})+E_L(R_M)-E_L(R_1)-E_L(R_0)]\)\f] * - eloc_yp = \f$E_L(R_{M+1})\f$ - eloc_xp = \f$E_L(R_{M})\f$ - eloc_y = \f$E_L(R_{1})\f$ - eloc_x = \f$E_L(R_{0})\f$ */ //Wpolymer *= exp(-oneover2tau*(eloc_yp+eloc_xp-eloc_x-eloc_y)); Wpolymer +=(eloc_yp+eloc_xp-eloc_x-eloc_y); //move the anchor and swap the local energies for Wpolymer anchor=head; //increment the index i++; if(i<NumCuts) { eloc_xp = eloc_yp; eloc_x = eloc_y; eloc_y = Reptile->tails[i+1]->Properties(LOCALENERGY); } } Wpolymer = std::exp(-tauover2*Wpolymer); double accept = std::min(1.0,Wpolymer); if(Random() < accept){//move accepted Reptile->updateEnds(); ++nAccept; } else { ++nReject; if(UseBounce) { NumTurns++; Reptile->flip(); } } //RealType Bounce = UseBounce ? 1.0-accept: 0.5; //if(Random()<Bounce) { // Reptile->flip(); // LogOut->getStream() << "Bounce = " << Bounce << " " << NumTurns << " " << polymer.MoveHead << endl; // NumTurns++;//increase the number of turns //} }
void DMCWOS::advanceWalkerByWalker(BRANCHER& Branch) { //Pooma::Clock timer; RealType oneovertau = 1.0/Tau; RealType oneover2tau = 0.5*oneovertau; RealType g = sqrt(Tau); RealType vwos; //MCWalkerConfiguration::PropertyContainer_t Properties; int nh = H.size()+1; // extract the WOS potential //WOSPotential* wos = dynamic_cast<WOSPotential*>(H.getHamiltonian("wos")); MCWalkerConfiguration::iterator it = W.begin(); MCWalkerConfiguration::iterator it_end = W.end(); while(it != it_end) { (*it)->Properties(WEIGHT) = 1.0; (*it)->Properties(MULTIPLICITY) = 1.0; //copy the properties of the working walker W.Properties = (*it)->Properties; //save old local energy ValueType eold = W.Properties(LOCALENERGY); /* If doing WOS then we have to redo the old calculation for ROld to calculate G(Rold,Rnew). So we extract wos from the Hamiltonian and re-evaluate it. */ ValueType emixed = eold; if(wos_ref) { W.R = (*it)->R; DistanceTable::update(W); ValueType psi(Psi.evaluateLog(W));//not used anyway eold = H.evaluate(W); emixed += 0.5*W.Properties(WOSVAR)*Tau_var; } //ValueType emixed = eold + 0.5*W.Properties(WOSVAR)*Tau; //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); W.R = g*deltaR + (*it)->R + (*it)->Drift; //update the distance table associated with W DistanceTable::update(W); //evaluate wave function ValueType logpsi(Psi.evaluateLog(W)); //update the properties W.Properties(LOCALENERGY) = H.evaluate(W); W.Properties(LOGPSI) =logpsi; W.Properties(SIGN) = Psi.getSign(); bool accepted=false; //deltaR = W.R - (*it)->R - (*it)->Drift; //RealType forwardGF = exp(-oneover2tau*Dot(deltaR,deltaR)); //RealType forwardGF = exp(-0.5*Dot(deltaR,deltaR)); RealType logGf = -0.5*Dot(deltaR,deltaR); //scale the drift term to prevent persistent cofigurations ValueType vsq = Dot(W.G,W.G); //converting gradients to drifts, D = tau*G (reuse G) // W.G *= Tau;//original implementation with bare drift ValueType scale = ((-1.0+sqrt(1.0+2.0*Tau*vsq))/vsq); drift = scale*W.G; deltaR = (*it)->R - W.R - drift; //RealType backwardGF = exp(-oneover2tau*Dot(deltaR,deltaR)); RealType logGb = -oneover2tau*Dot(deltaR,deltaR); //set acceptance probability //RealType prob= std::min(backwardGF/forwardGF*W.Properties(PSISQ)/(*it)->Properties(PSISQ),1.0); //RealType prob= std::min(exp(logGb-logGf)*W.Properties(PSISQ)/(*it)->Properties(PSISQ),1.0); RealType prob= std::min(exp(logGb-logGf +2.0*(W.Properties(LOGPSI)-(*it)->Properties(LOGPSI))),1.0); if(Random() > prob){ (*it)->Properties(AGE)++; emixed += emixed; } else { accepted=true; W.Properties(AGE) = 0; (*it)->R = W.R; (*it)->Drift = drift; (*it)->Properties = W.Properties; H.copy((*it)->getEnergyBase()); emixed += W.Properties(LOCALENERGY) + 0.5*W.Properties(WOSVAR)*Tau_var; } //calculate the weight and multiplicity ValueType M = Branch.branchGF(Tau,emixed*0.5,1.0-prob); if((*it)->Properties(AGE) > 3.0) M = min(0.5,M); if((*it)->Properties(AGE) > 0.9) M = min(1.0,M); (*it)->Properties(WEIGHT) = M; (*it)->Properties(MULTIPLICITY) = M + Random(); //node-crossing: kill it for the time being if(Branch(W.Properties(SIGN),(*it)->Properties(SIGN))) { accepted=false; (*it)->Properties(WEIGHT) = 0.0; (*it)->Properties(MULTIPLICITY) = 0.0; } /* if(Branch(W.Properties(PSI),(*it)->Properties(PSI))) { accepted=false; (*it)->Properties(WEIGHT) = 0.0; (*it)->Properties(MULTIPLICITY) = 0.0; } else { RealType logGf = -0.5*Dot(deltaR,deltaR); ValueType vsq = Dot(W.G,W.G); ValueType scale = ((-1.0+sqrt(1.0+2.0*Tau*vsq))/vsq); drift = scale*W.G; deltaR = (*it)->R - W.R - drift; RealType logGb = -oneover2tau*Dot(deltaR,deltaR); RealType prob = std::min(exp(logGb-logGf)*W.Properties(PSISQ)/(*it)->Properties(PSISQ),1.0); if(Random() > prob){ (*it)->Properties(AGE)++; emixed += emixed; } else { accepted=true; W.Properties(AGE) = 0; (*it)->R = W.R; (*it)->Drift = drift; (*it)->Properties = W.Properties; // H.update(W.Energy[(*it)->ID]); //H.get((*it)->E); H.copy((*it)->getEnergyBase()); emixed += W.Properties(LOCALENERGY) + 0.5*W.Properties(WOSVAR)*Tau_var; } //calculate the weight and multiplicity ValueType M = Branch.branchGF(Tau,emixed*0.5,1.0-prob); if((*it)->Properties(AGE) > 3.0) M = min(0.5,M); if((*it)->Properties(AGE) > 0.9) M = min(1.0,M); (*it)->Properties(WEIGHT) = M; (*it)->Properties(MULTIPLICITY) = M + Random(); } */ if(accepted) ++nAccept; else ++nReject; ++it; } }
void VMC::advanceAllWalkers() { static ParticleSet::ParticlePos_t deltaR(W.getTotalNum()); WalkerSetRef Wref(W); Wref.resize(W.getActiveWalkers(),W.getTotalNum()); //Pooma::Clock timer; RealType oneovertau = 1.0/Tau; RealType oneover2tau = 0.5*oneovertau; RealType g = sqrt(Tau); MCWalkerConfiguration::PropertyContainer_t Properties; makeGaussRandom(Wref.R); Wref.R *= g; int nptcl = W.getTotalNum(); int iw = 0; MCWalkerConfiguration::iterator it = W.begin(); while(it != W.end()) { const ParticleSet::ParticlePos_t& r = (*it)->R; for(int jat=0; jat<nptcl; jat++) { Wref.R(iw,jat) += r(jat) + (*it)->Drift(jat); } iw++; it++; } DistanceTable::update(Wref); OrbitalBase::ValueVectorType psi(iw), energy(iw); Psi.evaluate(Wref,psi); H.evaluate(Wref,energy); //multiply tau to convert gradient to drift term Wref.G *= Tau; iw = 0; it = W.begin(); while(it != W.end()) { ValueType eold = Properties(LocalEnergy); for(int iat=0; iat<nptcl; iat++) deltaR(iat) = Wref.R(iw,iat) - (*it)->R(iat) - (*it)->Drift(iat); RealType forwardGF = exp(-oneover2tau*Dot(deltaR,deltaR)); for(int iat=0; iat<nptcl; iat++) deltaR(iat) = (*it)->R(iat) - Wref.R(iw,iat) - Wref.G(iw,iat); RealType backwardGF = exp(-oneover2tau*Dot(deltaR,deltaR)); ValueType psisq = psi(iw)*psi(iw); if(Random() > backwardGF/forwardGF*psisq/(*it)->Properties(PsiSq)) { ++nReject; (*it)->Properties(Age) += 1; } else { (*it)->Properties(Age) = 0; for(int iat=0; iat<nptcl; iat++) (*it)->R(iat) = Wref.R(iw,iat); for(int iat=0; iat<nptcl; iat++) (*it)->Drift(iat) = Wref.G(iw,iat); (*it)->Properties(Sign) = psi(iw); (*it)->Properties(PsiSq) = psisq; (*it)->Properties(LocalEnergy) = energy(iw); ++nAccept; } iw++;it++; } }
void VMC::advanceWalkerByWalker() { //Pooma::Clock timer; RealType oneovertau = 1.0/Tau; RealType oneover2tau = 0.5*oneovertau; RealType g = sqrt(Tau); MCWalkerConfiguration::PropertyContainer_t Properties; static ParticleSet::ParticlePos_t deltaR(W.getTotalNum()); static ParticleSet::ParticlePos_t drift(W.getTotalNum()); int nh = H.size()+1; for (MCWalkerConfiguration::iterator it = W.begin(); it != W.end(); ++it) { //copy the properties of the working walker Properties = (*it)->Properties; //save old local energy ValueType eold = Properties(LocalEnergy); //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); W.R = g*deltaR + (*it)->R + (*it)->Drift; //update the distance table associated with W DistanceTable::update(W); //evaluate wave function ValueType psi = Psi.evaluate(W); //update the properties Properties(PsiSq) = psi*psi; Properties(Sign) = psi; Properties(LocalEnergy) = H.evaluate(W); // deltaR = W.R - (*it)->R - (*it)->Drift; // RealType forwardGF = exp(-oneover2tau*Dot(deltaR,deltaR)); RealType forwardGF = exp(-0.5*Dot(deltaR,deltaR)); //converting gradients to drifts, D = tau*G (reuse G) //W.G *= Tau;//original implementation with bare drift ValueType vsq = Dot(W.G,W.G); ValueType scale = ((-1.0+sqrt(1.0+2.0*Tau*vsq))/vsq); drift = scale*W.G; // W.G *= scale; deltaR = (*it)->R - W.R - drift; RealType backwardGF = exp(-oneover2tau*Dot(deltaR,deltaR)); //forwardGF/backwardGF*Properties(PsiSq)/(*it)->Properties(PsiSq) //RealType prob = min(Properties(PsiSq)/(*it)->Properties(PsiSq),1.0); //cout << "forward/backward " << forwardGF << " " << backwardGF << endl; if(Random() > backwardGF/forwardGF*Properties(PsiSq)/(*it)->Properties(PsiSq)) { (*it)->Properties(Age)++; ++nReject; } else { Properties(Age) = 0; (*it)->R = W.R; (*it)->Drift = drift; (*it)->Properties = Properties; H.get((*it)->E); ++nAccept; } } }
void VMCParticleByParticle::advanceWalkerByWalker() { MCWalkerConfiguration::iterator it(W.begin()),it_end(W.end()); while(it != it_end) { //MCWalkerConfiguration::WalkerData_t& w_buffer = *(W.DataSet[iwalker]); Walker_t& thisWalker(**it); Buffer_t& w_buffer(thisWalker.DataSet); W.R = thisWalker.R; w_buffer.rewind(); W.copyFromBuffer(w_buffer); Psi.copyFromBuffer(W,w_buffer); ValueType psi_old = thisWalker.Properties(SIGN); ValueType psi = psi_old; //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); bool moved = false; for(int iat=0; iat<W.getTotalNum(); iat++) { PosType dr = m_sqrttau*deltaR[iat]+thisWalker.Drift[iat]; PosType newpos = W.makeMove(iat,dr); //RealType ratio = Psi.ratio(W,iat); RealType ratio = Psi.ratio(W,iat,dG,dL); G = W.G+dG; //RealType forwardGF = exp(-0.5*dot(deltaR[iat],deltaR[iat])); //dr = (*it)->R[iat]-newpos-Tau*G[iat]; //RealType backwardGF = exp(-oneover2tau*dot(dr,dr)); RealType logGf = -0.5e0*dot(deltaR[iat],deltaR[iat]); ValueType vsq = Dot(G,G); ValueType scale = ((-1.0e0+sqrt(1.0e0+2.0e0*Tau*vsq))/vsq); dr = thisWalker.R[iat]-newpos-scale*G[iat]; //dr = (*it)->R[iat]-newpos-Tau*G[iat]; RealType logGb = -m_oneover2tau*dot(dr,dr); RealType prob = std::min(1.0e0,ratio*ratio*exp(logGb-logGf)); //alternatively if(Random() < prob) { moved = true; ++nAccept; W.acceptMove(iat); //Psi.update(W,iat); Psi.update2(W,iat); W.G = G; W.L += dL; //(*it)->Drift = Tau*G; thisWalker.Drift = scale*G; } else { ++nReject; Psi.restore(iat); } } if(moved) { w_buffer.rewind(); W.copyToBuffer(w_buffer); psi = Psi.evaluate(W,w_buffer); thisWalker.R = W.R; RealType eloc=H.evaluate(W); thisWalker.resetProperty(log(abs(psi)),psi,eloc); H.saveProperty(thisWalker.getPropertyBase()); } //Keep until everything is tested: debug routine // if(moved) { // //(*it)->Data.rewind(); // //W.copyToBuffer((*it)->Data); // //psi = Psi.evaluate(W,(*it)->Data); // DataSet[iwalker]->rewind(); // W.copyToBuffer(*DataSet[iwalker]); // psi = Psi.evaluate(W,*DataSet[iwalker]); // cout << "What is the ratio " << psi/psi_old << endl; // // std::cout << "Are they the same ratio=" << psi << endl; // RealType please = H.evaluate(W); // /**@note Debug the particle-by-particle move */ // G = W.G; // L = W.L; // DistanceTable::update(W); // ValueType psinew = Psi.evaluate(W); // ValueType dsum=pow(psi-psinew,2); // std::cout<< "Diff of updated and calculated gradients/laplacians" << endl; // for(int iat=0; iat<W.getTotalNum(); iat++) { // dsum += pow(G[iat][0]-W.G[iat][0],2)+pow(G[iat][1]-W.G[iat][1],2) // +pow(G[iat][2]-W.G[iat][2],2)+pow(L[iat]-W.L[iat],2); // std::cout << G[iat]-W.G[iat] << " " << L[iat]-W.L[iat] <<endl; // } // std::cout << "Are they the same ratio=" << psi << " eval=" << psinew << " " << sqrt(dsum) << endl; // std::cout << "============================ " << endl; // for(int i=0; i<W.getLocalNum(); i++) { // cout << W.G[i] << " " << W.L[i] << endl; // } // (*it)->Properties(PSISQ) = psi*psi; // (*it)->Properties(PSI) = psi; // (*it)->Properties(LOCALENERGY) = H.evaluate(W); // std::cout << "Energy " << please << " " << (*it)->Properties(LOCALENERGY) // << endl; // H.copy((*it)->getEnergyBase()); // ValueType vsq = Dot(W.G,W.G); // ValueType scale = ((-1.0+sqrt(1.0+2.0*Tau*vsq))/vsq); // (*it)->Drift = scale*W.G; // (*it)->R =W.R; // } else { ++nAllRejected; } ++it; } }
bool VMCPbyPMultiple::run() { Estimators->reportHeader(); //going to add routines to calculate how much we need bool require_register = W.createAuxDataSet(); multiEstimator->initialize(W,H1,Psi1,Tau,require_register); Estimators->reset(); //create an output engine HDFWalkerOutput WO(RootName); IndexType block = 0; Pooma::Clock timer; RealType oneovertau = 1.0/Tau; RealType oneover2tau = 0.5*oneovertau; RealType g = sqrt(Tau); RealType nPsi_minus_one = nPsi-1; MCWalkerConfiguration::iterator it; MCWalkerConfiguration::iterator it_end(W.end()); MCWalkerConfiguration::PropertyContainer_t Properties; ParticleSet::ParticleGradient_t dG(W.getTotalNum()); IndexType accstep=0; IndexType nAcceptTot = 0; IndexType nRejectTot = 0; do { //Blocks loop IndexType step = 0; timer.start(); nAccept = 0; nReject=0; IndexType nAllRejected = 0; do { //Steps loop it = W.begin(); int iwalker=0; while(it != it_end) { //Walkers loop MCWalkerConfiguration::WalkerData_t& w_buffer = *(W.DataSet[iwalker]); W.R = (*it)->R; w_buffer.rewind(); // Copy walker info in W W.copyFromBuffer(w_buffer); for(int ipsi=0; ipsi<nPsi; ipsi++){ // Copy wave function info in W and Psi1 Psi1[ipsi]->copyFromBuffer(W,w_buffer); Psi1[ipsi]->G=W.G; Psi1[ipsi]->L=W.L; } // Point to the correct walker in the ratioij buffer RealType *ratioijPtr=multiEstimator->RatioIJ[iwalker]; ValueType psi_old = (*it)->Properties(SIGN); ValueType psi = psi_old; //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); bool moved = false; for(int iat=0; iat<W.getTotalNum(); iat++) { //Particles loop PosType dr = g*deltaR[iat]+(*it)->Drift[iat]; PosType newpos = W.makeMove(iat,dr); for(int ipsi=0; ipsi<nPsi; ipsi++){ // Compute ratios before and after the move ratio[ipsi] = Psi1[ipsi]->ratio(W,iat,dG,*dL[ipsi]); // Compute Gradient in new position *G[ipsi]=Psi1[ipsi]->G + dG; // Initialize: sumratio[i]=(Psi[i]/Psi[i])^2=1.0 sumratio[ipsi]=1.0; } // Compute new (Psi[i]/Psi[j])^2 and their sum int indexij(0); for(int ipsi=0; ipsi< nPsi_minus_one; ipsi++){ for(int jpsi=ipsi+1; jpsi < nPsi; jpsi++){ RealType rji=ratio[jpsi]/ratio[ipsi]; rji = rji*rji*ratioijPtr[indexij]; ratioij[indexij++]=rji; sumratio[ipsi] += rji; sumratio[jpsi] += 1.0/rji; } } RealType logGf = -0.5*dot(deltaR[iat],deltaR[iat]); ValueType scale = Tau; drift=0.0; // Evaluate new Umbrella Weight and new drift for(int ipsi=0; ipsi< nPsi; ipsi++){ invsumratio[ipsi]=1.0/sumratio[ipsi]; drift += invsumratio[ipsi]*(*G[ipsi]); } drift *= scale; dr = (*it)->R[iat]-newpos-drift[iat]; RealType logGb = -oneover2tau*dot(dr,dr); // td = Target Density ratio RealType td=pow(ratio[0],2) *sumratio[0]/(*it)->Properties(SUMRATIO); RealType prob = std::min(1.0,td*exp(logGb-logGf)); if(Random() < prob) { /* Electron move is accepted. Update: -ratio (Psi[i]/Psi[j])^2 for this walker -Gradient and laplacian for each Psi1[i] -Drift -buffered info for each Psi1[i]*/ moved = true; ++nAccept; W.acceptMove(iat); // Update Buffer for (Psi[i]/Psi[j])^2 std::copy(ratioij.begin(),ratioij.end(),ratioijPtr); // Update Umbrella weight UmbrellaWeight=invsumratio; // Store sumratio for next Accept/Reject step (*it)->Properties(SUMRATIO)=sumratio[0]; for(int ipsi=0; ipsi< nPsi; ipsi++){ ////Update local Psi1[i] buffer for the next move Psi1[ipsi]->update2(W,iat); // Update G and L in Psi1[i] Psi1[ipsi]->G = *G[ipsi]; Psi1[ipsi]->L += *dL[ipsi]; } // Update Drift (*it)->Drift = drift; } else { ++nReject; for(int ipsi=0; ipsi< nPsi; ipsi++) Psi1[ipsi]->restore(iat); } } if(moved) { /* The walker moved: Info are copied back to buffers: -copy (Psi[i]/Psi[j])^2 to ratioijBuffer -Gradient and laplacian for each Psi1[i] -Drift -buffered info for each Psi1[i] Physical properties are updated */ (*it)->R = W.R; w_buffer.rewind(); W.copyToBuffer(w_buffer); for(int ipsi=0; ipsi< nPsi; ipsi++){ W.G=Psi1[ipsi]->G; W.L=Psi1[ipsi]->L; psi = Psi1[ipsi]->evaluate(W,w_buffer); RealType et = H1[ipsi]->evaluate(W); H1[ipsi]->copy((*it)->getEnergyBase(ipsi)); multiEstimator->updateSample(iwalker,ipsi,et,UmbrellaWeight[ipsi]); } } else { ++nAllRejected; } ++it; ++iwalker; } ++step;++accstep; Estimators->accumulate(W); } while(step<nSteps); timer.stop(); nAcceptTot += nAccept; nRejectTot += nReject; Estimators->flush(); Estimators->setColumn(AcceptIndex, static_cast<RealType>(nAccept)/static_cast<RealType>(nAccept+nReject)); Estimators->report(accstep); LogOut->getStream() << "Block " << block << " " << timer.cpu_time() << " Fixed_configs " << static_cast<RealType>(nAllRejected)/static_cast<RealType>(step*W.getActiveWalkers()) << " nPsi " << nPsi << endl; if(pStride) WO.get(W); nAccept = 0; nReject = 0; ++block; } while(block<nBlocks); LogOut->getStream() << "Ratio = " << static_cast<RealType>(nAcceptTot)/static_cast<RealType>(nAcceptTot+nRejectTot) << endl; if(!pStride) WO.get(W); Estimators->finalize(); return true; }
/** Advance all the walkers one timstep. */ void VMCMultiple::advanceWalkerByWalker() { m_oneover2tau = 0.5/Tau; m_sqrttau = std::sqrt(Tau); //MCWalkerConfiguration::PropertyContainer_t Properties; MCWalkerConfiguration::iterator it(W.begin()); MCWalkerConfiguration::iterator it_end(W.end()); int iwlk(0); int nPsi_minus_one(nPsi-1); while(it != it_end) { MCWalkerConfiguration::Walker_t &thisWalker(**it); //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); W.R = m_sqrttau*deltaR + thisWalker.R + thisWalker.Drift; //update the distance table associated with W //DistanceTable::update(W); W.update(); //Evaluate Psi and graidients and laplacians //\f$\sum_i \ln(|psi_i|)\f$ and catching the sign separately for(int ipsi=0; ipsi< nPsi; ipsi++) { logpsi[ipsi]=Psi1[ipsi]->evaluateLog(W); Psi1[ipsi]->L=W.L; Psi1[ipsi]->G=W.G; sumratio[ipsi]=1.0; } // Compute the sum over j of Psi^2[j]/Psi^2[i] for each i for(int ipsi=0; ipsi< nPsi_minus_one; ipsi++) { for(int jpsi=ipsi+1; jpsi< nPsi; jpsi++) { RealType ratioij=Norm[ipsi]/Norm[jpsi]*std::exp(2.0*(logpsi[jpsi]-logpsi[ipsi])); sumratio[ipsi] += ratioij; sumratio[jpsi] += 1.0/ratioij; } } for(int ipsi=0; ipsi< nPsi; ipsi++) invsumratio[ipsi]=1.0/sumratio[ipsi]; // Only these properties need to be updated // Using the sum of the ratio Psi^2[j]/Psi^2[iwref] // because these are number of order 1. Potentially // the sum of Psi^2[j] can get very big //Properties(LOGPSI) =logpsi[0]; //Properties(SUMRATIO) = sumratio[0]; RealType logGf = -0.5*Dot(deltaR,deltaR); //RealType scale = Tau; // ((-1.0+sqrt(1.0+2.0*Tau*vsq))/vsq); //accumulate the weighted drift: using operators in ParticleBase/ParticleAttribOps.h //to handle complex-to-real assignment and axpy operations. //drift = invsumratio[0]*Psi1[0]->G; //for(int ipsi=1; ipsi< nPsi ;ipsi++) { // drift += invsumratio[ipsi]*Psi1[ipsi]->G; //} PAOps<RealType,DIM>::scale(invsumratio[0],Psi1[0]->G,drift); for(int ipsi=1; ipsi< nPsi ; ipsi++) { PAOps<RealType,DIM>::axpy(invsumratio[ipsi],Psi1[ipsi]->G,drift); } //drift *= scale; setScaledDrift(Tau,drift); deltaR = thisWalker.R - W.R - drift; RealType logGb = -m_oneover2tau*Dot(deltaR,deltaR); //Original //RealType g = Properties(SUMRATIO)/thisWalker.Properties(SUMRATIO)* // exp(logGb-logGf+2.0*(Properties(LOGPSI)-thisWalker.Properties(LOGPSI))); //Reuse Multiplicity to store the sumratio[0] RealType g = sumratio[0]/thisWalker.Multiplicity* std::exp(logGb-logGf+2.0*(logpsi[0]-thisWalker.Properties(LOGPSI))); if(Random() > g) { thisWalker.Age++; ++nReject; } else { thisWalker.Age=0; thisWalker.Multiplicity=sumratio[0]; thisWalker.R = W.R; thisWalker.Drift = drift; for(int ipsi=0; ipsi<nPsi; ipsi++) { W.L=Psi1[ipsi]->L; W.G=Psi1[ipsi]->G; RealType et = H1[ipsi]->evaluate(W); thisWalker.Properties(ipsi,LOGPSI)=logpsi[ipsi]; thisWalker.Properties(ipsi,SIGN) =Psi1[ipsi]->getPhase(); thisWalker.Properties(ipsi,UMBRELLAWEIGHT)=invsumratio[ipsi]; thisWalker.Properties(ipsi,LOCALENERGY)=et; //multiEstimator->updateSample(iwlk,ipsi,et,invsumratio[ipsi]); H1[ipsi]->saveProperty(thisWalker.getPropertyBase(ipsi)); } ++nAccept; } ++it; ++iwlk; } }
bool VMCPbyPMultiple::run() { useDrift = (useDriftOpt=="yes"); if(useDrift) app_log() << " VMCPbyPMultiple::run useDrift=yes" << endl; else app_log() << " VMCPbyPMultiple::run useDrift=no" << endl; //TEST CACHE //Estimators->reportHeader(AppendRun); //going to add routines to calculate how much we need bool require_register = W.createAuxDataSet(); vector<RealType>Norm(nPsi),tmpNorm(nPsi); if(equilBlocks > 0) { for(int ipsi=0; ipsi< nPsi; ipsi++) { Norm[ipsi]=1.0; tmpNorm[ipsi]=0.0; } } else { for(int ipsi=0; ipsi< nPsi; ipsi++) Norm[ipsi]=std::exp(branchEngine->LogNorm[ipsi]); } multiEstimator->initialize(W,H1,Psi1,Tau,Norm,require_register); //TEST CACHE //Estimators->reset(); Estimators->start(nBlocks); //TEST CACHE IndexType block = 0; m_oneover2tau = 0.5/Tau; m_sqrttau = std::sqrt(Tau); RealType nPsi_minus_one = nPsi-1; ParticleSet::ParticleGradient_t dG(W.getTotalNum()); IndexType nAcceptTot = 0; IndexType nRejectTot = 0; MCWalkerConfiguration::iterator it; MCWalkerConfiguration::iterator it_end(W.end()); do //Blocks loop { IndexType step = 0; nAccept = 0; nReject=0; IndexType nAllRejected = 0; Estimators->startBlock(nSteps); do //Steps loop { it = W.begin(); int iwalker=0; while(it != it_end) //Walkers loop { Walker_t& thisWalker(**it); Walker_t::Buffer_t& w_buffer(thisWalker.DataSet); W.R = thisWalker.R; w_buffer.rewind(); // Copy walker info in W W.copyFromBuffer(w_buffer); for(int ipsi=0; ipsi<nPsi; ipsi++) { // Copy wave function info in W and Psi1 Psi1[ipsi]->copyFromBuffer(W,w_buffer); Psi1[ipsi]->G=W.G; Psi1[ipsi]->L=W.L; } // Point to the correct walker in the ratioij buffer RealType *ratioijPtr=multiEstimator->RatioIJ[iwalker]; //This is not used //ValueType psi_old = thisWalker.Properties(SIGN); //ValueType psi = psi_old; //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); bool moved = false; for(int iat=0; iat<W.getTotalNum(); iat++) //Particles loop { PosType dr = m_sqrttau*deltaR[iat]; if(useDrift) dr += thisWalker.Drift[iat]; PosType newpos = W.makeMove(iat,dr); for(int ipsi=0; ipsi<nPsi; ipsi++) { // Compute ratios before and after the move ratio[ipsi] = Psi1[ipsi]->ratio(W,iat,dG,*dL[ipsi]); logpsi2[ipsi]=std::log(ratio[ipsi]*ratio[ipsi]); // Compute Gradient in new position *G[ipsi]=Psi1[ipsi]->G + dG; // Initialize: sumratio[i]=(Psi[i]/Psi[i])^2=1.0 sumratio[ipsi]=1.0; } // Compute new (Psi[i]/Psi[j])^2 and their sum int indexij(0); for(int ipsi=0; ipsi< nPsi_minus_one; ipsi++) { for(int jpsi=ipsi+1; jpsi < nPsi; jpsi++, indexij++) { // Ratio between norms is already included in ratioijPtr from initialize. RealType rji=std::exp(logpsi2[jpsi]-logpsi2[ipsi])*ratioijPtr[indexij]; ratioij[indexij]=rji; sumratio[ipsi] += rji; sumratio[jpsi] += 1.0/rji; } } // Evaluate new Umbrella Weight for(int ipsi=0; ipsi< nPsi ; ipsi++) { invsumratio[ipsi]=1.0/sumratio[ipsi]; } RealType td=ratio[0]*ratio[0]*sumratio[0]/(*it)->Multiplicity; if(useDrift) { // Evaluate new drift PAOps<RealType,DIM>::scale(invsumratio[0],Psi1[0]->G,drift); for(int ipsi=1; ipsi< nPsi ; ipsi++) { PAOps<RealType,DIM>::axpy(invsumratio[ipsi],Psi1[ipsi]->G,drift); } setScaledDrift(Tau,drift); RealType logGf = -0.5*dot(deltaR[iat],deltaR[iat]); dr = thisWalker.R[iat]-newpos-drift[iat]; RealType logGb = -m_oneover2tau*dot(dr,dr); td *=std::exp(logGb-logGf); } // td = Target Density ratio //RealType td=pow(ratio[0],2)*sumratio[0]/(*it)->Properties(SUMRATIO); //td=ratio[0]*ratio[0]*sumratio[0]/(*it)->Multiplicity; //RealType prob = std::min(1.0,td*exp(logGb-logGf)); RealType prob = std::min(1.0,td); if(Random() < prob) { /* Electron move is accepted. Update: -ratio (Psi[i]/Psi[j])^2 for this walker -Gradient and laplacian for each Psi1[i] -Drift -buffered info for each Psi1[i]*/ moved = true; ++nAccept; W.acceptMove(iat); // Update Buffer for (Psi[i]/Psi[j])^2 std::copy(ratioij.begin(),ratioij.end(),ratioijPtr); // Update Umbrella weight UmbrellaWeight=invsumratio; // Store sumratio for next Accept/Reject step to Multiplicity //thisWalker.Properties(SUMRATIO)=sumratio[0]; thisWalker.Multiplicity=sumratio[0]; for(int ipsi=0; ipsi< nPsi; ipsi++) { //Update local Psi1[i] buffer for the next move Psi1[ipsi]->acceptMove(W,iat); //Update G and L in Psi1[i] Psi1[ipsi]->G = *G[ipsi]; Psi1[ipsi]->L += *dL[ipsi]; thisWalker.Properties(ipsi,LOGPSI)+=std::log(abs(ratio[ipsi])); } // Update Drift if(useDrift) (*it)->Drift = drift; } else { ++nReject; W.rejectMove(iat); for(int ipsi=0; ipsi< nPsi; ipsi++) Psi1[ipsi]->rejectMove(iat); } } if(moved) { /* The walker moved: Info are copied back to buffers: -copy (Psi[i]/Psi[j])^2 to ratioijBuffer -Gradient and laplacian for each Psi1[i] -Drift -buffered info for each Psi1[i] Physical properties are updated */ (*it)->Age=0; (*it)->R = W.R; w_buffer.rewind(); W.copyToBuffer(w_buffer); for(int ipsi=0; ipsi< nPsi; ipsi++) { W.G=Psi1[ipsi]->G; W.L=Psi1[ipsi]->L; ValueType psi = Psi1[ipsi]->evaluate(W,w_buffer); RealType et = H1[ipsi]->evaluate(W); //multiEstimator->updateSample(iwalker,ipsi,et,UmbrellaWeight[ipsi]); //Properties is used for UmbrellaWeight and UmbrellaEnergy thisWalker.Properties(ipsi,UMBRELLAWEIGHT)=UmbrellaWeight[ipsi]; thisWalker.Properties(ipsi,LOCALENERGY)=et; H1[ipsi]->auxHevaluate(W); H1[ipsi]->saveProperty(thisWalker.getPropertyBase(ipsi)); } } else { ++nAllRejected; } ++it; ++iwalker; } ++step; ++CurrentStep; Estimators->accumulate(W); } while(step<nSteps); //Modify Norm. if(block < equilBlocks) { for(int ipsi=0; ipsi< nPsi; ipsi++) { tmpNorm[ipsi]+=multiEstimator->esum(ipsi,MultipleEnergyEstimator::WEIGHT_INDEX); } if(block==(equilBlocks-1) || block==(nBlocks-1)) { RealType SumNorm(0.e0); for(int ipsi=0; ipsi< nPsi; ipsi++) SumNorm+=tmpNorm[ipsi]; for(int ipsi=0; ipsi< nPsi; ipsi++) { Norm[ipsi]=tmpNorm[ipsi]/SumNorm; branchEngine->LogNorm[ipsi]=std::log(Norm[ipsi]); } } } Estimators->stopBlock(static_cast<RealType>(nAccept)/static_cast<RealType>(nAccept+nReject)); nAcceptTot += nAccept; nRejectTot += nReject; nAccept = 0; nReject = 0; ++block; //record the current configuration recordBlock(block); //re-evaluate the ratio and update the Norm multiEstimator->initialize(W,H1,Psi1,Tau,Norm,false); } while(block<nBlocks); ////Need MPI-IO //app_log() // << "Ratio = " // << static_cast<RealType>(nAcceptTot)/static_cast<RealType>(nAcceptTot+nRejectTot) // << endl; return finalize(block); }
/** advance all the walkers with killnode==no * @param nat number of particles to move * * When killnode==no, any move resulting in node-crossing is treated * as a normal rejection. */ void DMCNonLocalUpdatePbyP::advanceWalkers(WalkerIter_t it, WalkerIter_t it_end, bool measure) { for(; it!=it_end; ++it) { Walker_t& thisWalker(**it); Walker_t::Buffer_t& w_buffer(thisWalker.DataSet); RealType eold(thisWalker.Properties(LOCALENERGY)); RealType vqold(thisWalker.Properties(DRIFTSCALE)); //RealType emixed(eold), enew(eold); RealType enew(eold); W.R = thisWalker.R; w_buffer.rewind(); W.copyFromBuffer(w_buffer); Psi.copyFromBuffer(W,w_buffer); //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); bool notcrossed(true); int nAcceptTemp(0); int nRejectTemp(0); RealType rr_proposed=0.0; RealType rr_accepted=0.0; for(int iat=0; iat<NumPtcl; ++iat) { //PosType dr(m_sqrttau*deltaR[iat]+thisWalker.Drift[iat]); RealType sc=getDriftScale(Tau,W.G[iat]); PosType dr(m_sqrttau*deltaR[iat]+sc*real(W.G[iat])); //RealType rr=dot(dr,dr); RealType rr=Tau*dot(deltaR[iat],deltaR[iat]); rr_proposed+=rr; if(rr>m_r2max)//too big move { ++nRejectTemp; continue; } PosType newpos(W.makeMove(iat,dr)); RealType ratio=Psi.ratio(W,iat,dG,dL); //node is crossed reject the move if(Psi.getPhase() > numeric_limits<RealType>::epsilon()) { ++nRejectTemp; ++nNodeCrossing; W.rejectMove(iat); Psi.rejectMove(iat); } else { G = W.G+dG; RealType logGf = -0.5*dot(deltaR[iat],deltaR[iat]); //RealType scale=getDriftScale(Tau,G); RealType scale=getDriftScale(Tau,G[iat]); dr = thisWalker.R[iat]-newpos-scale*real(G[iat]); RealType logGb = -m_oneover2tau*dot(dr,dr); RealType prob = std::min(1.0,ratio*ratio*std::exp(logGb-logGf)); if(RandomGen() < prob) { ++nAcceptTemp; W.acceptMove(iat); Psi.acceptMove(W,iat); W.G = G; W.L += dL; //assignDrift(scale,G,thisWalker.Drift); rr_accepted+=rr; } else { ++nRejectTemp; W.rejectMove(iat); Psi.rejectMove(iat); } } }//end of drift+diffusion for all the particles of a walker nonLocalOps.reset(); if(nAcceptTemp>0) {//need to overwrite the walker properties thisWalker.R = W.R; w_buffer.rewind(); W.copyToBuffer(w_buffer); RealType psi = Psi.evaluate(W,w_buffer); enew= H.evaluate(W,nonLocalOps.Txy); thisWalker.resetProperty(std::log(abs(psi)),psi,enew,rr_accepted,rr_proposed,1.0); H.saveProperty(thisWalker.getPropertyBase()); thisWalker.Drift=W.G; } else { thisWalker.Age++; thisWalker.Properties(R2ACCEPTED)=0.0; ++nAllRejected; enew=eold;//copy back old energy } int ibar = nonLocalOps.selectMove(RandomGen()); //make a non-local move if(ibar) { int iat=nonLocalOps.id(ibar); PosType newpos(W.makeMove(iat, nonLocalOps.delta(ibar))); RealType ratio=Psi.ratio(W,iat,dG,dL); W.acceptMove(iat); Psi.acceptMove(W,iat); W.G += dG; W.L += dL; PAOps<RealType,OHMMS_DIM>::copy(W.G,thisWalker.Drift); thisWalker.R[iat]=W.R[iat]; w_buffer.rewind(); W.copyToBuffer(w_buffer); RealType psi = Psi.evaluate(W,w_buffer); ++NonLocalMoveAccepted; } thisWalker.Weight *= branchEngine->branchWeight(eold,enew); nAccept += nAcceptTemp; nReject += nRejectTemp; } }
bool DMCParticleByParticle::run() { //add columns IndexType PopIndex = Estimators->addColumn("Population"); IndexType EtrialIndex = Estimators->addColumn("E_T"); //write the header Estimators->reportHeader(); MolecuFixedNodeBranch<RealType> brancher(Tau,W.getActiveWalkers()); //initialize parameters for fixed-node branching brancher.put(qmcNode,LogOut); if(BranchInfo != "default") brancher.read(BranchInfo); else { /*if VMC/DMC directly preceded DMC (Counter > 0) then use the average value of the energy estimator for the reference energy of the brancher*/ if(Counter) { RealType e_ref = W.getLocalEnergy(); LOGMSG("Overwriting the reference energy by the local energy " << e_ref) brancher.setEguess(e_ref); } } MCWalkerConfiguration::iterator it(W.begin()); MCWalkerConfiguration::iterator it_end(W.end()); while(it!=it_end) { (*it)->Properties(WEIGHT) = 1.0; (*it)->Properties(MULTIPLICITY) = 1.0; ++it; } //going to add routines to calculate how much we need bool require_register = W.createAuxDataSet(); int iwalker=0; it = W.begin(); it_end = W.end(); if(require_register) { while(it != it_end) { W.DataSet[iwalker]->rewind(); W.registerData(**it,*(W.DataSet[iwalker])); Psi.registerData(W,*(W.DataSet[iwalker])); ++it;++iwalker; } } Estimators->reset(); IndexType block = 0; Pooma::Clock timer; int Population = W.getActiveWalkers(); int tPopulation = W.getActiveWalkers(); RealType Eest = brancher.E_T; RealType E_T = Eest; RealType oneovertau = 1.0/Tau; RealType oneover2tau = 0.5*oneovertau; RealType g = sqrt(Tau); MCWalkerConfiguration::PropertyContainer_t Properties; ParticleSet::ParticlePos_t deltaR(W.getTotalNum()); ParticleSet::ParticleGradient_t G(W.getTotalNum()), dG(W.getTotalNum()); ParticleSet::ParticleLaplacian_t L(W.getTotalNum()), dL(W.getTotalNum()); IndexType accstep=0; IndexType nAcceptTot = 0; IndexType nRejectTot = 0; IndexType nat = W.getTotalNum(); int ncross = 0; do { IndexType step = 0; timer.start(); nAccept = 0; nReject=0; IndexType nAllRejected = 0; do { Population = W.getActiveWalkers(); it = W.begin(); it_end = W.end(); iwalker=0; while(it != it_end) { MCWalkerConfiguration::WalkerData_t& w_buffer = *(W.DataSet[iwalker]); (*it)->Properties(WEIGHT) = 1.0; (*it)->Properties(MULTIPLICITY) = 1.0; //save old local energy ValueType eold((*it)->Properties(LOCALENERGY)); ValueType emixed(eold); W.R = (*it)->R; w_buffer.rewind(); W.copyFromBuffer(w_buffer); Psi.copyFromBuffer(W,w_buffer); ValueType psi_old((*it)->Properties(SIGN)); ValueType psi(psi_old); //create a 3N-Dimensional Gaussian with variance=1 makeGaussRandom(deltaR); bool notcrossed(true); int nAcceptTemp(0); int nRejectTemp(0); int iat=0; while(notcrossed && iat<nat){ PosType dr(g*deltaR[iat]+(*it)->Drift[iat]); PosType newpos(W.makeMove(iat,dr)); RealType ratio(Psi.ratio(W,iat,dG,dL)); if(ratio < 0.0) {//node is crossed, stop here notcrossed = false; } else { G = W.G+dG; RealType logGf = -0.5*dot(deltaR[iat],deltaR[iat]); ValueType vsq = Dot(G,G); ValueType scale = ((-1.0+sqrt(1.0+2.0*Tau*vsq))/vsq); dr = (*it)->R[iat]-newpos-scale*G[iat]; //dr = (*it)->R[iat]-newpos-Tau*G[iat]; RealType logGb = -oneover2tau*dot(dr,dr); //RealType ratio2 = pow(ratio,2) RealType prob = std::min(1.0,pow(ratio,2)*exp(logGb-logGf)); if(Random() < prob) { ++nAcceptTemp; W.acceptMove(iat); Psi.update2(W,iat); W.G = G; W.L += dL; // (*it)->Drift = Tau*G; (*it)->Drift = scale*G; } else { ++nRejectTemp; Psi.restore(iat); } } ++iat; } if(notcrossed) { if(nAcceptTemp) {//need to overwrite the walker properties w_buffer.rewind(); W.copyToBuffer(w_buffer); psi = Psi.evaluate(W,w_buffer); (*it)->R = W.R; (*it)->Properties(AGE) = 0; //This is not so useful: allow overflow/underflow (*it)->Properties(LOGPSI) = log(fabs(psi)); (*it)->Properties(SIGN) = psi; (*it)->Properties(LOCALENERGY) = H.evaluate(W); H.copy((*it)->getEnergyBase()); (*it)->Properties(LOCALPOTENTIAL) = H.getLocalPotential(); emixed += (*it)->Properties(LOCALENERGY); } else { //WARNMSG("All the particle moves are rejected.") (*it)->Properties(AGE)++; ++nAllRejected; emixed += eold; } ValueType M = brancher.branchGF(Tau,emixed*0.5,0.0); // if((*it)->Properties(AGE) > 3.0) M = min(0.5,M); //persistent configurations if((*it)->Properties(AGE) > 1.9) M = std::min(0.5,M); if((*it)->Properties(AGE) > 0.9) M = std::min(1.0,M); (*it)->Properties(WEIGHT) = M; (*it)->Properties(MULTIPLICITY) = M + Random(); nAccept += nAcceptTemp; nReject += nRejectTemp; } else {//set the weight and multiplicity to zero (*it)->Properties(WEIGHT) = 0.0; (*it)->Properties(MULTIPLICITY) = 0.0; nReject += W.getTotalNum();//not make sense } ++it; ++iwalker; } ++step;++accstep; Estimators->accumulate(W); Eest = brancher.update(Population,Eest); //E_T = brancher.update(Population,Eest); brancher.branch(accstep,W); } while(step<nSteps); //WARNMSG("The number of a complete rejectoin " << nAllRejected) timer.stop(); nAcceptTot += nAccept; nRejectTot += nReject; Estimators->flush(); Estimators->setColumn(PopIndex,static_cast<RealType>(Population)); Estimators->setColumn(EtrialIndex,Eest); //E_T); Estimators->setColumn(AcceptIndex, static_cast<RealType>(nAccept)/static_cast<RealType>(nAccept+nReject)); Estimators->report(accstep); Eest = Estimators->average(0); LogOut->getStream() << "Block " << block << " " << timer.cpu_time() << " Fixed_configs " << static_cast<RealType>(nAllRejected)/static_cast<RealType>(step*W.getActiveWalkers()) << endl; if(pStride) { //create an output engine HDFWalkerOutput WO(RootName); WO.get(W); brancher.write(WO.getGroupID()); } nAccept = 0; nReject = 0; block++; } while(block<nBlocks); LogOut->getStream() << "Ratio = " << static_cast<double>(nAcceptTot)/static_cast<double>(nAcceptTot+nRejectTot) << endl; if(!pStride) { //create an output engine HDFWalkerOutput WO(RootName); WO.get(W); brancher.write(WO.getGroupID()); } Estimators->finalize(); return true; }
/** Make Metropolis move to the walkers and save in a temporary array. * @param it the iterator of the first walker to work on * @param tauinv inverse of the time step * * R + D + X */ void MCWalkerConfiguration::sample(iterator it, RealType tauinv) { makeGaussRandom(R); R *= tauinv; R += (*it)->R + (*it)->Drift; }
void ParticleSet::randomizeFromSource (ParticleSet &src) { SpeciesSet& srcSpSet(src.getSpeciesSet()); SpeciesSet& spSet(getSpeciesSet()); int srcChargeIndx = srcSpSet.addAttribute("charge"); int srcMemberIndx = srcSpSet.addAttribute("membersize"); int ChargeIndex = spSet.addAttribute("charge"); int MemberIndx = spSet.addAttribute("membersize"); int Nsrc = src.getTotalNum(); int Nptcl = getTotalNum(); int NumSpecies = spSet.TotalNum; int NumSrcSpecies = srcSpSet.TotalNum; //Store information about charges and number of each species vector<int> Zat, Zspec, NofSpecies, NofSrcSpecies, CurElec; Zat.resize(Nsrc); Zspec.resize(NumSrcSpecies); NofSpecies.resize(NumSpecies); CurElec.resize(NumSpecies); NofSrcSpecies.resize(NumSrcSpecies); for(int spec=0; spec<NumSrcSpecies; spec++) { Zspec[spec] = (int)round(srcSpSet(srcChargeIndx,spec)); NofSrcSpecies[spec] = (int)round(srcSpSet(srcMemberIndx,spec)); } for(int spec=0; spec<NumSpecies; spec++) { NofSpecies[spec] = (int)round(spSet(MemberIndx,spec)); CurElec[spec] = first(spec); } int totQ=0; for(int iat=0; iat<Nsrc; iat++) totQ+=Zat[iat] = Zspec[src.GroupID[iat]]; app_log() << " Total ion charge = " << totQ << endl; totQ -= Nptcl; app_log() << " Total system charge = " << totQ << endl; // Now, loop over ions, attaching electrons to them to neutralize // charge int spToken = 0; // This is decremented when we run out of electrons in each species int spLeft = NumSpecies; vector<PosType> gaussRand (Nptcl); makeGaussRandom (gaussRand); for (int iat=0; iat<Nsrc; iat++) { // Loop over electrons to add, selecting round-robin from the // electron species int z = Zat[iat]; while (z > 0 && spLeft) { int sp = spToken++ % NumSpecies; if (NofSpecies[sp]) { NofSpecies[sp]--; z--; int elec = CurElec[sp]++; app_log() << " Assigning " << (sp ? "down" : "up ") << " electron " << elec << " to ion " << iat << " with charge " << z << endl; double radius = 0.5* std::sqrt((double)Zat[iat]); R[elec] = src.R[iat] + radius * gaussRand[elec]; } else spLeft--; } } // Assign remaining electrons int ion=0; for (int sp=0; sp < NumSpecies; sp++) { for (int ie=0; ie<NofSpecies[sp]; ie++) { int iat = ion++ % Nsrc; double radius = std::sqrt((double)Zat[iat]); int elec = CurElec[sp]++; R[elec] = src.R[iat] + radius * gaussRand[elec]; } } }