QMCPACK
QMCUpdateBase.h
Go to the documentation of this file.
1 //////////////////////////////////////////////////////////////////////////////////////
2 // This file is distributed under the University of Illinois/NCSA Open Source License.
3 // See LICENSE file in top directory for details.
4 //
5 // Copyright (c) 2016 Jeongnim Kim and QMCPACK developers.
6 //
7 // File developed by: Jeremy McMinnis, jmcminis@gmail.com, University of Illinois at Urbana-Champaign
8 // Raymond Clay III, j.k.rofling@gmail.com, Lawrence Livermore National Laboratory
9 // Jeongnim Kim, jeongnim.kim@gmail.com, University of Illinois at Urbana-Champaign
10 // Jaron T. Krogel, krogeljt@ornl.gov, Oak Ridge National Laboratory
11 // Mark A. Berrill, berrillma@ornl.gov, Oak Ridge National Laboratory
12 //
13 // File created by: Jeongnim Kim, jeongnim.kim@gmail.com, University of Illinois at Urbana-Champaign
14 //////////////////////////////////////////////////////////////////////////////////////
15 
16 
17 /** @file
18  * @brief Declare QMCUpdateBase class
19  */
20 #ifndef QMCPLUSPLUS_QMCUPDATE_BASE_H
21 #define QMCPLUSPLUS_QMCUPDATE_BASE_H
27 #include "SimpleFixedNodeBranch.h"
28 #include "DriverDebugChecks.h"
30 
31 namespace qmcplusplus
32 {
33 class TraceManager;
34 class WalkerLogManager;
35 /** @ingroup QMC
36  * @brief Base class for update methods for each step
37  *
38  * QMCUpdateBase provides the common functions to update all the walkers for each time step.
39  * Derived classes should implement advanceWalkers to complete a step.
40  */
41 class QMCUpdateBase : public QMCTraits
42 {
43 public:
47 #ifdef MIXED_PRECISION
50 #else
51  using mPosType = PosType;
53 #endif
54 
55  ///number of steps per measurement
56  int nSubSteps;
57  /// determine additional checks for debugging purpose
59  std::string debug_checks_str_;
60  ///MaxAge>0 indicates branch is done
62  ///counter for number of moves accepted
64  ///counter for number of moves rejected
66  ///Total number of the steps when all the particle moves are rejected.
68  ///Total number of node crossings per block
70  ///Total numer of non-local moves accepted
72  ///timestep
74  ///spin mass
76  ///use Drift
77  bool UseDrift;
78 
79  /// Constructor.
81  ///Alt Constructor.
83  TrialWaveFunction& psi,
84  TrialWaveFunction& guide,
85  QMCHamiltonian& h,
87  ///destructor
88  virtual ~QMCUpdateBase();
89 
90  inline RealType acceptRatio() const
91  {
92  return static_cast<RealType>(nAccept) / static_cast<RealType>(nAccept + nReject);
93  }
94 
95  /** reset the QMCUpdateBase parameters
96  * @param brancher engine which handles branching
97  *
98  * Update time-step variables to move walkers
99  */
100  void resetRun(BranchEngineType* brancher,
102  TraceManager* traces,
103  const DriftModifierBase* driftmodifer);
104 
105  void resetRun2(BranchEngineType* brancher,
107  TraceManager* traces,
108  WalkerLogCollector* wlog_collector_,
109  const DriftModifierBase* driftmodifer);
110 
111  inline RealType getTau()
112  {
113  //SpeciesSet tspecies(W.getSpeciesSet());
114  //int massind=tspecies.addAttribute("mass");
115  //RealType mass = tspecies(massind,0);
116  //return m_tauovermass*mass;
117  return Tau;
118  }
119 
120  inline void setTau(RealType t)
121  {
122  //SpeciesSet tspecies(W.getSpeciesSet());
123  //int massind=tspecies.addAttribute("mass");
124  //RealType mass = tspecies(massind,0);
125  //RealType oneovermass = 1.0/mass;
126  //RealType oneoversqrtmass = std::sqrt(oneovermass);
127  // // Tau=brancher->getTau();
128  // // assert (Tau==i);
129  //m_tauovermass = i/mass;
130  Tau = t;
131  m_tauovermass = t * MassInvS[0];
132  m_oneover2tau = 0.5 / (m_tauovermass);
134  }
135 
136  inline RealType getSpinMass() { return spinMass; }
137 
138  inline void setSpinMass(RealType m) { spinMass = m; }
139 
140  inline void getLogs(std::vector<RealType>& logs) { Psi.getLogs(logs); }
141 
142  inline void set_step(int step) { W.current_step = step; }
143 
144 
145  ///** start a run */
146  void startRun(int blocks, bool record);
147  /** stop a run */
148  void stopRun();
149  void stopRun2();
150 
151  /** prepare to start a block
152  * @param steps number of steps within the block
153  */
154  void startBlock(int steps);
155 
156  /** stop a block
157  */
158  void stopBlock(bool collectall = true);
159 
160  /** set the multiplicity of the walkers to branch */
161  void setMultiplicity(WalkerIter_t it, WalkerIter_t it_end);
162 
163  inline void setMultiplicity(Walker_t& awalker) const
164  {
165  constexpr RealType onehalf(0.5);
166  constexpr RealType cone(1);
167  RealType M = awalker.Weight;
168  if (awalker.Age > MaxAge)
169  M = std::min(onehalf, M);
170  else if (awalker.Age > 0)
171  M = std::min(cone, M);
172  awalker.Multiplicity = M + RandomGen();
173  }
174 
175  /** initialize Walker buffers for PbyP update
176  */
177  virtual void initWalkersForPbyP(WalkerIter_t it, WalkerIter_t it_end);
178 
179  /** initialize Walker for walker update
180  */
181  virtual void initWalkers(WalkerIter_t it, WalkerIter_t it_end);
182 
183  /** process options
184  */
185  virtual bool put(xmlNodePtr cur);
186 
187  inline void accumulate(WalkerIter_t it, WalkerIter_t it_end) { Estimators->accumulate(W, it, it_end); }
188 
189  /** advance walkers executed at each step
190  *
191  * Derived classes implement how to move walkers and accept/reject
192  * moves.
193  */
194  virtual void advanceWalkers(WalkerIter_t it, WalkerIter_t it_end, bool recompute);
195 
196  ///move a walker
197  virtual void advanceWalker(Walker_t& thisWalker, bool recompute) = 0;
198 
200  std::vector<PosType>& dR,
201  std::vector<int>& iats,
202  std::vector<int>& rs,
203  std::vector<RealType>& ratios)
204  {
205  return 0.0;
206  };
207 
208  ///normalization offset for cs type runs.
210 
211  // virtual void estimateNormWalkers(std::vector<TrialWaveFunction*>& pclone
212  // , std::vector<MCWalkerConfiguration*>& wclone
213  // , std::vector<QMCHamiltonian*>& hclone
214  // , std::vector<RandomGenerator*>& rng
215  // , std::vector<RealType>& ratio_i_0){};
216  int RMC_checkIndex(int N, int NMax)
217  {
218  if (N < 0)
219  return N + NMax;
220  else if (N >= NMax)
221  return N - NMax;
222  else
223  return N;
224  }
225 
227  {
228  if (it >= last)
229  it -= (last - first);
230  else if (it < first)
231  it += (last - first);
232  }
233 
235  {
236  RealType logGb = 0.0;
237  for (int iat = 0; iat < W.getTotalNum(); ++iat)
238  {
239  RealType mass_over_tau = 1.0 / (SqrtTauOverMass[iat] * SqrtTauOverMass[iat]);
240  logGb += 0.5 * dot(displ[iat], displ[iat]) * mass_over_tau;
241  }
242  return -logGb;
243  }
244 
245 public:
246  ///traces
249 
250 protected:
251  ///update particle-by-particle
253  ///number of particles
255  ///Time-step factor \f$ 1/(2\tau)\f$
257  ///Time-step factor \f$ \sqrt{\tau}\f$
259  ///tau/mass
261  ///maximum displacement^2
263  ///walker ensemble
265  ///trial function
267  ///guide function
269  ///Hamiltonian
271  ///random number generator
273  ///branch engine, stateless reference to the one in QMCDriver
275  ///drift modifer, stateless reference to the one in QMCDriver
277  ///estimator
279  ///parameters
281  ///1/Mass per species
282  std::vector<RealType> MassInvS;
283  ///1/Mass per particle
284  std::vector<RealType> MassInvP;
285  ///sqrt(tau/Mass) per particle
286  std::vector<RealType> SqrtTauOverMass;
287  ///temporary storage for drift
289  ///temporary storage for random displacement
291  ///temporart storage for spin displacement
293  ///storage for differential gradients for PbyP update
295  ///storage for differential laplacians for PbyP update
297 
298  /** evaluate the ratio of scaled velocity and velocity
299  * @param g gradient
300  * @param gscaled scaled gradient
301  * @return the ratio
302  */
304 
305  ///copy constructor (disabled)
306  QMCUpdateBase(const QMCUpdateBase&) = delete;
307 
308  /// check logpsi and grad and lap against values computed from scratch
309  static void checkLogAndGL(ParticleSet& pset, TrialWaveFunction& twf, const std::string_view location);
310 
311 private:
312  ///set default parameters
313  void setDefaults();
314  /// Copy operator (disabled).
315  QMCUpdateBase& operator=(const QMCUpdateBase&) { return *this; }
316  ///
318 };
319 } // namespace qmcplusplus
320 
321 #endif
TraceManager * Traces
traces
const BranchEngineType * branchEngine
branch engine, stateless reference to the one in QMCDriver
ParticleSet::ParticlePos drift
temporary storage for drift
WalkerLogCollector * wlog_collector
RealType m_tauovermass
tau/mass
Base class for update methods for each step.
Definition: QMCUpdateBase.h:41
A set of walkers that are to be advanced by Metropolis Monte Carlo.
MCWalkerConfiguration::iterator WalkerIter_t
Definition: QMCUpdateBase.h:45
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
void getLogs(std::vector< RealType > &lvals)
void accumulate(MCWalkerConfiguration &W)
accumulate the measurements
void setSpinMass(RealType m)
TrialWaveFunction & Psi
trial function
RealType spinMass
spin mass
Definition: QMCUpdateBase.h:75
size_t getTotalNum() const
Definition: ParticleSet.h:493
ParticleSet::ParticleScalar deltaS
temporart storage for spin displacement
virtual bool put(xmlNodePtr cur)
process options
constexpr std::complex< float > cone
Definition: BLAS.hpp:50
Collection of Local Energy Operators.
ParticleAttrib< QTFull::ValueType > ParticleLaplacian
Definition: Configuration.h:96
Declaration of NonLocalTOperator.
this class implements drift modification
MCWalkerConfiguration::Walker_t Walker_t
Definition: QMCUpdateBase.h:44
Crowd-level resource for walker log collection.
Timer accumulates time and call counts.
Definition: NewTimer.h:135
RealType m_sqrttau
Time-step factor .
int current_step
current MC step
Definition: ParticleSet.h:134
ParticleSet::ParticlePos deltaR
temporary storage for random displacement
int Age
Age of this walker age is incremented when a walker is not moved after a sweep.
Definition: Walker.h:100
T min(T a, T b)
DriverDebugChecks debug_checks_
determine additional checks for debugging purpose
Definition: QMCUpdateBase.h:58
void startBlock(int steps)
prepare to start a block
IndexType nNodeCrossing
Total number of node crossings per block.
Definition: QMCUpdateBase.h:69
void resetRun(BranchEngineType *brancher, EstimatorManagerBase *est, TraceManager *traces, const DriftModifierBase *driftmodifer)
reset the QMCUpdateBase parameters
Specialized paritlce class for atomistic simulations.
Definition: ParticleSet.h:55
WalkerList_t::iterator iterator
FIX: a type alias of iterator for an object should not be for just one of many objects it holds...
Manages the state of QMC sections and handles population control for DMCs.
class to handle a set of parameters
Definition: ParameterSet.h:27
EstimatorManagerBase * Estimators
estimator
virtual void initWalkers(WalkerIter_t it, WalkerIter_t it_end)
initialize Walker for walker update
std::vector< RealType > SqrtTauOverMass
sqrt(tau/Mass) per particle
virtual void advanceWalker(Walker_t &thisWalker, bool recompute)=0
move a walker
RealType Tau
timestep
Definition: QMCUpdateBase.h:73
RealType logBackwardGF(const ParticleSet::ParticlePos &displ)
TrialWaveFunction & Guide
guide function
QTBase::PosType PosType
Definition: Configuration.h:61
void stopBlock(bool collectall=true)
stop a block
ParticleSet::ParticleLaplacian dL
MCWalkerConfiguration & W
walker ensemble
void setDefaults()
set default parameters
std::vector< RealType > MassInvP
1/Mass per particle
void getLogs(std::vector< RealType > &logs)
static void checkLogAndGL(ParticleSet &pset, TrialWaveFunction &twf, const std::string_view location)
check logpsi and grad and lap against values computed from scratch
Manager class of scalar estimators.
virtual RealType advanceWalkerForEE(Walker_t &w1, std::vector< PosType > &dR, std::vector< int > &iats, std::vector< int > &rs, std::vector< RealType > &ratios)
OHMMS_INDEXTYPE IndexType
define other types
Definition: Configuration.h:65
RealType getNodeCorrection(const ParticleSet::ParticleGradient &g, ParticleSet::ParticlePos &gscaled)
evaluate the ratio of scaled velocity and velocity
QMCUpdateBase(MCWalkerConfiguration &w, TrialWaveFunction &psi, QMCHamiltonian &h, RandomBase< FullPrecRealType > &rg)
Constructor.
int RMC_checkIndex(int N, int NMax)
WalkerConfigurations::Walker_t Walker_t
Class to manage a set of ScalarEstimators.
Declaration of a TrialWaveFunction.
int nSubSteps
number of steps per measurement
Definition: QMCUpdateBase.h:56
IndexType MaxAge
MaxAge>0 indicates branch is done.
Definition: QMCUpdateBase.h:61
void RMC_checkWalkerBounds(WalkerIter_t &it, WalkerIter_t first, WalkerIter_t last)
MakeReturn< UnaryNode< FnSqrt, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t sqrt(const Vector< T1, C1 > &l)
virtual ~QMCUpdateBase()
destructor
QMCUpdateBase & operator=(const QMCUpdateBase &)
Copy operator (disabled).
void setMultiplicity(Walker_t &awalker) const
Class to represent a many-body trial wave function.
RealType csoffset
normalization offset for cs type runs.
RandomBase< FullPrecRealType > & RandomGen
random number generator
ParticleAttrib< SingleParticlePos > ParticlePos
Definition: Configuration.h:92
virtual void initWalkersForPbyP(WalkerIter_t it, WalkerIter_t it_end)
initialize Walker buffers for PbyP update
FullPrecRealType Weight
Weight of the walker.
Definition: Walker.h:102
Tensor< typename BinaryReturn< T1, T2, OpMultiply >::Type_t, D > dot(const AntiSymTensor< T1, D > &lhs, const AntiSymTensor< T2, D > &rhs)
ParticleAttrib< QTFull::GradType > ParticleGradient
Definition: Configuration.h:95
declare a handler of DMC branching
IndexType nAccept
counter for number of moves accepted
Definition: QMCUpdateBase.h:63
RealType m_oneover2tau
Time-step factor .
IndexType nReject
counter for number of moves rejected
Definition: QMCUpdateBase.h:65
virtual void advanceWalkers(WalkerIter_t it, WalkerIter_t it_end, bool recompute)
advance walkers executed at each step
void accumulate(WalkerIter_t it, WalkerIter_t it_end)
traits for QMC variables
Definition: Configuration.h:49
ParticleSet::ParticleLaplacian L
storage for differential laplacians for PbyP update
ParticleSet::ParticleGradient dG
ParticleSet::ParticleGradient G
storage for differential gradients for PbyP update
FullPrecRealType Multiplicity
Number of copies for branching When Multiplicity = 0, this walker will be destroyed.
Definition: Walker.h:106
std::vector< RealType > MassInvS
1/Mass per species
void setMultiplicity(WalkerIter_t it, WalkerIter_t it_end)
set the multiplicity of the walkers to branch
IndexType NumPtcl
number of particles
QTBase::TensorType TensorType
Definition: Configuration.h:63
ParticleAttrib< Scalar_t > ParticleScalar
Definition: Configuration.h:91
Declaration of a MCWalkerConfiguration.
A container class to represent a walker.
Definition: Walker.h:49
QMCHamiltonian & H
Hamiltonian.
SimpleFixedNodeBranch BranchEngineType
Definition: QMCUpdateBase.h:46
ParameterSet myParams
parameters
RealType m_r2max
maximum displacement^2
RealType acceptRatio() const
Definition: QMCUpdateBase.h:90
void startRun(int blocks, bool record)
start a run
IndexType NonLocalMoveAccepted
Total numer of non-local moves accepted.
Definition: QMCUpdateBase.h:71
Declaration of QMCHamiltonian.
IndexType nAllRejected
Total number of the steps when all the particle moves are rejected.
Definition: QMCUpdateBase.h:67
void resetRun2(BranchEngineType *brancher, EstimatorManagerBase *est, TraceManager *traces, WalkerLogCollector *wlog_collector_, const DriftModifierBase *driftmodifer)
bool UpdatePbyP
update particle-by-particle
const DriftModifierBase * DriftModifier
drift modifer, stateless reference to the one in QMCDriver