QMCPACK
QMCDriver.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 // Jeongnim Kim, jeongnim.kim@gmail.com, University of Illinois at Urbana-Champaign
9 // Cynthia Gu, zg1@ornl.gov, Oak Ridge National Laboratory
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 /**
18  * @file QMCDriver.h
19  * @brief Declaration of QMCDriver
20  */
21 #ifndef QMCPLUSPLUS_QMCDRIVER_H
22 #define QMCPLUSPLUS_QMCDRIVER_H
23 
24 #include "Configuration.h"
25 #include "OhmmsData/ParameterSet.h"
26 #include "Pools/PooledData.h"
27 #include "Utilities/TimerManager.h"
29 #include "Utilities/ProjectData.h"
38 class Communicate;
39 
40 namespace qmcplusplus
41 {
42 /** @defgroup QMCDrivers QMC Driver group
43  * QMC drivers that implement QMC algorithms
44  */
45 
46 /** @defgroup WalkerByWalker QMC Drivers using walker-by-walker update
47  * @brief Move all the particles for each walker
48  */
49 
50 /** @defgroup ParticleByParticle QMC Drivers using particle-by-particle update
51  * @brief Move particle by particle
52  */
53 
54 /** @defgroup MultiplePsi QMC Drivers for energy differences
55  * @brief Umbrella sampling over multiple H/Psi
56  *
57  * This class of QMC drivers are suitable to evaluate
58  * the energy differences of multiple H-Psi pairs.
59  */
60 
61 //forward declarations: Do not include headers if not needed
62 class MCWalkerConfiguration;
63 class HDFWalkerOutput;
64 class TraceManager;
65 class WalkerLogManager;
66 
67 /** @ingroup QMCDrivers
68  * @{
69  * @brief abstract base class for QMC engines
70  */
71 class QMCDriver : public QMCDriverInterface, public QMCTraits, public MPIObjectBase
72 {
73 public:
74  /** enumeration coupled with QMCMode */
75  enum
76  {
81  };
82 
86 
87  /** bits to classify QMCDriver
88  *
89  * - qmc_driver_mode[QMC_UPDATE_MODE]? particle-by-particle: walker-by-walker
90  * - qmc_driver_mode[QMC_MULTIPLE]? multiple H/Psi : single H/Psi
91  * - qmc_driver_mode[QMC_OPTIMIZE]? optimization : vmc/dmc/rmc
92  */
93  std::bitset<QMC_MODE_MAX> qmc_driver_mode;
94 
95  /// whether to allow traces
97  /// traces xml
98  xmlNodePtr traces_xml;
99 
100  /// whether to allow traces
102  /// traces xml
103  xmlNodePtr walker_logs_xml;
104 
105  /// Constructor.
106  QMCDriver(const ProjectData& project_data,
108  TrialWaveFunction& psi,
109  QMCHamiltonian& h,
110  Communicate* comm,
111  const std::string& QMC_driver_type,
112  bool enable_profiling = false);
113 
114  ~QMCDriver() override;
115 
116  ///return current step
117  inline int current() const { return CurrentStep; }
118 
119  /** set the update mode
120  * @param pbyp if true, use particle-by-particle update
121  */
122  inline void setUpdateMode(bool pbyp) override { qmc_driver_mode[QMC_UPDATE_MODE] = pbyp; }
123 
124  /** Set the status of the QMCDriver
125  * @param aname the root file name
126  * @param h5name root name of the master hdf5 file containing previous qmcrun
127  * @param append if true, the run is a continuation of the previous qmc
128  *
129  * All output files will be of
130  * the form "aname.s00X.suffix", where "X" is number
131  * of previous QMC runs for the simulation and "suffix"
132  * is the suffix for the output file.
133  */
134  void setStatus(const std::string& aname, const std::string& h5name, bool append) override;
135 
136  /** add QMCHamiltonian/TrialWaveFunction pair for multiple
137  * @param h QMCHamiltonian
138  * @param psi TrialWaveFunction
139  *
140  * *Multiple* drivers use multiple H/Psi pairs to perform correlated sampling
141  * for energy difference evaluations.
142  */
143  void add_H_and_Psi(QMCHamiltonian* h, TrialWaveFunction* psi) override;
144 
145  /** initialize with xmlNode
146  */
147  void process(xmlNodePtr cur) override;
148 
149  /** return a xmlnode with update **/
150  xmlNodePtr getQMCNode();
151 
152  void putWalkers(std::vector<xmlNodePtr>& wset) override;
153 
154  inline void putTraces(xmlNodePtr txml) override { traces_xml = txml; }
155 
156  inline void requestTraces(bool traces) override { allow_traces = traces; }
157 
158  inline void putWalkerLogs(xmlNodePtr wlxml) override { walker_logs_xml = wlxml; }
159 
160  inline void requestWalkerLogs(bool allow_walker_logs_) override { allow_walker_logs = allow_walker_logs_; }
161 
162  std::string getEngineName() override { return QMCType; }
163 
164  template<class PDT>
165  void setValue(const std::string& aname, PDT x)
166  {
167  m_param.setValue(aname, x);
168  }
169 
170  ///set the BranchEngineType
171  void setBranchEngine(std::unique_ptr<BranchEngineType>&& be) override { branchEngine = std::move(be); }
172 
173  ///return BranchEngineType*
174  std::unique_ptr<BranchEngineType> getBranchEngine() override { return std::move(branchEngine); }
175 
176  int addObservable(const std::string& aname)
177  {
178  if (Estimators)
179  return Estimators->addObservable(aname.c_str());
180  else
181  return -1;
182  }
183 
185 
186  void setTau(RealType i) { Tau = i; }
187 
188  ///set global offsets of the walkers
189  void setWalkerOffsets();
190 
191  ///Observables manager
193 
194  ///Traces manager
195  std::unique_ptr<TraceManager> Traces;
196 
197  ///Traces manager
198  std::unique_ptr<WalkerLogManager> wlog_manager_;
199 
200  ///return the random generators
202  {
204  for (int i = 0; i < Rng.size(); ++i)
205  RngRefs.push_back(*Rng[i]);
206  return RngRefs;
207  }
208 
209  ///return the i-th random generator
210  inline RandomBase<FullPrecRealType>& getRng(int i) override { return (*Rng[i]); }
211 
212  unsigned long getDriverMode() override { return qmc_driver_mode.to_ulong(); }
213 
214 protected:
215  /// @brief top-level project data information
217  ///branch engine
218  std::unique_ptr<BranchEngineType> branchEngine;
219  ///drift modifer
221  ///randomize it
223  ///flag to append or restart the run
224  bool AppendRun;
225  ///flag to turn off dumping configurations
227  ///true, if it is a real QMC engine
229  /** the number of times this QMCDriver is executed
230  *
231  * MyCounter is initialized to zero by the constructor and is incremented
232  * whenever a run is completed by calling finalize(int block) or
233  * using MyCounter++ as in RQMC.
234  */
236  ///the number to delay updates by
237  int kDelay;
238  /** period of dumping walker configurations and everything else for restart
239  *
240  * The unit is a block.
241  */
243  /** period of dumping walker positions and IDs for Forward Walking
244  *
245  * The unit is in steps.
246  */
247 
248  ///Period to recalculate the walker properties from scratch.
250 
251  /** period of recording walker configurations
252  *
253  * Default is 0 indicating that only the last configuration will be saved.
254  */
256 
257  /** period of recording walker positions and IDs for forward walking afterwards
258  *
259  */
261 
262  ///current step
264 
265  ///maximum number of blocks
267 
268  ///maximum number of steps
270 
271  ///number of steps between a step: VMCs do not evaluate energies
273 
274  ///number of warmup steps
276 
277  ///counter for number of moves accepted
279 
280  ///counter for number of moves /rejected
282 
283  /// the number of blocks between recomptePsi
285 
286  ///the number of walkers
288  ///the number of saved samples
290  ///alternate method of setting QMC run parameters
292  ///samples per thread
294  ///target population
296 
297 
298  ///timestep
300 
301  ///maximum cpu in secs
303 
304  ///Time-step factor \f$ 1/(2\tau)\f$
306  ///Time-step factor \f$ \sqrt{\tau}\f$
308 
309  ///pointer to qmc node in xml file
310  xmlNodePtr qmcNode;
311 
312  ///type of QMC driver
313  const std::string QMCType;
314  ///the root of h5File
315  std::string h5FileRoot;
316  ///root of all the output files
317  std::string RootName;
318 
319  ///store any parameter that has to be read from a file
321 
322  ///walker ensemble
324 
325  ///trial function
327 
328  ///Hamiltonian
330 
331  ///record engine for walkers
332  std::unique_ptr<HDFWalkerOutput> wOut;
333 
334  ///a list of TrialWaveFunctions for multiple method
335  std::vector<TrialWaveFunction*> Psi1;
336 
337  ///a list of QMCHamiltonians for multiple method
338  std::vector<QMCHamiltonian*> H1;
339 
340  ///Random number generators
342 
343  ///a list of mcwalkerset element
344  std::vector<xmlNodePtr> mcwalkerNodePtr;
345 
346  ///temporary storage for drift
348 
349  ///temporary storage for random displacement
351 
352  ///spin mass for spinor calcs
354 
355  bool putQMCInfo(xmlNodePtr cur);
356 
357  void addWalkers(int nwalkers);
358 
359  /** record the state of the block
360  * @param block current block
361  *
362  * virtual function with a default implementation
363  */
364  void recordBlock(int block) override;
365 
366  /** finalize a qmc section
367  * @param block current block
368  * @param dumpwalkers if true, dump walkers
369  *
370  * Accumulate energy and weight is written to a hdf5 file.
371  * Finialize the estimators
372  */
373  bool finalize(int block, bool dumpwalkers = true);
374 
375  int rotation;
376  std::string getRotationName(std::string RootName);
377  std::string getLastRotationName(std::string RootName);
378  const std::string& get_root_name() const override { return RootName; }
379 
380 private:
382  ///profile the driver lifetime
384 };
385 /**@}*/
386 } // namespace qmcplusplus
387 
388 #endif
NewTimer & checkpoint_timer_
Definition: QMCDriver.h:381
std::string getRotationName(std::string RootName)
Definition: QMCDriver.cpp:278
std::string h5FileRoot
the root of h5File
Definition: QMCDriver.h:315
int MaxCPUSecs
maximum cpu in secs
Definition: QMCDriver.h:302
A set of walkers that are to be advanced by Metropolis Monte Carlo.
Base class for any object which needs to know about a MPI communicator.
Definition: MPIObjectBase.h:26
helper functions for EinsplineSetBuilder
Definition: Configuration.h:43
timer_manager class.
QTBase::RealType RealType
Definition: Configuration.h:58
QMCDriver(const ProjectData &project_data, MCWalkerConfiguration &w, TrialWaveFunction &psi, QMCHamiltonian &h, Communicate *comm, const std::string &QMC_driver_type, bool enable_profiling=false)
Constructor.
Definition: QMCDriver.cpp:44
std::string getLastRotationName(std::string RootName)
Definition: QMCDriver.cpp:293
std::vector< TrialWaveFunction * > Psi1
a list of TrialWaveFunctions for multiple method
Definition: QMCDriver.h:335
void recordBlock(int block) override
record the state of the block
Definition: QMCDriver.cpp:307
bool DumpConfig
flag to turn off dumping configurations
Definition: QMCDriver.h:226
std::string RootName
root of all the output files
Definition: QMCDriver.h:317
void setValue(const std::string &aname, PDT x)
Definition: QMCDriver.h:165
RealType m_oneover2tau
Time-step factor .
Definition: QMCDriver.h:305
RealType nSamplesPerThread
samples per thread
Definition: QMCDriver.h:293
xmlNodePtr qmcNode
pointer to qmc node in xml file
Definition: QMCDriver.h:310
void putWalkerLogs(xmlNodePtr wlxml) override
Definition: QMCDriver.h:158
class ProjectData
Definition: ProjectData.h:36
int addObservable(const std::string &aname)
Definition: QMCDriver.h:176
Collection of Local Energy Operators.
this class implements drift modification
ParticleSet::ParticlePos drift
temporary storage for drift
Definition: QMCDriver.h:347
IndexType CurrentStep
current step
Definition: QMCDriver.h:263
Timer accumulates time and call counts.
Definition: NewTimer.h:135
std::vector< std::unique_ptr< T > > UPtrVector
RandomBase< FullPrecRealType > & getRng(int i) override
return the i-th random generator
Definition: QMCDriver.h:210
abstract base class for QMC engines
Definition: QMCDriver.h:71
bool AppendRun
flag to append or restart the run
Definition: QMCDriver.h:224
RealType SpinMass
spin mass for spinor calcs
Definition: QMCDriver.h:353
IndexType nReject
counter for number of moves /rejected
Definition: QMCDriver.h:281
IndexType nTargetWalkers
the number of walkers
Definition: QMCDriver.h:287
std::unique_ptr< TraceManager > Traces
Traces manager.
Definition: QMCDriver.h:195
unsigned long getDriverMode() override
Definition: QMCDriver.h:212
void setBranchEngine(std::unique_ptr< BranchEngineType > &&be) override
set the BranchEngineType
Definition: QMCDriver.h:171
RealType getObservable(int i)
Definition: QMCDriver.h:184
int Period4ConfigDump
period of recording walker positions and IDs for forward walking afterwards
Definition: QMCDriver.h:260
ParameterSet m_param
store any parameter that has to be read from a file
Definition: QMCDriver.h:320
IndexType nSubSteps
number of steps between a step: VMCs do not evaluate energies
Definition: QMCDriver.h:272
Creates a common base class pointer for QMCDriver and QMCDriverNew to share.
Wrapping information on parallelism.
Definition: Communicate.h:68
IndexType nSteps
maximum number of steps
Definition: QMCDriver.h:269
std::unique_ptr< HDFWalkerOutput > wOut
record engine for walkers
Definition: QMCDriver.h:332
void setUpdateMode(bool pbyp) override
set the update mode
Definition: QMCDriver.h:122
bool finalize(int block, bool dumpwalkers=true)
finalize a qmc section
Definition: QMCDriver.cpp:318
const std::string & get_root_name() const override
Definition: QMCDriver.h:378
const std::string QMCType
type of QMC driver
Definition: QMCDriver.h:313
int kDelay
the number to delay updates by
Definition: QMCDriver.h:237
class to handle a set of parameters
Definition: ParameterSet.h:27
std::bitset< QMC_MODE_MAX > qmc_driver_mode
bits to classify QMCDriver
Definition: QMCDriver.h:93
Declaration of WaveFunctionPool.
std::unique_ptr< BranchEngineType > branchEngine
branch engine
Definition: QMCDriver.h:218
void putTraces(xmlNodePtr txml) override
Definition: QMCDriver.h:154
bool putQMCInfo(xmlNodePtr cur)
Parses the xml input file for parameter definitions for a single qmc simulation.
Definition: QMCDriver.cpp:399
void setValue(const std::string &aname_in, PDT aval)
QMCTraits::FullPrecRealType FullPrecRealType
Definition: QMCDriver.h:85
QMCHamiltonian & H
Hamiltonian.
Definition: QMCDriver.h:329
void setStatus(const std::string &aname, const std::string &h5name, bool append) override
Set the status of the QMCDriver.
Definition: QMCDriver.cpp:234
DriftModifierBase * DriftModifier
drift modifer
Definition: QMCDriver.h:220
MCWalkerConfiguration & W
walker ensemble
Definition: QMCDriver.h:323
RefVector< RandomBase< FullPrecRealType > > getRngRefs() const
return the random generators
Definition: QMCDriver.h:201
Manager class of scalar estimators.
OHMMS_INDEXTYPE IndexType
define other types
Definition: Configuration.h:65
std::unique_ptr< BranchEngineType > getBranchEngine() override
return BranchEngineType*
Definition: QMCDriver.h:174
WalkerConfigurations::Walker_t Walker_t
void setTau(RealType i)
Definition: QMCDriver.h:186
bool allow_walker_logs
whether to allow traces
Definition: QMCDriver.h:101
Class to manage a set of ScalarEstimators.
xmlNodePtr walker_logs_xml
traces xml
Definition: QMCDriver.h:103
Declaration of a TrialWaveFunction.
void add_H_and_Psi(QMCHamiltonian *h, TrialWaveFunction *psi) override
add QMCHamiltonian/TrialWaveFunction pair for multiple
Definition: QMCDriver.cpp:151
bool IsQMCDriver
true, if it is a real QMC engine
Definition: QMCDriver.h:228
void requestTraces(bool traces) override
Definition: QMCDriver.h:156
std::vector< std::reference_wrapper< T > > RefVector
std::string getEngineName() override
Definition: QMCDriver.h:162
Class to represent a many-body trial wave function.
std::vector< QMCHamiltonian * > H1
a list of QMCHamiltonians for multiple method
Definition: QMCDriver.h:338
IndexType nTargetSamples
the number of saved samples
Definition: QMCDriver.h:289
EstimatorManagerBase * Estimators
Observables manager.
Definition: QMCDriver.h:192
int Period4CheckPoint
period of dumping walker configurations and everything else for restart
Definition: QMCDriver.h:242
RealType nTargetPopulation
target population
Definition: QMCDriver.h:295
ParticleAttrib< SingleParticlePos > ParticlePos
Definition: Configuration.h:92
TrialWaveFunction & Psi
trial function
Definition: QMCDriver.h:326
std::vector< xmlNodePtr > mcwalkerNodePtr
a list of mcwalkerset element
Definition: QMCDriver.h:344
IndexType nAccept
counter for number of moves accepted
Definition: QMCDriver.h:278
RealType Tau
timestep
Definition: QMCDriver.h:299
declare a handler of DMC branching
const ProjectData & project_data_
top-level project data information
Definition: QMCDriver.h:216
IndexType nStepsBetweenSamples
alternate method of setting QMC run parameters
Definition: QMCDriver.h:291
IndexType nBlocks
maximum number of blocks
Definition: QMCDriver.h:266
ScopedProfiler driver_scope_profiler_
profile the driver lifetime
Definition: QMCDriver.h:383
QTFull::RealType FullPrecRealType
Definition: Configuration.h:66
int Period4WalkerDump
period of recording walker configurations
Definition: QMCDriver.h:255
int Period4CheckProperties
period of dumping walker positions and IDs for Forward Walking
Definition: QMCDriver.h:249
RealType m_sqrttau
Time-step factor .
Definition: QMCDriver.h:307
void addWalkers(int nwalkers)
Add walkers to the end of the ensemble of walkers.
Definition: QMCDriver.cpp:338
void setWalkerOffsets()
set global offsets of the walkers
Definition: QMCDriver.cpp:370
traits for QMC variables
Definition: Configuration.h:49
IndexType nWarmupSteps
number of warmup steps
Definition: QMCDriver.h:275
UPtrVector< RandomBase< FullPrecRealType > > Rng
Random number generators.
Definition: QMCDriver.h:341
std::unique_ptr< WalkerLogManager > wlog_manager_
Traces manager.
Definition: QMCDriver.h:198
IndexType nBlocksBetweenRecompute
the number of blocks between recomptePsi
Definition: QMCDriver.h:284
xmlNodePtr traces_xml
traces xml
Definition: QMCDriver.h:98
int current() const
return current step
Definition: QMCDriver.h:117
PooledData< RealType > Buffer_t
Definition: Walker.h:81
int MyCounter
the number of times this QMCDriver is executed
Definition: QMCDriver.h:235
ParticleSet::ParticlePos deltaR
temporary storage for random displacement
Definition: QMCDriver.h:350
void requestWalkerLogs(bool allow_walker_logs_) override
Definition: QMCDriver.h:160
A container class to represent a walker.
Definition: Walker.h:49
void putWalkers(std::vector< xmlNodePtr > &wset) override
Read walker configurations from *.config.h5 files.
Definition: QMCDriver.cpp:253
void process(xmlNodePtr cur) override
initialize with xmlNode
Definition: QMCDriver.cpp:172
Define a serialized buffer to store anonymous data.
bool allow_traces
whether to allow traces
Definition: QMCDriver.h:96
Declaration of QMCHamiltonian.
xmlNodePtr getQMCNode()
return a xmlnode with update
Definition: QMCDriver.cpp:516
bool ResetRandom
randomize it
Definition: QMCDriver.h:222