QMCPACK
QMCHamiltonian Class Reference

Collection of Local Energy Operators. More...

+ Collaboration diagram for QMCHamiltonian:

Classes

struct  QMCHamiltonianMultiWalkerResource
 

Public Types

enum  { DIM = OHMMS_DIM }
 
using Return_t = OperatorBase::Return_t
 
using PosType = OperatorBase::PosType
 
using TensorType = OperatorBase::TensorType
 
using RealType = OperatorBase::RealType
 
using ValueType = OperatorBase::ValueType
 
using FullPrecRealType = QMCTraits::FullPrecRealType
 
using PropertySetType = OperatorBase::PropertySetType
 
using BufferType = OperatorBase::BufferType
 
using Walker_t = OperatorBase::Walker_t
 
using WP = WalkerProperties::Indexes
 
using ValueMatrix = SPOSet::ValueMatrix
 

Public Member Functions

 QMCHamiltonian (const std::string &aname="psi0")
 constructor More...
 
 ~QMCHamiltonian ()
 destructor More...
 
void addOperator (std::unique_ptr< OperatorBase > &&h, const std::string &aname, bool physical=true)
 add an operator More...
 
void addOperatorType (const std::string &name, const std::string &type)
 record the name-type pair of an operator More...
 
const std::string & getOperatorType (const std::string &name)
 return type of named H element or fail More...
 
int size () const
 return the number of Hamiltonians More...
 
int total_size () const
 return the total number of Hamiltonians (physical + aux) More...
 
OperatorBasegetHamiltonian (const std::string &aname)
 return OperatorBase with the name aname More...
 
OperatorBasegetHamiltonian (int i)
 return i-th OperatorBase More...
 
RefVector< OperatorBasegetTWFDependentComponents ()
 return components, auxH not included, depending on TWF. More...
 
void initialize_traces (TraceManager &tm, ParticleSet &P)
 initialize trace data More...
 
void collect_walker_traces (Walker_t &walker, int step)
 collect walker trace data More...
 
void finalize_traces ()
 finalize trace data More...
 
int addObservables (ParticleSet &P)
 add each term to the PropertyList for averages More...
 
void registerObservables (std::vector< ObservableHelper > &h5desc, hdf_archive &file) const
 register obsevables so that their averages can be dumped to hdf5 More...
 
void registerCollectables (std::vector< ObservableHelper > &h5desc, hdf_archive &file) const
 register collectables so that their averages can be dumped to hdf5 More...
 
void informOperatorsOfListener ()
 Some Hamiltonian components need to be informed that they are in a per particle reporting situation so additional state can be added either to them or the objects they are strongly coupled with. More...
 
int startIndex () const
 retrun the starting index More...
 
int sizeOfObservables () const
 return the size of observables More...
 
int sizeOfCollectables () const
 return the size of collectables More...
 
RealType getObservable (int i) const
 return the value of the i-th observable More...
 
int getObservable (std::string Oname) const
 return the value of the observable with a set name if it exists More...
 
std::string getObservableName (int i) const
 return the name of the i-th observable More...
 
template<class IT , typename = std::enable_if_t<std::is_same<std::add_pointer<FullPrecRealType>::type, IT>::value>>
void saveProperty (IT first)
 save the values of Hamiltonian elements to the Properties More...
 
template<class IT , typename = std::enable_if_t<std::is_same<std::add_pointer<FullPrecRealType>::type, IT>::value>>
void setProperty (IT first)
 
void updateSource (ParticleSet &s)
 remove a named Hamiltonian from the list More...
 
FullPrecRealType getLocalEnergy ()
 
FullPrecRealType getLocalPotential ()
 
FullPrecRealType getKineticEnergy ()
 
void auxHevaluate (ParticleSet &P)
 
void auxHevaluate (ParticleSet &P, Walker_t &ThisWalker)
 This is more efficient. Only calculate auxH elements if moves are accepted. More...
 
void auxHevaluate (ParticleSet &P, Walker_t &ThisWalker, bool do_properties, bool do_collectables)
 Evaluate properties only. More...
 
void rejectedMove (ParticleSet &P, Walker_t &ThisWalker)
 Looks like a hack see DMCBatched.cpp and DMC.cpp weight is used like temporary flag from DMC. More...
 
void setPrimary (bool primary)
 set PRIMARY bit of all the components More...
 
FullPrecRealType evaluate (ParticleSet &P)
 evaluate Local Energy More...
 
FullPrecRealType evaluateDeterministic (ParticleSet &P)
 evaluate Local Energy deterministically. More...
 
FullPrecRealType evaluateWithToperator (ParticleSet &P)
 evaluate Local energy with Toperators updated. More...
 
FullPrecRealType evaluateValueAndDerivatives (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, Vector< ValueType > &dhpsioverpsi)
 evaluate energy and derivatives wrt to the variables More...
 
FullPrecRealType evaluateIonDerivsDeterministicFast (ParticleSet &P, ParticleSet &ions, TrialWaveFunction &psi_in, TWFFastDerivWrapper &psi_wrapper, ParticleSet::ParticlePos &dedr, ParticleSet::ParticlePos &wf_grad)
 evaluate local energy and derivatives w.r.t ionic coordinates, but deterministically. More...
 
void evaluateElecGrad (ParticleSet &P, TrialWaveFunction &psi, ParticleSet::ParticlePos &EGrad, RealType delta=1e-5)
 Evaluate the electron gradient of the local energy. More...
 
FullPrecRealType evaluateIonDerivs (ParticleSet &P, ParticleSet &ions, TrialWaveFunction &psi, ParticleSet::ParticlePos &hf_terms, ParticleSet::ParticlePos &pulay_terms, ParticleSet::ParticlePos &wf_grad)
 evaluate local energy and derivatives w.r.t ionic coordinates. More...
 
FullPrecRealType evaluateIonDerivsDeterministic (ParticleSet &P, ParticleSet &ions, TrialWaveFunction &psi, ParticleSet::ParticlePos &hf_terms, ParticleSet::ParticlePos &pulay_terms, ParticleSet::ParticlePos &wf_grad)
 evaluate local energy and derivatives w.r.t ionic coordinates, but deterministically. More...
 
void setNonLocalMoves (xmlNodePtr cur)
 set non local moves options More...
 
void setNonLocalMoves (const std::string &non_local_move_option, const double tau, const double alpha, const double gamma)
 
int makeNonLocalMoves (ParticleSet &P)
 make non local moves More...
 
bool has_L2 ()
 determine if L2 potential is present More...
 
void computeL2DK (ParticleSet &P, int iel, TensorType &D, PosType &K)
 compute D matrix and K vector for L2 potential propagator More...
 
void computeL2D (ParticleSet &P, int iel, TensorType &D)
 compute D matrix for L2 potential propagator More...
 
FullPrecRealType evaluateVariableEnergy (ParticleSet &P, bool free_nlpp)
 evaluate energy More...
 
FullPrecRealType getEnsembleAverage ()
 return an average value of the LocalEnergy More...
 
void resetTargetParticleSet (ParticleSet &P)
 
const std::string & getName () const
 
bool get (std::ostream &os) const
 
void setRandomGenerator (RandomBase< FullPrecRealType > *rng)
 
void createResource (ResourceCollection &collection) const
 initialize a shared resource and hand it to a collection More...
 
std::unique_ptr< QMCHamiltonianmakeClone (ParticleSet &qp, TrialWaveFunction &psi) const
 return a clone More...
 

Static Public Member Functions

static void mw_registerKineticListener (QMCHamiltonian &ham_leader, ListenerVector< RealType > listener)
 Listener Registration This must be called on a QMCHamiltonian that has acquired multiwalker resources. More...
 
static void mw_registerLocalEnergyListener (QMCHamiltonian &ham_leader, ListenerVector< RealType > listener)
 
static void mw_registerLocalPotentialListener (QMCHamiltonian &ham_leader, ListenerVector< RealType > listener)
 
static void mw_registerLocalIonPotentialListener (QMCHamiltonian &ham_leader, ListenerVector< RealType > listener)
 
static std::vector< QMCHamiltonian::FullPrecRealTypemw_evaluate (const RefVectorWithLeader< QMCHamiltonian > &ham_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list)
 batched version of evaluate for LocalEnergy More...
 
static std::vector< QMCHamiltonian::FullPrecRealTypemw_evaluateWithToperator (const RefVectorWithLeader< QMCHamiltonian > &ham_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list)
 batched version of evaluate Local energy with Toperators updated. More...
 
static std::vector< QMCHamiltonian::FullPrecRealTypemw_evaluateValueAndDerivatives (const RefVectorWithLeader< QMCHamiltonian > &ham_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list, const opt_variables_type &optvars, RecordArray< ValueType > &dlogpsi, RecordArray< ValueType > &dhpsioverpsi)
 
static std::vector< int > mw_makeNonLocalMoves (const RefVectorWithLeader< QMCHamiltonian > &ham_list, const RefVectorWithLeader< TrialWaveFunction > &wf_list, const RefVectorWithLeader< ParticleSet > &p_list)
 
static void updateComponent (OperatorBase &op, QMCHamiltonian &ham, ParticleSet &pset)
 accumulate local energy and update Observables and PropertyList More...
 
static void updateKinetic (QMCHamiltonian &ham, ParticleSet &pset)
 extract kinetic and potential energies. More...
 
static void acquireResource (ResourceCollection &collection, const RefVectorWithLeader< QMCHamiltonian > &ham_list)
 acquire external resource Note: use RAII ResourceCollectionLock whenever possible More...
 
static void releaseResource (ResourceCollection &collection, const RefVectorWithLeader< QMCHamiltonian > &ham_list)
 release external resource Note: use RAII ResourceCollectionLock whenever possible More...
 

Private Member Functions

void resetObservables (int start, int ncollects)
 reset Observables and counters More...
 
void reportToListeners ()
 

Static Private Member Functions

static RefVectorWithLeader< OperatorBaseextract_HC_list (const RefVectorWithLeader< QMCHamiltonian > &ham_list, int id)
 

Private Attributes

int myIndex
 starting index More...
 
int numCollectables
 starting index More...
 
FullPrecRealType LocalEnergy
 Current Local Energy. More...
 
FullPrecRealType KineticEnergy
 Current Kinetic Energy. More...
 
FullPrecRealType NewLocalEnergy
 Current Local Energy for the proposed move. More...
 
const std::string myName
 getName is in the way More...
 
std::vector< std::unique_ptr< OperatorBase > > H
 vector of Hamiltonians More...
 
NonLocalECPotentialnlpp_ptr
 pointer to NonLocalECP More...
 
L2Potentiall2_ptr
 pointer to L2Potential More...
 
std::vector< std::unique_ptr< OperatorBase > > auxH
 vector of Hamiltonians More...
 
NewTimerham_timer_
 Total timer for H evaluation. More...
 
NewTimereval_vals_derivs_timer_
 Total timer for H evaluation. More...
 
NewTimereval_ion_derivs_fast_timer_
 Total timer for H ion deriv evaluation;. More...
 
std::vector< std::reference_wrapper< NewTimer > > my_timers_
 timers for H components More...
 
std::map< std::string, std::string > operator_types
 types of component operators More...
 
PropertySetType Observables
 data More...
 
TraceRequest request
 traces variables More...
 
bool streaming_position
 
Array< TraceInt, 1 > * id_sample
 
Array< TraceInt, 1 > * pid_sample
 
Array< TraceInt, 1 > * step_sample
 
Array< TraceInt, 1 > * gen_sample
 
Array< TraceInt, 1 > * age_sample
 
Array< TraceInt, 1 > * mult_sample
 
Array< TraceReal, 1 > * weight_sample
 
Array< TraceReal, 2 > * position_sample
 
ResourceHandle< QMCHamiltonianMultiWalkerResourcemw_res_handle_
 

Static Private Attributes

static constexpr std::array< std::string_view, 8 > available_quantities_
 

Friends

class HamiltonianFactory
 

Detailed Description

Collection of Local Energy Operators.

Note that QMCHamiltonian is not derived from QMCHmailtonianBase.

Definition at line 49 of file QMCHamiltonian.h.

Member Typedef Documentation

◆ BufferType

Definition at line 61 of file QMCHamiltonian.h.

◆ FullPrecRealType

Definition at line 59 of file QMCHamiltonian.h.

◆ PosType

Definition at line 55 of file QMCHamiltonian.h.

◆ PropertySetType

using PropertySetType = OperatorBase::PropertySetType

Definition at line 60 of file QMCHamiltonian.h.

◆ RealType

Definition at line 57 of file QMCHamiltonian.h.

◆ Return_t

Definition at line 54 of file QMCHamiltonian.h.

◆ TensorType

Definition at line 56 of file QMCHamiltonian.h.

◆ ValueMatrix

Definition at line 64 of file QMCHamiltonian.h.

◆ ValueType

Definition at line 58 of file QMCHamiltonian.h.

◆ Walker_t

Definition at line 62 of file QMCHamiltonian.h.

◆ WP

Definition at line 63 of file QMCHamiltonian.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
DIM 

Definition at line 65 of file QMCHamiltonian.h.

66  {
67  DIM = OHMMS_DIM
68  };
#define OHMMS_DIM
Definition: config.h:64

Constructor & Destructor Documentation

◆ QMCHamiltonian()

QMCHamiltonian ( const std::string &  aname = "psi0")

constructor

Definition at line 49 of file QMCHamiltonian.cpp.

50  : myIndex(0),
51  numCollectables(0),
52  myName(aname),
53  nlpp_ptr(nullptr),
54  l2_ptr(nullptr),
55  ham_timer_(createGlobalTimer("Hamiltonian:" + aname + "::evaluate", timer_level_medium)),
56  eval_vals_derivs_timer_(createGlobalTimer("Hamiltonian:" + aname + "::ValueParamDerivs", timer_level_medium)),
58  createGlobalTimer("Hamiltonian:" + aname + ":::evaluateIonDerivsDeterministicFast", timer_level_medium))
59 #if !defined(REMOVE_TRACEMANAGER)
60  ,
61  streaming_position(false),
62  id_sample(nullptr),
63  pid_sample(nullptr),
64  step_sample(nullptr),
65  gen_sample(nullptr),
66  age_sample(nullptr),
67  mult_sample(nullptr),
68  weight_sample(nullptr),
69  position_sample(nullptr)
70 #endif
71 {}
Array< TraceInt, 1 > * gen_sample
NewTimer & ham_timer_
Total timer for H evaluation.
Array< TraceInt, 1 > * step_sample
Array< TraceInt, 1 > * pid_sample
NewTimer & eval_vals_derivs_timer_
Total timer for H evaluation.
NonLocalECPotential * nlpp_ptr
pointer to NonLocalECP
NewTimer & eval_ion_derivs_fast_timer_
Total timer for H ion deriv evaluation;.
Array< TraceReal, 1 > * weight_sample
int myIndex
starting index
int numCollectables
starting index
NewTimer & createGlobalTimer(const std::string &myname, timer_levels mylevel)
Array< TraceInt, 1 > * age_sample
L2Potential * l2_ptr
pointer to L2Potential
const std::string myName
getName is in the way
Array< TraceReal, 2 > * position_sample
Array< TraceInt, 1 > * mult_sample
Array< TraceInt, 1 > * id_sample

◆ ~QMCHamiltonian()

destructor

Definition at line 78 of file QMCHamiltonian.cpp.

79 {
80  //@todo clean up H and auxH
81 }

Member Function Documentation

◆ acquireResource()

void acquireResource ( ResourceCollection collection,
const RefVectorWithLeader< QMCHamiltonian > &  ham_list 
)
static

acquire external resource Note: use RAII ResourceCollectionLock whenever possible

Definition at line 1042 of file QMCHamiltonian.cpp.

References QMCHamiltonian::extract_HC_list(), RefVectorWithLeader< T >::getLeader(), and ResourceCollection::lendResource().

1044 {
1045  auto& ham_leader = ham_list.getLeader();
1046  ham_leader.mw_res_handle_ = collection.lendResource<QMCHamiltonianMultiWalkerResource>();
1047  for (int i_ham_op = 0; i_ham_op < ham_leader.H.size(); ++i_ham_op)
1048  {
1049  const auto HC_list(extract_HC_list(ham_list, i_ham_op));
1050  ham_leader.H[i_ham_op]->acquireResource(collection, HC_list);
1051  }
1052 }
static RefVectorWithLeader< OperatorBase > extract_HC_list(const RefVectorWithLeader< QMCHamiltonian > &ham_list, int id)

◆ addOperator()

void addOperator ( std::unique_ptr< OperatorBase > &&  h,
const std::string &  aname,
bool  physical = true 
)

add an operator

add a new Hamiltonian the the list of Hamiltonians.

Parameters
han operator
anamename of h
physicalif true, a physical operator

Definition at line 99 of file QMCHamiltonian.cpp.

References APP_ABORT, qmcplusplus::app_log(), qmcplusplus::app_warning(), QMCHamiltonian::auxH, qmcplusplus::createGlobalTimer(), QMCHamiltonian::getName(), QMCHamiltonian::H, QMCHamiltonian::l2_ptr, QMCHamiltonian::my_timers_, QMCHamiltonian::nlpp_ptr, OperatorBase::PHYSICAL, and qmcplusplus::timer_level_fine.

Referenced by ACForce::add2Hamiltonian(), OperatorBase::add2Hamiltonian(), CostFunctionCrowdData::CostFunctionCrowdData(), ECPotentialBuilder::put(), and qmcplusplus::TEST_CASE().

100 {
101  //change UpdateMode[PHYSICAL] of h so that cloning can be done correctly
102  h->getUpdateMode()[OperatorBase::PHYSICAL] = physical;
103  if (physical)
104  {
105  for (int i = 0; i < H.size(); ++i)
106  {
107  if (H[i]->getName() == aname)
108  {
109  app_warning() << "QMCHamiltonian::addOperator cannot " << aname << ". The name is already used" << std::endl;
110  return;
111  }
112  }
113  app_log() << " QMCHamiltonian::addOperator " << aname << " to H, physical Hamiltonian " << std::endl;
114  h->setName(aname);
115  H.push_back(std::move(h));
116  std::string tname = "Hamiltonian:" + aname;
117  my_timers_.push_back(createGlobalTimer(tname, timer_level_fine));
118  }
119  else
120  {
121  //ignore timers for now
122  for (int i = 0; i < auxH.size(); ++i)
123  {
124  if (auxH[i]->getName() == aname)
125  {
126  app_warning() << "QMCHamiltonian::addOperator cannot " << aname << ". The name is already used" << std::endl;
127  return;
128  }
129  }
130  app_log() << " QMCHamiltonian::addOperator " << aname << " to auxH " << std::endl;
131  h->setName(aname);
132  auxH.push_back(std::move(h));
133  }
134 
135  //assign save NLPP if found
136  // name is fixed in ECPotentialBuilder::put()
137  if (aname == "NonLocalECP")
138  {
139  if (nlpp_ptr == nullptr)
140  {
141  // original h arguments moved to either H or auxH
142  nlpp_ptr = physical ? dynamic_cast<NonLocalECPotential*>(H.back().get())
143  : dynamic_cast<NonLocalECPotential*>(auxH.back().get());
144  }
145  else
146  {
147  APP_ABORT("QMCHamiltonian::addOperator nlpp_ptr is supposed to be null. Something went wrong!");
148  }
149  }
150 
151  //save L2 potential if found
152  // name is fixed in ECPotentialBuilder::put()
153  if (aname == "L2")
154  {
155  if (l2_ptr == nullptr)
156  {
157  l2_ptr = physical ? dynamic_cast<L2Potential*>(H.back().get()) : dynamic_cast<L2Potential*>(auxH.back().get());
158  }
159  else
160  {
161  APP_ABORT("QMCHamiltonian::addOperator l2_ptr is supposed to be null. Something went wrong!");
162  }
163  }
164 }
std::ostream & app_warning()
Definition: OutputManager.h:69
std::ostream & app_log()
Definition: OutputManager.h:65
std::vector< std::reference_wrapper< NewTimer > > my_timers_
timers for H components
NonLocalECPotential * nlpp_ptr
pointer to NonLocalECP
NewTimer & createGlobalTimer(const std::string &myname, timer_levels mylevel)
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
L2Potential * l2_ptr
pointer to L2Potential
const std::string & getName() const
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ addOperatorType()

void addOperatorType ( const std::string &  name,
const std::string &  type 
)

record the name-type pair of an operator

Definition at line 167 of file QMCHamiltonian.cpp.

References qmcplusplus::app_log(), and QMCHamiltonian::operator_types.

168 {
169  app_log() << "QMCHamiltonian::addOperatorType added type " << type << " named " << name << std::endl;
170  operator_types[name] = type;
171 }
std::ostream & app_log()
Definition: OutputManager.h:65
std::map< std::string, std::string > operator_types
types of component operators

◆ auxHevaluate() [1/3]

void auxHevaluate ( ParticleSet P)

Definition at line 723 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, QMCHamiltonian::myIndex, QMCHamiltonian::Observables, and ParticleSet::PropertyList.

Referenced by SOVMCUpdateAll::advanceWalker(), SODMCUpdatePbyPWithRejectionFast::advanceWalker(), DMCUpdateAllWithRejection::advanceWalker(), VMCUpdateAll::advanceWalker(), VMCUpdatePbyP::advanceWalker(), DMCUpdatePbyPWithRejectionFast::advanceWalker(), SOVMCUpdatePbyP::advanceWalker(), DMCUpdatePbyPL2::advanceWalker(), DMCUpdateAllWithKill::advanceWalker(), VMCBatched::advanceWalkers(), RMCUpdatePbyPWithDrift::advanceWalkersRMC(), RMCUpdateAllWithDrift::advanceWalkersRMC(), RMCUpdatePbyPWithDrift::advanceWalkersVMC(), RMCUpdateAllWithDrift::advanceWalkersVMC(), QMCDriverNew::initialLogEvaluation(), RMCUpdateAllWithDrift::initWalkers(), QMCUpdateBase::initWalkersForPbyP(), and WaveFunctionTester::runGradSourceTest().

724 {
725  for (int i = 0; i < auxH.size(); ++i)
726  {
727  RealType sink = auxH[i]->evaluate(P);
728  auxH[i]->setObservables(Observables);
729 #if !defined(REMOVE_TRACEMANAGER)
730  auxH[i]->collectScalarTraces();
731 #endif
732  auxH[i]->setParticlePropertyList(P.PropertyList, myIndex);
733  //H[i]->setParticlePropertyList(P.PropertyList,myIndex);
734  }
735 }
int myIndex
starting index
QMCTraits::RealType RealType
PropertySetType Observables
data
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ auxHevaluate() [2/3]

void auxHevaluate ( ParticleSet P,
Walker_t ThisWalker 
)

This is more efficient. Only calculate auxH elements if moves are accepted.

Definition at line 738 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, QMCHamiltonian::collect_walker_traces(), ParticleSet::current_step, QMCHamiltonian::myIndex, QMCHamiltonian::Observables, and ParticleSet::PropertyList.

739 {
740 #if !defined(REMOVE_TRACEMANAGER)
741  collect_walker_traces(ThisWalker, P.current_step);
742 #endif
743  for (int i = 0; i < auxH.size(); ++i)
744  {
745  auxH[i]->setHistories(ThisWalker);
746  RealType sink = auxH[i]->evaluate(P);
747  auxH[i]->setObservables(Observables);
748 #if !defined(REMOVE_TRACEMANAGER)
749  auxH[i]->collectScalarTraces();
750 #endif
751  auxH[i]->setParticlePropertyList(P.PropertyList, myIndex);
752  }
753 }
int myIndex
starting index
void collect_walker_traces(Walker_t &walker, int step)
collect walker trace data
QMCTraits::RealType RealType
PropertySetType Observables
data
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ auxHevaluate() [3/3]

void auxHevaluate ( ParticleSet P,
Walker_t ThisWalker,
bool  do_properties,
bool  do_collectables 
)

Evaluate properties only.

Definition at line 755 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, QMCHamiltonian::collect_walker_traces(), OperatorBase::COLLECTABLE, ParticleSet::current_step, QMCHamiltonian::myIndex, QMCHamiltonian::Observables, and ParticleSet::PropertyList.

756 {
757 #if !defined(REMOVE_TRACEMANAGER)
758  collect_walker_traces(ThisWalker, P.current_step);
759 #endif
760  for (int i = 0; i < auxH.size(); ++i)
761  {
762  bool is_property = !(auxH[i]->getMode(OperatorBase::COLLECTABLE));
763  bool is_collectable = (auxH[i]->getMode(OperatorBase::COLLECTABLE));
764  if ((is_property && do_properties) || (is_collectable && do_collectables))
765  {
766  auxH[i]->setHistories(ThisWalker);
767  RealType sink = auxH[i]->evaluate(P);
768  auxH[i]->setObservables(Observables);
769 #if !defined(REMOVE_TRACEMANAGER)
770  auxH[i]->collectScalarTraces();
771 #endif
772  auxH[i]->setParticlePropertyList(P.PropertyList, myIndex);
773  }
774  }
775 }
int myIndex
starting index
void collect_walker_traces(Walker_t &walker, int step)
collect walker trace data
QMCTraits::RealType RealType
PropertySetType Observables
data
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ collect_walker_traces()

void collect_walker_traces ( Walker_t walker,
int  step 
)

collect walker trace data

Definition at line 491 of file QMCHamiltonian.cpp.

References QMCHamiltonian::age_sample, QMCHamiltonian::DIM, QMCHamiltonian::mult_sample, QMCHamiltonian::position_sample, QMCHamiltonian::request, TraceRequest::streaming_default_scalars, QMCHamiltonian::streaming_position, qmcplusplus::walker, and QMCHamiltonian::weight_sample.

Referenced by QMCHamiltonian::auxHevaluate(), and QMCHamiltonian::rejectedMove().

492 {
494  {
495  (*id_sample)(0) = walker.getWalkerID();
496  (*pid_sample)(0) = walker.getParentID();
497  (*step_sample)(0) = step;
498  (*gen_sample)(0) = walker.Generation;
499  (*age_sample)(0) = walker.Age;
500  (*mult_sample)(0) = walker.Multiplicity;
501  (*weight_sample)(0) = walker.Weight;
502  }
503  if (streaming_position)
504  for (int i = 0; i < walker.R.size(); ++i)
505  for (int d = 0; d < DIM; ++d)
506  (*position_sample)(i, d) = walker.R[i][d];
507 }
Array< TraceReal, 1 > * weight_sample
TraceRequest request
traces variables
Array< TraceInt, 1 > * age_sample
Array< TraceReal, 2 > * position_sample
Array< TraceInt, 1 > * mult_sample

◆ computeL2D()

void computeL2D ( ParticleSet P,
int  iel,
TensorType D 
)
inline

compute D matrix for L2 potential propagator

Parameters
rsingle particle coordinate
Ddiffusion matrix (outputted)

Definition at line 381 of file QMCHamiltonian.h.

References L2Potential::evaluateD(), and QMCHamiltonian::l2_ptr.

Referenced by DMCUpdatePbyPL2::advanceWalker().

382  {
383  if (l2_ptr != nullptr)
384  l2_ptr->evaluateD(P, iel, D);
385  }
void evaluateD(ParticleSet &P, int iel, TensorType &D)
L2Potential * l2_ptr
pointer to L2Potential

◆ computeL2DK()

void computeL2DK ( ParticleSet P,
int  iel,
TensorType D,
PosType K 
)
inline

compute D matrix and K vector for L2 potential propagator

Parameters
rsingle particle coordinate
Ddiffusion matrix (outputted)
Kdrift modification vector (outputted)

Definition at line 371 of file QMCHamiltonian.h.

References L2Potential::evaluateDK(), qmcplusplus::Units::energy::K, and QMCHamiltonian::l2_ptr.

Referenced by DMCUpdatePbyPL2::advanceWalker().

372  {
373  if (l2_ptr != nullptr)
374  l2_ptr->evaluateDK(P, iel, D, K);
375  }
L2Potential * l2_ptr
pointer to L2Potential
void evaluateDK(ParticleSet &P, int iel, TensorType &D, PosType &K)
Definition: L2Potential.cpp:97

◆ createResource()

void createResource ( ResourceCollection collection) const

initialize a shared resource and hand it to a collection

Definition at line 1035 of file QMCHamiltonian.cpp.

References ResourceCollection::addResource(), and QMCHamiltonian::H.

Referenced by CostFunctionCrowdData::CostFunctionCrowdData(), and QMCDriverNew::initializeQMC().

1036 {
1037  auto resource_index = collection.addResource(std::make_unique<QMCHamiltonianMultiWalkerResource>());
1038  for (int i = 0; i < H.size(); ++i)
1039  H[i]->createResource(collection);
1040 }
void createResource(ResourceCollection &collection) const
initialize a shared resource and hand it to a collection
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ evaluate()

evaluate Local Energy

Evaluate all the Hamiltonians for the N-particle configuration.

Parameters
PParticleSet
Returns
the local energy

P.R, P.G and P.L are used to evaluate the LocalEnergy.

Parameters
Pinput configuration containing N particles
Returns
the local energy

Definition at line 540 of file QMCHamiltonian.cpp.

References QMCHamiltonian::H, QMCHamiltonian::ham_timer_, QMCHamiltonian::LocalEnergy, QMCHamiltonian::my_timers_, QMCHamiltonian::updateComponent(), and QMCHamiltonian::updateKinetic().

Referenced by SOVMCUpdateAll::advanceWalker(), VMCUpdateAll::advanceWalker(), SOVMCUpdatePbyP::advanceWalker(), VMCUpdatePbyP::advanceWalker(), DMCUpdateAllWithKill::advanceWalker(), RMCUpdatePbyPWithDrift::advanceWalkersRMC(), RMCUpdateAllWithDrift::advanceWalkersRMC(), RMCUpdatePbyPWithDrift::advanceWalkersVMC(), RMCUpdateAllWithDrift::advanceWalkersVMC(), QMCHamiltonian::evaluateVariableEnergy(), QMCMain::executeCMCSection(), RMCUpdateAllWithDrift::initWalkers(), QMCUpdateBase::initWalkers(), RMCUpdatePbyPWithDrift::initWalkersForPbyP(), QMCUpdateBase::initWalkersForPbyP(), WaveFunctionTester::printEloc(), WaveFunctionTester::runCloneTest(), WaveFunctionTester::runDerivCloneTest(), WaveFunctionTester::runDerivNLPPTest(), WaveFunctionTester::runDerivTest(), WaveFunctionTester::runGradSourceTest(), WaveFunctionTester::runRatioTest(), WaveFunctionTester::runRatioTest2(), WaveFunctionTester::runZeroVarianceTest(), and qmcplusplus::test_hcpBe_rotation().

541 {
542  ScopedTimer local_timer(ham_timer_);
543  LocalEnergy = 0.0;
544  for (int i = 0; i < H.size(); ++i)
545  {
546  ScopedTimer h_timer(my_timers_[i]);
547  H[i]->evaluate(P);
548  updateComponent(*H[i], *this, P);
549 #if !defined(REMOVE_TRACEMANAGER)
550  H[i]->collectScalarTraces();
551 #endif
552  }
553  updateKinetic(*this, P);
554  return LocalEnergy;
555 }
NewTimer & ham_timer_
Total timer for H evaluation.
std::vector< std::reference_wrapper< NewTimer > > my_timers_
timers for H components
static void updateComponent(OperatorBase &op, QMCHamiltonian &ham, ParticleSet &pset)
accumulate local energy and update Observables and PropertyList
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
FullPrecRealType LocalEnergy
Current Local Energy.
static void updateKinetic(QMCHamiltonian &ham, ParticleSet &pset)
extract kinetic and potential energies.
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ evaluateDeterministic()

QMCHamiltonian::FullPrecRealType evaluateDeterministic ( ParticleSet P)

evaluate Local Energy deterministically.

Defaults to evaluate(P) for operators without a stochastic component. For the nonlocal PP, the quadrature grid is not rerandomized.

Parameters
PParticleSet
Returns
Local energy.

Definition at line 557 of file QMCHamiltonian.cpp.

References QMCHamiltonian::H, QMCHamiltonian::ham_timer_, QMCHamiltonian::LocalEnergy, QMCHamiltonian::my_timers_, QMCHamiltonian::updateComponent(), and QMCHamiltonian::updateKinetic().

Referenced by QMCHamiltonian::evaluateElecGrad(), and qmcplusplus::TEST_CASE().

558 {
559  ScopedTimer local_timer(ham_timer_);
560  LocalEnergy = 0.0;
561  for (int i = 0; i < H.size(); ++i)
562  {
563  ScopedTimer h_timer(my_timers_[i]);
564  H[i]->evaluateDeterministic(P);
565  updateComponent(*H[i], *this, P);
566 #if !defined(REMOVE_TRACEMANAGER)
567  H[i]->collectScalarTraces();
568 #endif
569  }
570  updateKinetic(*this, P);
571  return LocalEnergy;
572 }
NewTimer & ham_timer_
Total timer for H evaluation.
std::vector< std::reference_wrapper< NewTimer > > my_timers_
timers for H components
static void updateComponent(OperatorBase &op, QMCHamiltonian &ham, ParticleSet &pset)
accumulate local energy and update Observables and PropertyList
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
FullPrecRealType LocalEnergy
Current Local Energy.
static void updateKinetic(QMCHamiltonian &ham, ParticleSet &pset)
extract kinetic and potential energies.
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ evaluateElecGrad()

void evaluateElecGrad ( ParticleSet P,
TrialWaveFunction psi,
ParticleSet::ParticlePos EGrad,
RealType  delta = 1e-5 
)

Evaluate the electron gradient of the local energy.

Parameters
psiTrial Wave Function
Pelectron particle set
EGradan Nelec x 3 real array which corresponds to d/d[r_i]_j E_L
Afinite difference step size if applicable. Default is to use finite diff with delta=1e-5.
Returns
EGrad. Function itself returns nothing.

Definition at line 863 of file QMCHamiltonian.cpp.

References QMCHamiltonian::evaluateDeterministic(), TrialWaveFunction::evaluateLog(), ParticleSet::getTotalNum(), OHMMS_DIM, ParticleSet::R, and ParticleSet::update().

Referenced by ACForce::evaluate().

867 {
868  int nelec = P.getTotalNum();
869  RealType ep(0.0);
870  RealType em(0.0);
871  RealType e0(0.0);
872  for (int iel = 0; iel < nelec; iel++)
873  {
874  for (int dim = 0; dim < OHMMS_DIM; dim++)
875  {
876  RealType r0 = P.R[iel][dim];
877  ep = 0;
878  em = 0;
879  //Plus
880  RealType rp = r0 + delta;
881  P.R[iel][dim] = rp;
882  P.update();
883  psi.evaluateLog(P);
884  ep = evaluateDeterministic(P);
885 
886  //minus
887  RealType rm = r0 - delta;
888  P.R[iel][dim] = rm;
889  P.update();
890  psi.evaluateLog(P);
891  em = evaluateDeterministic(P);
892 
893  Egrad[iel][dim] = (ep - em) / (2.0 * delta);
894  P.R[iel][dim] = r0;
895  P.update();
896  psi.evaluateLog(P);
897  }
898  }
899 }
#define OHMMS_DIM
Definition: config.h:64
FullPrecRealType evaluateDeterministic(ParticleSet &P)
evaluate Local Energy deterministically.
QMCTraits::RealType RealType

◆ evaluateIonDerivs()

QMCHamiltonian::FullPrecRealType evaluateIonDerivs ( ParticleSet P,
ParticleSet ions,
TrialWaveFunction psi,
ParticleSet::ParticlePos hf_terms,
ParticleSet::ParticlePos pulay_terms,
ParticleSet::ParticlePos wf_grad 
)

evaluate local energy and derivatives w.r.t ionic coordinates.

Parameters
Ptarget particle set (electrons)
ionssource particle set (ions)
psiTrial wave function
hf_termsRe [(dH)Psi]/Psi
pulay_termsRe [(H-E_L)dPsi]/Psi
wf_gradRe (dPsi/Psi)
Returns
Local Energy.

Definition at line 900 of file QMCHamiltonian.cpp.

References qmcplusplus::convertToReal(), TrialWaveFunction::evalGradSource(), ParticleSet::getTotalNum(), and QMCHamiltonian::H.

Referenced by qmcplusplus::TEST_CASE().

906 {
907  ParticleSet::ParticleGradient wfgradraw_(ions.getTotalNum());
908  wfgradraw_ = 0.0;
909  RealType localEnergy = 0.0;
910 
911  for (int i = 0; i < H.size(); ++i)
912  localEnergy += H[i]->evaluateWithIonDerivs(P, ions, psi, hf_term, pulay_terms);
913 
914  for (int iat = 0; iat < ions.getTotalNum(); iat++)
915  {
916  wfgradraw_[iat] = psi.evalGradSource(P, ions, iat);
917  convertToReal(wfgradraw_[iat], wf_grad[iat]);
918  }
919  return localEnergy;
920 }
void convertToReal(const T1 &in, T2 &out)
generic conversion from type T1 to type T2 using implicit conversion
Definition: ConvertToReal.h:32
QMCTraits::RealType RealType
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
ParticleAttrib< QTFull::GradType > ParticleGradient
Definition: Configuration.h:95

◆ evaluateIonDerivsDeterministic()

QMCHamiltonian::FullPrecRealType evaluateIonDerivsDeterministic ( ParticleSet P,
ParticleSet ions,
TrialWaveFunction psi,
ParticleSet::ParticlePos hf_terms,
ParticleSet::ParticlePos pulay_terms,
ParticleSet::ParticlePos wf_grad 
)

evaluate local energy and derivatives w.r.t ionic coordinates, but deterministically.

Parameters
Ptarget particle set (electrons)
ionssource particle set (ions)
psiTrial wave function
hf_termsRe [(dH)Psi]/Psi
pulay_termsRe [(H-E_L)dPsi]/Psi
wf_gradRe (dPsi/Psi)
Returns
Local Energy.

Definition at line 922 of file QMCHamiltonian.cpp.

References qmcplusplus::convertToReal(), TrialWaveFunction::evalGradSource(), ParticleSet::getTotalNum(), and QMCHamiltonian::H.

Referenced by ACForce::evaluate(), and qmcplusplus::TEST_CASE().

928 {
929  ParticleSet::ParticleGradient wfgradraw_(ions.getTotalNum());
930  wfgradraw_ = 0.0;
931  RealType localEnergy = 0.0;
932 
933  for (int i = 0; i < H.size(); ++i)
934  localEnergy += H[i]->evaluateWithIonDerivsDeterministic(P, ions, psi, hf_term, pulay_terms);
935 
936  for (int iat = 0; iat < ions.getTotalNum(); iat++)
937  {
938  wfgradraw_[iat] = psi.evalGradSource(P, ions, iat);
939  convertToReal(wfgradraw_[iat], wf_grad[iat]);
940  }
941  return localEnergy;
942 }
void convertToReal(const T1 &in, T2 &out)
generic conversion from type T1 to type T2 using implicit conversion
Definition: ConvertToReal.h:32
QMCTraits::RealType RealType
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
ParticleAttrib< QTFull::GradType > ParticleGradient
Definition: Configuration.h:95

◆ evaluateIonDerivsDeterministicFast()

QMCHamiltonian::FullPrecRealType evaluateIonDerivsDeterministicFast ( ParticleSet P,
ParticleSet ions,
TrialWaveFunction psi_in,
TWFFastDerivWrapper psi_wrapper,
ParticleSet::ParticlePos dedr,
ParticleSet::ParticlePos wf_grad 
)

evaluate local energy and derivatives w.r.t ionic coordinates, but deterministically.

Parameters
Ptarget particle set (electrons)
ionssource particle set (ions)
psiTrial wave function
hf_termsRe [(dH)Psi]/Psi
pulay_termsRe [(H-E_L)dPsi]/Psi
wf_gradRe (dPsi/Psi)
Returns
Local Energy.

Definition at line 1097 of file QMCHamiltonian.cpp.

References TWFFastDerivWrapper::buildX(), TWFFastDerivWrapper::computeGSDerivative(), qmcplusplus::convertToReal(), QMCHamiltonian::eval_ion_derivs_fast_timer_, TWFFastDerivWrapper::evaluateJastrowGradSource(), ParticleSet::first(), TWFFastDerivWrapper::getGSMatrices(), TWFFastDerivWrapper::getIonGradM(), TWFFastDerivWrapper::getM(), ParticleSet::getTotalNum(), TWFFastDerivWrapper::getTWFGroupIndex(), QMCHamiltonian::H, TWFFastDerivWrapper::invertMatrices(), ParticleSet::last(), TWFFastDerivWrapper::numGroups(), TWFFastDerivWrapper::numOrbitals(), OHMMS_DIM, TWFFastDerivWrapper::trAB(), ParticleSet::update(), and TWFFastDerivWrapper::wipeMatrices().

Referenced by ACForce::evaluate(), and qmcplusplus::TEST_CASE().

1103 {
1105  P.update();
1106  //resize everything;
1107  const int ngroups = psi_wrapper_in.numGroups();
1108 
1109  std::vector<ValueMatrix> X_; //Working arrays for derivatives
1110  std::vector<ValueMatrix> Minv_; //Working array for derivatives.
1111  std::vector<ValueMatrix> B_;
1112  std::vector<ValueMatrix> B_gs_;
1113  std::vector<ValueMatrix> M_;
1114  std::vector<ValueMatrix> M_gs_;
1115 
1116  std::vector<std::vector<ValueMatrix>> dM_;
1117  std::vector<std::vector<ValueMatrix>> dM_gs_;
1118  std::vector<std::vector<ValueMatrix>> dB_;
1119  std::vector<std::vector<ValueMatrix>> dB_gs_;
1120 
1121  {
1122  M_.resize(ngroups);
1123  M_gs_.resize(ngroups);
1124  X_.resize(ngroups);
1125  B_.resize(ngroups);
1126  B_gs_.resize(ngroups);
1127  Minv_.resize(ngroups);
1128 
1129  for (int gid = 0; gid < ngroups; gid++)
1130  {
1131  const int sid = psi_wrapper_in.getTWFGroupIndex(gid);
1132  const int norbs = psi_wrapper_in.numOrbitals(sid);
1133  const int first = P.first(gid);
1134  const int last = P.last(gid);
1135  const int nptcls = last - first;
1136 
1137  M_[sid].resize(nptcls, norbs);
1138  B_[sid].resize(nptcls, norbs);
1139 
1140  M_gs_[sid].resize(nptcls, nptcls);
1141  Minv_[sid].resize(nptcls, nptcls);
1142  B_gs_[sid].resize(nptcls, nptcls);
1143  X_[sid].resize(nptcls, nptcls);
1144  }
1145 
1146  dM_.resize(OHMMS_DIM);
1147  dM_gs_.resize(OHMMS_DIM);
1148  dB_.resize(OHMMS_DIM);
1149  dB_gs_.resize(OHMMS_DIM);
1150 
1151  for (int idim = 0; idim < OHMMS_DIM; idim++)
1152  {
1153  dM_[idim].resize(ngroups);
1154  dB_[idim].resize(ngroups);
1155  dM_gs_[idim].resize(ngroups);
1156  dB_gs_[idim].resize(ngroups);
1157 
1158  for (int gid = 0; gid < ngroups; gid++)
1159  {
1160  const int sid = psi_wrapper_in.getTWFGroupIndex(gid);
1161  const int norbs = psi_wrapper_in.numOrbitals(sid);
1162  const int first = P.first(gid);
1163  const int last = P.last(gid);
1164  const int nptcls = last - first;
1165 
1166  dM_[idim][sid].resize(nptcls, norbs);
1167  dB_[idim][sid].resize(nptcls, norbs);
1168  dM_gs_[idim][sid].resize(nptcls, nptcls);
1169  dB_gs_[idim][sid].resize(nptcls, nptcls);
1170  }
1171  }
1172  psi_wrapper_in.wipeMatrices(M_);
1173  psi_wrapper_in.wipeMatrices(M_gs_);
1174  psi_wrapper_in.wipeMatrices(X_);
1175  psi_wrapper_in.wipeMatrices(B_);
1176  psi_wrapper_in.wipeMatrices(Minv_);
1177  psi_wrapper_in.wipeMatrices(B_gs_);
1178 
1179  for (int idim = 0; idim < OHMMS_DIM; idim++)
1180  {
1181  psi_wrapper_in.wipeMatrices(dM_[idim]);
1182  psi_wrapper_in.wipeMatrices(dM_gs_[idim]);
1183  psi_wrapper_in.wipeMatrices(dB_[idim]);
1184  psi_wrapper_in.wipeMatrices(dB_gs_[idim]);
1185  }
1186  }
1187  ParticleSet::ParticleGradient wfgradraw_(ions.getTotalNum());
1188  ParticleSet::ParticleGradient pulay_(ions.getTotalNum());
1189  ParticleSet::ParticleGradient hf_(ions.getTotalNum());
1190  ParticleSet::ParticleGradient dedr_complex(ions.getTotalNum());
1191  ParticleSet::ParticlePos pulayterms_(ions.getTotalNum());
1192  ParticleSet::ParticlePos hfdiag_(ions.getTotalNum());
1193  wfgradraw_ = 0.0;
1194  RealType localEnergy = 0.0;
1195 
1196  {
1197  psi_wrapper_in.getM(P, M_);
1198  }
1199  {
1200  psi_wrapper_in.getGSMatrices(M_, M_gs_);
1201  psi_wrapper_in.invertMatrices(M_gs_, Minv_);
1202  }
1203  //Build B-matrices. Only for non-diagonal observables right now.
1204  for (int i = 0; i < H.size(); ++i)
1205  {
1206  if (H[i]->dependsOnWaveFunction())
1207  {
1208  H[i]->evaluateOneBodyOpMatrix(P, psi_wrapper_in, B_);
1209  }
1210  else
1211  {
1212  localEnergy += H[i]->evaluateWithIonDerivsDeterministic(P, ions, psi_in, hfdiag_, pulayterms_);
1213  }
1214  }
1215 
1216  ValueType nondiag_cont = 0.0;
1217  RealType nondiag_cont_re = 0.0;
1218 
1219  psi_wrapper_in.getGSMatrices(B_, B_gs_);
1220  nondiag_cont = psi_wrapper_in.trAB(Minv_, B_gs_);
1221  convertToReal(nondiag_cont, nondiag_cont_re);
1222  localEnergy += nondiag_cont_re;
1223 
1224  {
1225  psi_wrapper_in.buildX(Minv_, B_gs_, X_);
1226  }
1227  //And now we compute the 3N force derivatives. 3 at a time for each atom.
1228  for (int iat = 0; iat < ions.getTotalNum(); iat++)
1229  {
1230  //The total wavefunction derivative has two contributions. One from determinantal piece,
1231  //One from the Jastrow. Jastrow is easy, so we evaluate it here, then add on the
1232  //determinantal piece at the end of this block.
1233 
1234  wfgradraw_[iat] = psi_wrapper_in.evaluateJastrowGradSource(P, ions, iat);
1235  for (int idim = 0; idim < OHMMS_DIM; idim++)
1236  {
1237  psi_wrapper_in.wipeMatrices(dM_[idim]);
1238  psi_wrapper_in.wipeMatrices(dM_gs_[idim]);
1239  psi_wrapper_in.wipeMatrices(dB_[idim]);
1240  psi_wrapper_in.wipeMatrices(dB_gs_[idim]);
1241  }
1242 
1243  {
1244  //ion derivative of slater matrix.
1245  psi_wrapper_in.getIonGradM(P, ions, iat, dM_);
1246  }
1247  for (int i = 0; i < H.size(); ++i)
1248  {
1249  if (H[i]->dependsOnWaveFunction())
1250  {
1251  H[i]->evaluateOneBodyOpMatrixForceDeriv(P, ions, psi_wrapper_in, iat, dB_);
1252  }
1253  }
1254 
1255  for (int idim = 0; idim < OHMMS_DIM; idim++)
1256  {
1257  psi_wrapper_in.getGSMatrices(dB_[idim], dB_gs_[idim]);
1258  psi_wrapper_in.getGSMatrices(dM_[idim], dM_gs_[idim]);
1259 
1260  ValueType fval = 0.0;
1261  fval = psi_wrapper_in.computeGSDerivative(Minv_, X_, dM_gs_[idim], dB_gs_[idim]);
1262  dedr_complex[iat][idim] = fval;
1263 
1264  ValueType wfcomp = 0.0;
1265  wfcomp = psi_wrapper_in.trAB(Minv_, dM_gs_[idim]);
1266  wfgradraw_[iat][idim] += wfcomp; //The determinantal piece of the WF grad.
1267  }
1268  convertToReal(dedr_complex[iat], dEdR[iat]);
1269  convertToReal(wfgradraw_[iat], wf_grad[iat]);
1270  }
1271  dEdR += hfdiag_;
1272  return localEnergy;
1273 }
void convertToReal(const T1 &in, T2 &out)
generic conversion from type T1 to type T2 using implicit conversion
Definition: ConvertToReal.h:32
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
#define OHMMS_DIM
Definition: config.h:64
NewTimer & eval_ion_derivs_fast_timer_
Total timer for H ion deriv evaluation;.
QMCTraits::RealType RealType
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
ParticleAttrib< SingleParticlePos > ParticlePos
Definition: Configuration.h:92
ParticleAttrib< QTFull::GradType > ParticleGradient
Definition: Configuration.h:95
OperatorBase::ValueType ValueType

◆ evaluateValueAndDerivatives()

QMCHamiltonian::FullPrecRealType evaluateValueAndDerivatives ( ParticleSet P,
const opt_variables_type optvars,
Vector< ValueType > &  dlogpsi,
Vector< ValueType > &  dhpsioverpsi 
)

evaluate energy and derivatives wrt to the variables

Parameters
PParticleSet
optvarscurrent optimiable variables
dlogpsi$\partial \ln \Psi({\bf R})/\partial \alpha $
dhpsioverpsi$\partial(\hat{h}\Psi({\bf R})/\Psi({\bf R})) /\partial \alpha $
compute_derivif true, compute dhpsioverpsi of the non-local potential component

Definition at line 649 of file QMCHamiltonian.cpp.

References QMCHamiltonian::eval_vals_derivs_timer_, QMCHamiltonian::H, QMCHamiltonian::KineticEnergy, QMCHamiltonian::LocalEnergy, and QMCHamiltonian::my_timers_.

Referenced by WaveFunctionTester::runDerivNLPPTest(), and qmcplusplus::test_hcpBe_rotation().

653 {
654  // The first componennt must be BareKineticEnergy for both handling KineticEnergy and dlogpsi computation
655  // by calling TWF::evaluateDerivatives inside BareKineticEnergy::evaluateValueAndDerivatives
656  assert(dynamic_cast<BareKineticEnergy*>(H[0].get()) &&
657  "BUG: The first componennt in Hamiltonian must be BareKineticEnergy.");
659 
660  {
661  ScopedTimer h_timer(my_timers_[0]);
662  LocalEnergy = KineticEnergy = H[0]->evaluateValueAndDerivatives(P, optvars, dlogpsi, dhpsioverpsi);
663  }
664 
665  for (int i = 1; i < H.size(); ++i)
666  {
667  ScopedTimer h_timer(my_timers_[i]);
668  LocalEnergy += H[i]->evaluateValueAndDerivatives(P, optvars, dlogpsi, dhpsioverpsi);
669  }
670  return LocalEnergy;
671 }
NewTimer & eval_vals_derivs_timer_
Total timer for H evaluation.
FullPrecRealType KineticEnergy
Current Kinetic Energy.
std::vector< std::reference_wrapper< NewTimer > > my_timers_
timers for H components
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
FullPrecRealType LocalEnergy
Current Local Energy.
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ evaluateVariableEnergy()

QMCHamiltonian::FullPrecRealType evaluateVariableEnergy ( ParticleSet P,
bool  free_nlpp 
)

evaluate energy

Parameters
Pquantum particleset
free_nlppif true, non-local PP is a variable
Returns
KE + NonLocal potential

Definition at line 710 of file QMCHamiltonian.cpp.

References QMCHamiltonian::evaluate(), and QMCHamiltonian::H.

Referenced by WaveFunctionTester::runDerivNLPPTest().

711 {
712  RealType nlpp = 0.0;
713  RealType ke = H[0]->evaluate(P);
714  if (free_nlpp)
715  for (int i = 1; i < H.size(); ++i)
716  {
717  if (H[i]->isNonLocal())
718  nlpp += H[i]->evaluate(P);
719  }
720  return ke + nlpp;
721 }
QMCTraits::RealType RealType
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
FullPrecRealType evaluate(ParticleSet &P)
evaluate Local Energy

◆ evaluateWithToperator()

QMCHamiltonian::FullPrecRealType evaluateWithToperator ( ParticleSet P)

evaluate Local energy with Toperators updated.

Parameters
PParticleSEt
Returns
Local energy

Definition at line 796 of file QMCHamiltonian.cpp.

References QMCHamiltonian::H, QMCHamiltonian::ham_timer_, QMCHamiltonian::LocalEnergy, QMCHamiltonian::my_timers_, QMCHamiltonian::updateComponent(), and QMCHamiltonian::updateKinetic().

Referenced by DMCUpdateAllWithRejection::advanceWalker(), SODMCUpdatePbyPWithRejectionFast::advanceWalker(), DMCUpdatePbyPWithRejectionFast::advanceWalker(), and DMCUpdatePbyPL2::advanceWalker().

797 {
798  ScopedTimer local_timer(ham_timer_);
799  LocalEnergy = 0.0;
800  for (int i = 0; i < H.size(); ++i)
801  {
802  ScopedTimer h_timer(my_timers_[i]);
803  H[i]->evaluateWithToperator(P);
804  updateComponent(*H[i], *this, P);
805 #if !defined(REMOVE_TRACEMANAGER)
806  H[i]->collectScalarTraces();
807 #endif
808  }
809  updateKinetic(*this, P);
810  return LocalEnergy;
811 }
NewTimer & ham_timer_
Total timer for H evaluation.
std::vector< std::reference_wrapper< NewTimer > > my_timers_
timers for H components
static void updateComponent(OperatorBase &op, QMCHamiltonian &ham, ParticleSet &pset)
accumulate local energy and update Observables and PropertyList
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
FullPrecRealType LocalEnergy
Current Local Energy.
static void updateKinetic(QMCHamiltonian &ham, ParticleSet &pset)
extract kinetic and potential energies.
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ extract_HC_list()

RefVectorWithLeader< OperatorBase > extract_HC_list ( const RefVectorWithLeader< QMCHamiltonian > &  ham_list,
int  id 
)
staticprivate

Definition at line 1087 of file QMCHamiltonian.cpp.

References RefVectorWithLeader< T >::getLeader(), and QMCHamiltonian::H.

Referenced by QMCHamiltonian::acquireResource(), QMCHamiltonian::mw_evaluate(), QMCHamiltonian::mw_evaluateValueAndDerivatives(), QMCHamiltonian::mw_evaluateWithToperator(), and QMCHamiltonian::releaseResource().

1089 {
1090  RefVectorWithLeader<OperatorBase> HC_list(*ham_list.getLeader().H[id]);
1091  HC_list.reserve(ham_list.size());
1092  for (QMCHamiltonian& H : ham_list)
1093  HC_list.push_back(*(H.H[id]));
1094  return HC_list;
1095 }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
QMCHamiltonian(const std::string &aname="psi0")
constructor

◆ finalize_traces()

void finalize_traces ( )

finalize trace data

Definition at line 510 of file QMCHamiltonian.cpp.

References QMCHamiltonian::age_sample, QMCHamiltonian::auxH, QMCHamiltonian::gen_sample, QMCHamiltonian::H, QMCHamiltonian::id_sample, QMCHamiltonian::mult_sample, QMCHamiltonian::pid_sample, QMCHamiltonian::position_sample, QMCHamiltonian::request, TraceRequest::reset(), QMCHamiltonian::step_sample, TraceRequest::streaming(), TraceRequest::streaming_default_scalars, QMCHamiltonian::streaming_position, and QMCHamiltonian::weight_sample.

Referenced by QMCUpdateBase::stopRun2().

511 {
513  {
514  delete id_sample;
515  delete pid_sample;
516  delete step_sample;
517  delete gen_sample;
518  delete age_sample;
519  delete mult_sample;
520  delete weight_sample;
521  }
522  if (streaming_position)
523  delete position_sample;
524  if (request.streaming())
525  {
526  for (int i = 0; i < H.size(); ++i)
527  H[i]->deleteTraceQuantities();
528  for (int i = 0; i < auxH.size(); ++i)
529  auxH[i]->deleteTraceQuantities();
530  }
531  streaming_position = false;
532  request.reset();
533 }
Array< TraceInt, 1 > * gen_sample
Array< TraceInt, 1 > * step_sample
Array< TraceInt, 1 > * pid_sample
bool streaming(const std::string &name)
Definition: TraceManager.h:256
Array< TraceReal, 1 > * weight_sample
TraceRequest request
traces variables
Array< TraceInt, 1 > * age_sample
Array< TraceReal, 2 > * position_sample
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians
Array< TraceInt, 1 > * mult_sample
Array< TraceInt, 1 > * id_sample

◆ get()

bool get ( std::ostream &  os) const

Definition at line 83 of file QMCHamiltonian.cpp.

References QMCHamiltonian::H.

84 {
85  for (int i = 0; i < H.size(); i++)
86  {
87  os << " " << std::setw(16) << std::left << H[i]->getName();
88  H[i]->get(os);
89  os << "\n";
90  }
91  return true;
92 }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ getEnsembleAverage()

QMCHamiltonian::FullPrecRealType getEnsembleAverage ( )

return an average value of the LocalEnergy

Introduced to get a collective value

Definition at line 944 of file QMCHamiltonian.cpp.

References QMCHamiltonian::H.

945 {
946  FullPrecRealType sum = 0.0;
947  for (int i = 0; i < H.size(); i++)
948  sum += H[i]->getEnsembleAverage();
949  return sum;
950 }
QMCTraits::FullPrecRealType FullPrecRealType
FullPrecRealType getEnsembleAverage()
return an average value of the LocalEnergy
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ getHamiltonian() [1/2]

OperatorBase * getHamiltonian ( const std::string &  aname)

return OperatorBase with the name aname

return pointer to the QMCHamtiltonian with the name

Parameters
anamename of a OperatorBase
Returns
0 if aname is not found.
Parameters
anamethe name of Hamiltonian
Returns
the pointer to the named term.

If not found, return 0

Definition at line 958 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, QMCHamiltonian::getName(), and QMCHamiltonian::H.

Referenced by qmcplusplus::TEST_CASE(), and qmcplusplus::test_hcpBe_rotation().

959 {
960  for (int i = 0; i < H.size(); ++i)
961  if (H[i]->getName() == aname)
962  return H[i].get();
963  for (int i = 0; i < auxH.size(); ++i)
964  if (auxH[i]->getName() == aname)
965  return auxH[i].get();
966  return nullptr;
967 }
const std::string & getName() const
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ getHamiltonian() [2/2]

OperatorBase* getHamiltonian ( int  i)
inline

return i-th OperatorBase

Parameters
iindex of the OperatorBase
Returns
H[i]

Definition at line 101 of file QMCHamiltonian.h.

References QMCHamiltonian::H.

101 { return H[i].get(); }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ getKineticEnergy()

FullPrecRealType getKineticEnergy ( )
inline

Definition at line 216 of file QMCHamiltonian.h.

References QMCHamiltonian::KineticEnergy.

Referenced by qmcplusplus::test_hcpBe_rotation().

216 { return KineticEnergy; }
FullPrecRealType KineticEnergy
Current Kinetic Energy.

◆ getLocalEnergy()

◆ getLocalPotential()

FullPrecRealType getLocalPotential ( )
inline

Definition at line 215 of file QMCHamiltonian.h.

References QMCHamiltonian::KineticEnergy, and QMCHamiltonian::LocalEnergy.

Referenced by WaveFunctionTester::runDerivCloneTest(), and WaveFunctionTester::runDerivTest().

215 { return LocalEnergy - KineticEnergy; }
FullPrecRealType KineticEnergy
Current Kinetic Energy.
FullPrecRealType LocalEnergy
Current Local Energy.

◆ getName()

const std::string& getName ( ) const
inline

Definition at line 405 of file QMCHamiltonian.h.

References QMCHamiltonian::myName.

Referenced by QMCHamiltonian::addOperator(), and QMCHamiltonian::getHamiltonian().

405 { return myName; }
const std::string myName
getName is in the way

◆ getOperatorType()

const std::string & getOperatorType ( const std::string &  name)

return type of named H element or fail

Definition at line 174 of file QMCHamiltonian.cpp.

References APP_ABORT, and QMCHamiltonian::operator_types.

Referenced by qmcplusplus::TEST_CASE().

175 {
176  std::map<std::string, std::string>::iterator type = operator_types.find(name);
177  if (type == operator_types.end())
178  APP_ABORT("QMCHamiltonian::getOperatorType\n operator type not found for name " + name);
179  return type->second;
180 }
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
std::map< std::string, std::string > operator_types
types of component operators

◆ getTWFDependentComponents()

RefVector< OperatorBase > getTWFDependentComponents ( )

return components, auxH not included, depending on TWF.

Definition at line 969 of file QMCHamiltonian.cpp.

References QMCHamiltonian::H.

Referenced by QMCCostFunctionBatched::correlatedSampling(), CostFunctionCrowdData::CostFunctionCrowdData(), and QMCCostFunctionBatched::getConfigurations().

970 {
971  RefVector<OperatorBase> components;
972  for (int i = 0; i < H.size(); i++)
973  if (H[i]->dependsOnWaveFunction())
974  components.push_back(*H[i]);
975  return components;
976 }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ has_L2()

bool has_L2 ( )
inline

determine if L2 potential is present

Definition at line 364 of file QMCHamiltonian.h.

References QMCHamiltonian::l2_ptr.

Referenced by DMCUpdatePbyPL2::advanceWalker().

364 { return l2_ptr != nullptr; }
L2Potential * l2_ptr
pointer to L2Potential

◆ initialize_traces()

void initialize_traces ( TraceManager tm,
ParticleSet P 
)

initialize trace data

Definition at line 310 of file QMCHamiltonian.cpp.

References QMCHamiltonian::age_sample, APP_ABORT, qmcplusplus::app_log(), QMCHamiltonian::auxH, TraceManager::checkout_int(), TraceManager::checkout_real(), TraceRequest::contribute_array(), TraceRequest::contribute_combined(), TraceRequest::contribute_scalar(), TraceRequest::determine_stream_write(), QMCHamiltonian::DIM, QMCHamiltonian::gen_sample, OperatorBase::getName(), QMCHamiltonian::H, QMCHamiltonian::id_sample, TraceRequest::incorporate(), OperatorBase::isClassical(), OperatorBase::isClassicalClassical(), OperatorBase::isQuantum(), OperatorBase::isQuantumClassical(), OperatorBase::isQuantumQuantum(), TraceManager::make_combined_trace(), QMCHamiltonian::mult_sample, omp_get_thread_num(), QMCHamiltonian::pid_sample, QMCHamiltonian::position_sample, TraceRequest::relay_stream_info(), QMCHamiltonian::request, TraceManager::request, TraceManager::screen_writes(), QMCHamiltonian::step_sample, TraceRequest::streaming(), TraceRequest::streaming_array(), TraceRequest::streaming_default_scalars, QMCHamiltonian::streaming_position, TraceManager::streaming_traces, TraceManager::update_status(), TraceManager::user_report(), TraceManager::verbose, and QMCHamiltonian::weight_sample.

Referenced by QMCUpdateBase::startRun().

311 {
312  static bool first_init = true;
313  bool trace_log = first_init && tm.verbose && omp_get_thread_num() == 0;
314  if (trace_log)
315  app_log() << "\n Hamiltonian is initializing traces" << std::endl;
316 
317  //fill std::string vectors for combined trace quantities
318  std::vector<std::string> Eloc;
319  std::vector<std::string> Vloc;
320  std::vector<std::string> Vq, Vc, Vqq, Vqc, Vcc;
321  for (int i = 0; i < H.size(); ++i)
322  Eloc.push_back(H[i]->getName());
323  for (int i = 1; i < H.size(); ++i)
324  Vloc.push_back(H[i]->getName());
325 
326  // These contributions are based on the potential energy components.
327  // Loop starts at one to skip the kinetic energy component.
328  for (int i = 1; i < H.size(); ++i)
329  {
330  OperatorBase& h = *H[i];
331  if (h.isQuantum())
332  Vq.push_back(h.getName());
333  else if (h.isClassical())
334  Vc.push_back(h.getName());
335  else if (h.isQuantumQuantum())
336  Vqq.push_back(h.getName());
337  else if (h.isQuantumClassical())
338  Vqc.push_back(h.getName());
339  else if (h.isClassicalClassical())
340  Vcc.push_back(h.getName());
341  else if (omp_get_thread_num() == 0)
342  app_log() << " warning: potential named " << h.getName()
343  << " has not been classified according to its quantum domain (q,c,qq,qc,cc)\n estimators depending "
344  "on this classification may not function properly"
345  << std::endl;
346  }
347 
348 
349  //make trace quantities available
350  request.contribute_scalar("id", true); //default trace quantity
351  request.contribute_scalar("parent_id", true); //default trace quantity
352  request.contribute_scalar("step", true); //default trace quantity
353  request.contribute_scalar("generation", true); //default trace quantity
354  request.contribute_scalar("age", true); //default trace quantity
355  request.contribute_scalar("multiplicity", true); //default trace quantity
356  request.contribute_scalar("weight", true); //default trace quantity
357  request.contribute_array("position");
358  for (int i = 0; i < H.size(); ++i)
359  H[i]->contributeTraceQuantities();
360  for (int i = 0; i < auxH.size(); ++i)
361  auxH[i]->contributeTraceQuantities();
362 
363 
364  //note availability of combined quantities
365  request.contribute_combined("LocalEnergy", Eloc, true);
366  request.contribute_combined("LocalPotential", Vloc, true, true);
367  if (Vq.size() > 0)
368  request.contribute_combined("Vq", Vq, true, true);
369  if (Vc.size() > 0)
370  request.contribute_combined("Vc", Vc, true, true);
371  if (Vqq.size() > 0)
372  request.contribute_combined("Vqq", Vqq, true, true);
373  if (Vqc.size() > 0)
374  request.contribute_combined("Vqc", Vqc, true, true);
375  if (Vcc.size() > 0)
376  request.contribute_combined("Vcc", Vcc, true, true);
377 
378 
379  //collect trace requests
380  std::vector<TraceRequest*> requests;
381  // Hamiltonian request (id, step, weight, positions)
382  requests.push_back(&request);
383  // requests from Hamiltonian components
384  for (int i = 0; i < H.size(); ++i)
385  requests.push_back(&H[i]->getRequest());
386  // requests from other observables
387  for (int i = 0; i < auxH.size(); ++i)
388  requests.push_back(&auxH[i]->getRequest());
389 
390  //collect trace quantity availability/requests from contributors/requestors
391  for (int i = 0; i < requests.size(); ++i)
392  tm.request.incorporate(*requests[i]);
393 
394  //balance requests with availability, mark quantities as streaming/writing
395  tm.request.determine_stream_write();
396 
397  //relay updated streaming information to all contributors/requestors
398  for (int i = 0; i < requests.size(); ++i)
399  tm.request.relay_stream_info(*requests[i]);
400 
401  //set streaming/writing traces in general
402  tm.update_status();
403 
404  // setup traces, if any quantities should be streaming
405 
406  // tracing
407  bool tracing = request.streaming();
408  if (tracing != tm.streaming_traces)
409  APP_ABORT("QMCHamiltonian::initialize_traces trace request failed to initialize properly");
410  if (!tracing)
411  {
412  // Empty. Do not log if nothing will be done
413 
414  if (trace_log)
415  app_log() << " no traces streaming" << std::endl;
416  }
417  else
418  {
419  if (trace_log)
420  app_log() << " traces streaming" << std::endl;
421  //checkout trace quantities
422  //(requested sources checkout arrays to place samples in for streaming)
423  // checkout walker trace quantities
426  {
427  id_sample = tm.checkout_int<1>("id");
428  pid_sample = tm.checkout_int<1>("parent_id");
429  step_sample = tm.checkout_int<1>("step");
430  gen_sample = tm.checkout_int<1>("generation");
431  age_sample = tm.checkout_int<1>("age");
432  mult_sample = tm.checkout_int<1>("multiplicity");
433  weight_sample = tm.checkout_real<1>("weight");
434  }
435  if (streaming_position)
436  position_sample = tm.checkout_real<2>("position", P, DIM);
437  // checkout observable trace quantities
438  for (int i = 0; i < H.size(); ++i)
439  {
440  if (trace_log)
441  app_log() << " OperatorBase::checkoutTraceQuantities " << H[i]->getName() << std::endl;
442  H[i]->checkoutTraceQuantities(tm);
443  }
444  for (int i = 0; i < auxH.size(); ++i)
445  {
446  if (trace_log)
447  app_log() << " OperatorBase::checkoutTraceQuantities " << auxH[i]->getName() << std::endl;
448  auxH[i]->checkoutTraceQuantities(tm);
449  }
450  //setup combined traces that depend on H information
451  // LocalEnergy, LocalPotential, Vq, Vc, Vqq, Vqc, Vcc
452  if (Vloc.size() > 0 && request.streaming("LocalPotential"))
453  tm.make_combined_trace("LocalPotential", Vloc);
454  if (Eloc.size() > 0 && request.streaming("LocalEnergy"))
455  tm.make_combined_trace("LocalEnergy", Eloc);
456  if (Vq.size() > 0 && request.streaming("Vq"))
457  tm.make_combined_trace("Vq", Eloc);
458  if (Vc.size() > 0 && request.streaming("Vc"))
459  tm.make_combined_trace("Vc", Eloc);
460  if (Vqq.size() > 0 && request.streaming("Vqq"))
461  tm.make_combined_trace("Vqq", Eloc);
462  if (Vqc.size() > 0 && request.streaming("Vqc"))
463  tm.make_combined_trace("Vqc", Eloc);
464  if (Vcc.size() > 0 && request.streaming("Vcc"))
465  tm.make_combined_trace("Vcc", Eloc);
466 
467  //all trace samples have been created ( streaming instances)
468  // mark the ones that will be writing also
469  tm.screen_writes();
470 
471  //observables that depend on traces check them out
472  if (trace_log)
473  app_log() << "\n Hamiltonian is fulfilling trace requests from observables" << std::endl;
474  for (int i = 0; i < auxH.size(); ++i)
475  {
476  if (trace_log)
477  app_log() << " OperatorBase::getRequiredTraces " << auxH[i]->getName() << std::endl;
478  auxH[i]->getRequiredTraces(tm);
479  }
480  //report
481 
482  //write traces status to the log
483  if (trace_log)
484  tm.user_report();
485 
486  first_init = false;
487  }
488 }
Array< TraceInt, 1 > * gen_sample
Array< TraceInt, 1 > * step_sample
Array< TraceInt, 1 > * pid_sample
bool streaming_array(const std::string &name)
Definition: TraceManager.h:249
std::ostream & app_log()
Definition: OutputManager.h:65
bool streaming(const std::string &name)
Definition: TraceManager.h:256
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
Array< TraceReal, 1 > * weight_sample
TraceRequest request
traces variables
void contribute_scalar(const std::string &name, bool default_quantity=false)
Definition: TraceManager.h:189
Array< TraceInt, 1 > * age_sample
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
Array< TraceReal, 2 > * position_sample
void contribute_combined(const std::string &name, std::vector< std::string > &deps, bool scalar=false, bool array=false, bool default_quantity=false)
Definition: TraceManager.h:207
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
void contribute_array(const std::string &name, bool default_quantity=false)
Definition: TraceManager.h:198
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians
Array< TraceInt, 1 > * mult_sample
Array< TraceInt, 1 > * id_sample

◆ makeClone()

std::unique_ptr< QMCHamiltonian > makeClone ( ParticleSet qp,
TrialWaveFunction psi 
) const

return a clone

Definition at line 1066 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, PooledData< T >::clear(), ParticleSet::Collectables, QMCHamiltonian::H, QMCHamiltonian::myIndex, QMCHamiltonian::myName, QMCHamiltonian::numCollectables, and PooledData< T >::resize().

Referenced by CostFunctionCrowdData::CostFunctionCrowdData(), MCPopulation::createWalkers(), Crowd::Crowd(), CloneManager::makeClones(), WaveFunctionTester::runCloneTest(), WaveFunctionTester::runDerivCloneTest(), and MCPopulation::spawnWalker().

1067 {
1068  auto myclone = std::make_unique<QMCHamiltonian>(myName);
1069  for (int i = 0; i < H.size(); ++i)
1070  H[i]->add2Hamiltonian(qp, psi, *myclone);
1071  for (int i = 0; i < auxH.size(); ++i)
1072  auxH[i]->add2Hamiltonian(qp, psi, *myclone);
1073  //sync indices
1074  myclone->resetObservables(myIndex, numCollectables);
1075  //Hamiltonian needs to make sure qp.Collectables are the same as defined by the original Hamiltonian
1076  if (numCollectables)
1077  {
1078  qp.Collectables.clear();
1079  qp.Collectables.resize(numCollectables);
1080  }
1081  //Assume tau is correct for the Kinetic energy operator and assign to the rest of the clones
1082  //Return_t tau = H[0]->Tau;
1083  //myclone->setTau(tau);
1084  return myclone;
1085 }
int myIndex
starting index
int numCollectables
starting index
const std::string myName
getName is in the way
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ makeNonLocalMoves()

int makeNonLocalMoves ( ParticleSet P)

make non local moves

Parameters
Pparticle set
Returns
the number of accepted moves

Definition at line 1011 of file QMCHamiltonian.cpp.

References NonLocalECPotential::makeNonLocalMovesPbyP(), and QMCHamiltonian::nlpp_ptr.

Referenced by DMCUpdateAllWithRejection::advanceWalker(), SODMCUpdatePbyPWithRejectionFast::advanceWalker(), DMCUpdatePbyPWithRejectionFast::advanceWalker(), and DMCUpdatePbyPL2::advanceWalker().

1012 {
1013  if (nlpp_ptr == nullptr)
1014  return 0;
1015  else
1016  return nlpp_ptr->makeNonLocalMovesPbyP(P);
1017 }
NonLocalECPotential * nlpp_ptr
pointer to NonLocalECP
int makeNonLocalMovesPbyP(ParticleSet &P)
make non local moves with particle-by-particle moves

◆ mw_evaluate()

std::vector< QMCHamiltonian::FullPrecRealType > mw_evaluate ( const RefVectorWithLeader< QMCHamiltonian > &  ham_list,
const RefVectorWithLeader< TrialWaveFunction > &  wf_list,
const RefVectorWithLeader< ParticleSet > &  p_list 
)
static

batched version of evaluate for LocalEnergy

Encapsulation is ignored for ham_list hamiltonians method uses its status as QMCHamiltonian to break encapsulation. ParticleSet is also updated. Bugs could easily be created by accessing this scope. This should be set to static and fixed.

Definition at line 596 of file QMCHamiltonian.cpp.

References QMCHamiltonian::extract_HC_list(), RefVectorWithLeader< T >::getLeader(), QMCHamiltonian::getLocalEnergy(), qmcplusplus::ham, QMCHamiltonian::LocalEnergy, QMCHamiltonian::updateComponent(), and QMCHamiltonian::updateKinetic().

Referenced by QMCCostFunctionBatched::checkConfigurations(), QMCCostFunctionBatched::correlatedSampling(), Hdispatcher::flex_evaluate(), and qmcplusplus::TEST_CASE().

600 {
601  auto& ham_leader = ham_list.getLeader();
602  ScopedTimer local_timer(ham_leader.ham_timer_);
603  for (QMCHamiltonian& ham : ham_list)
604  ham.LocalEnergy = 0.0;
605 
606  const int num_ham_operators = ham_leader.H.size();
607 
608  // It is an invariant of the class that H[0] be the kinetic energy operator.
609  // This is enforced by HamiltonianFactory's constuctor and not QMCHamiltonians
610  int kinetic_index = 0;
611  {
612  ScopedTimer h_timer(ham_leader.my_timers_[kinetic_index]);
613  const auto HC_list(extract_HC_list(ham_list, kinetic_index));
614  if (ham_leader.mw_res_handle_.getResource().kinetic_listeners_.size() > 0)
615  ham_leader.H[kinetic_index]
616  ->mw_evaluatePerParticle(HC_list, wf_list, p_list, ham_leader.mw_res_handle_.getResource().kinetic_listeners_,
617  ham_leader.mw_res_handle_.getResource().ion_kinetic_listeners_);
618  else
619  ham_leader.H[kinetic_index]->mw_evaluate(HC_list, wf_list, p_list);
620  for (int iw = 0; iw < ham_list.size(); iw++)
621  updateComponent(HC_list[iw], ham_list[iw], p_list[iw]);
622  }
623 
624  for (int i_ham_op = 1; i_ham_op < num_ham_operators; ++i_ham_op)
625  {
626  ScopedTimer h_timer(ham_leader.my_timers_[i_ham_op]);
627  const auto HC_list(extract_HC_list(ham_list, i_ham_op));
628 
629  if (ham_leader.mw_res_handle_.getResource().potential_listeners_.size() > 0)
630  ham_leader.H[i_ham_op]->mw_evaluatePerParticle(HC_list, wf_list, p_list,
631  ham_leader.mw_res_handle_.getResource().potential_listeners_,
632  ham_leader.mw_res_handle_.getResource().ion_potential_listeners_);
633  else
634  ham_leader.H[i_ham_op]->mw_evaluate(HC_list, wf_list, p_list);
635  for (int iw = 0; iw < ham_list.size(); iw++)
636  updateComponent(HC_list[iw], ham_list[iw], p_list[iw]);
637  }
638 
639  for (int iw = 0; iw < ham_list.size(); iw++)
640  updateKinetic(ham_list[iw], p_list[iw]);
641 
642  std::vector<FullPrecRealType> local_energies(ham_list.size(), 0.0);
643  for (int iw = 0; iw < ham_list.size(); ++iw)
644  local_energies[iw] = ham_list[iw].getLocalEnergy();
645 
646  return local_energies;
647 }
static void updateComponent(OperatorBase &op, QMCHamiltonian &ham, ParticleSet &pset)
accumulate local energy and update Observables and PropertyList
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
FullPrecRealType LocalEnergy
Current Local Energy.
FullPrecRealType getLocalEnergy()
static void updateKinetic(QMCHamiltonian &ham, ParticleSet &pset)
extract kinetic and potential energies.
static RefVectorWithLeader< OperatorBase > extract_HC_list(const RefVectorWithLeader< QMCHamiltonian > &ham_list, int id)
QMCHamiltonian(const std::string &aname="psi0")
constructor

◆ mw_evaluateValueAndDerivatives()

std::vector< QMCHamiltonian::FullPrecRealType > mw_evaluateValueAndDerivatives ( const RefVectorWithLeader< QMCHamiltonian > &  ham_list,
const RefVectorWithLeader< TrialWaveFunction > &  wf_list,
const RefVectorWithLeader< ParticleSet > &  p_list,
const opt_variables_type optvars,
RecordArray< ValueType > &  dlogpsi,
RecordArray< ValueType > &  dhpsioverpsi 
)
static

Definition at line 673 of file QMCHamiltonian.cpp.

References QMCHamiltonian::extract_HC_list(), RefVectorWithLeader< T >::getLeader(), QMCHamiltonian::getLocalEnergy(), QMCHamiltonian::LocalEnergy, QMCHamiltonian::updateComponent(), and QMCHamiltonian::updateKinetic().

Referenced by QMCCostFunctionBatched::checkConfigurations(), QMCCostFunctionBatched::correlatedSampling(), and qmcplusplus::test_hcpBe_rotation().

680 {
681  std::vector<FullPrecRealType> local_energies(ham_list.size(), 0.0);
682  for (int iw = 0; iw < ham_list.size(); iw++)
683  ham_list[iw].LocalEnergy = 0.0;
684 
685  if (ham_list.size() > 0)
686  {
687  auto& ham_leader = ham_list.getLeader();
688  const int num_ham_operators = ham_leader.H.size();
689  for (int i_ham_op = 0; i_ham_op < num_ham_operators; ++i_ham_op)
690  {
691  ScopedTimer local_timer(ham_leader.my_timers_[i_ham_op]);
692  const auto HC_list(extract_HC_list(ham_list, i_ham_op));
693 
694  ham_leader.H[i_ham_op]->mw_evaluateWithParameterDerivatives(HC_list, p_list, optvars, dlogpsi, dhpsioverpsi);
695 
696  for (int iw = 0; iw < ham_list.size(); iw++)
697  updateComponent(HC_list[iw], ham_list[iw], p_list[iw]);
698  }
699 
700  for (int iw = 0; iw < ham_list.size(); iw++)
701  updateKinetic(ham_list[iw], p_list[iw]);
702 
703  for (int iw = 0; iw < ham_list.size(); ++iw)
704  local_energies[iw] = ham_list[iw].getLocalEnergy();
705  }
706 
707  return local_energies;
708 }
static void updateComponent(OperatorBase &op, QMCHamiltonian &ham, ParticleSet &pset)
accumulate local energy and update Observables and PropertyList
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
FullPrecRealType LocalEnergy
Current Local Energy.
FullPrecRealType getLocalEnergy()
static void updateKinetic(QMCHamiltonian &ham, ParticleSet &pset)
extract kinetic and potential energies.
static RefVectorWithLeader< OperatorBase > extract_HC_list(const RefVectorWithLeader< QMCHamiltonian > &ham_list, int id)

◆ mw_evaluateWithToperator()

std::vector< QMCHamiltonian::FullPrecRealType > mw_evaluateWithToperator ( const RefVectorWithLeader< QMCHamiltonian > &  ham_list,
const RefVectorWithLeader< TrialWaveFunction > &  wf_list,
const RefVectorWithLeader< ParticleSet > &  p_list 
)
static

batched version of evaluate Local energy with Toperators updated.

Definition at line 813 of file QMCHamiltonian.cpp.

References QMCHamiltonian::extract_HC_list(), QMCHamiltonian::getLocalEnergy(), qmcplusplus::ham, QMCHamiltonian::LocalEnergy, QMCHamiltonian::updateComponent(), and QMCHamiltonian::updateKinetic().

Referenced by Hdispatcher::flex_evaluateWithToperator().

817 {
818  for (QMCHamiltonian& ham : ham_list)
819  ham.LocalEnergy = 0.0;
820 
821  auto& ham_leader = ham_list.getLeader();
822  const int num_ham_operators = ham_leader.H.size();
823 
824  const int kinetic_index = 0;
825  {
826  ScopedTimer h_timer(ham_leader.my_timers_[kinetic_index]);
827  const auto HC_list(extract_HC_list(ham_list, kinetic_index));
828  if (ham_leader.mw_res_handle_.getResource().kinetic_listeners_.size() > 0)
829  ham_leader.H[kinetic_index]
830  ->mw_evaluatePerParticleWithToperator(HC_list, wf_list, p_list,
831  ham_leader.mw_res_handle_.getResource().kinetic_listeners_,
832  ham_leader.mw_res_handle_.getResource().ion_kinetic_listeners_);
833  else
834  ham_leader.H[kinetic_index]->mw_evaluateWithToperator(HC_list, wf_list, p_list);
835  for (int iw = 0; iw < ham_list.size(); iw++)
836  updateComponent(HC_list[iw], ham_list[iw], p_list[iw]);
837  }
838 
839  for (int i_ham_op = 1; i_ham_op < num_ham_operators; ++i_ham_op)
840  {
841  ScopedTimer local_timer(ham_leader.my_timers_[i_ham_op]);
842  const auto HC_list(extract_HC_list(ham_list, i_ham_op));
843  if (ham_leader.mw_res_handle_.getResource().potential_listeners_.size() > 0)
844  ham_leader.H[i_ham_op]
845  ->mw_evaluatePerParticleWithToperator(HC_list, wf_list, p_list,
846  ham_leader.mw_res_handle_.getResource().potential_listeners_,
847  ham_leader.mw_res_handle_.getResource().ion_potential_listeners_);
848  else
849  ham_leader.H[i_ham_op]->mw_evaluateWithToperator(HC_list, wf_list, p_list);
850  for (int iw = 0; iw < ham_list.size(); ++iw)
851  updateComponent(HC_list[iw], ham_list[iw], p_list[iw]);
852  }
853 
854  for (int iw = 0; iw < ham_list.size(); iw++)
855  updateKinetic(ham_list[iw], p_list[iw]);
856 
857  std::vector<FullPrecRealType> local_energies(ham_list.size());
858  for (int iw = 0; iw < ham_list.size(); ++iw)
859  local_energies[iw] = ham_list[iw].getLocalEnergy();
860 
861  return local_energies;
862 }
static void updateComponent(OperatorBase &op, QMCHamiltonian &ham, ParticleSet &pset)
accumulate local energy and update Observables and PropertyList
ScopeGuard< NewTimer > ScopedTimer
Definition: NewTimer.h:257
FullPrecRealType LocalEnergy
Current Local Energy.
FullPrecRealType getLocalEnergy()
static void updateKinetic(QMCHamiltonian &ham, ParticleSet &pset)
extract kinetic and potential energies.
static RefVectorWithLeader< OperatorBase > extract_HC_list(const RefVectorWithLeader< QMCHamiltonian > &ham_list, int id)
QMCHamiltonian(const std::string &aname="psi0")
constructor

◆ mw_makeNonLocalMoves()

std::vector< int > mw_makeNonLocalMoves ( const RefVectorWithLeader< QMCHamiltonian > &  ham_list,
const RefVectorWithLeader< TrialWaveFunction > &  wf_list,
const RefVectorWithLeader< ParticleSet > &  p_list 
)
static

Definition at line 1020 of file QMCHamiltonian.cpp.

References RefVectorWithLeader< T >::getLeader(), NonLocalECPotential::makeNonLocalMovesPbyP(), and QMCHamiltonian::nlpp_ptr.

Referenced by Hdispatcher::flex_makeNonLocalMoves().

1023 {
1024  auto& ham_leader = ham_list.getLeader();
1025 
1026  std::vector<int> num_accepts(ham_list.size(), 0);
1027  if (ham_list.getLeader().nlpp_ptr)
1028  {
1029  for (int iw = 0; iw < ham_list.size(); ++iw)
1030  num_accepts[iw] = ham_list[iw].nlpp_ptr->makeNonLocalMovesPbyP(p_list[iw]);
1031  }
1032  return num_accepts;
1033 }
NonLocalECPotential * nlpp_ptr
pointer to NonLocalECP
int makeNonLocalMovesPbyP(ParticleSet &P)
make non local moves with particle-by-particle moves

◆ rejectedMove()

void rejectedMove ( ParticleSet P,
Walker_t ThisWalker 
)

Looks like a hack see DMCBatched.cpp and DMC.cpp weight is used like temporary flag from DMC.

Definition at line 780 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, QMCHamiltonian::collect_walker_traces(), and ParticleSet::current_step.

Referenced by SODMCUpdatePbyPWithRejectionFast::advanceWalker(), DMCUpdatePbyPWithRejectionFast::advanceWalker(), DMCUpdatePbyPL2::advanceWalker(), DMCUpdateAllWithKill::advanceWalker(), RMCUpdatePbyPWithDrift::advanceWalkersRMC(), RMCUpdateAllWithDrift::advanceWalkersRMC(), RMCUpdatePbyPWithDrift::advanceWalkersVMC(), and RMCUpdateAllWithDrift::advanceWalkersVMC().

781 {
782  // weight should be 0 from DMC
783  // since other traced properties will be invalid
784  // (they will be from the walker moved before this one)
785 #if !defined(REMOVE_TRACEMANAGER)
786  collect_walker_traces(ThisWalker, P.current_step);
787 #endif
788  // ThisWalker.rejectedMove();
789  for (int i = 0; i < auxH.size(); ++i)
790  {
791  auxH[i]->setHistories(ThisWalker);
792  RealType sink = auxH[i]->rejectedMove(P);
793  }
794 }
void collect_walker_traces(Walker_t &walker, int step)
collect walker trace data
QMCTraits::RealType RealType
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ releaseResource()

void releaseResource ( ResourceCollection collection,
const RefVectorWithLeader< QMCHamiltonian > &  ham_list 
)
static

release external resource Note: use RAII ResourceCollectionLock whenever possible

Definition at line 1054 of file QMCHamiltonian.cpp.

References QMCHamiltonian::extract_HC_list(), RefVectorWithLeader< T >::getLeader(), and ResourceCollection::takebackResource().

1056 {
1057  auto& ham_leader = ham_list.getLeader();
1058  collection.takebackResource(ham_leader.mw_res_handle_);
1059  for (int i_ham_op = 0; i_ham_op < ham_leader.H.size(); ++i_ham_op)
1060  {
1061  const auto HC_list(extract_HC_list(ham_list, i_ham_op));
1062  ham_leader.H[i_ham_op]->releaseResource(collection, HC_list);
1063  }
1064 }
static RefVectorWithLeader< OperatorBase > extract_HC_list(const RefVectorWithLeader< QMCHamiltonian > &ham_list, int id)

◆ reportToListeners()

void reportToListeners ( )
private

◆ resetObservables()

void resetObservables ( int  start,
int  ncollects 
)
private

reset Observables and counters

Parameters
startstarting index within PropertyList
ncollectsnumber of collectables

Definition at line 245 of file QMCHamiltonian.cpp.

References QMCHamiltonian::addObservables(), APP_ABORT, QMCHamiltonian::auxH, QMCHamiltonian::H, QMCHamiltonian::myIndex, QMCHamiltonian::numCollectables, QMCHamiltonian::Observables, PooledData< T >::rewind(), and PooledData< T >::size().

246 {
247  Observables.clear();
248  BufferType collectables;
249  collectables.rewind();
250  for (int i = 0; i < H.size(); ++i)
251  H[i]->addObservables(Observables, collectables);
252  for (int i = 0; i < auxH.size(); ++i)
253  auxH[i]->addObservables(Observables, collectables);
254  if (collectables.size() != ncollects)
255  {
256  APP_ABORT(" QMCHamiltonian::resetObservables numCollectables != ncollects");
257  }
258  myIndex = start;
259  numCollectables = ncollects;
260 }
int addObservables(ParticleSet &P)
add each term to the PropertyList for averages
OperatorBase::BufferType BufferType
void rewind(size_type cur=0)
set the Current to a cursor
Definition: PooledData.h:56
int myIndex
starting index
int numCollectables
starting index
#define APP_ABORT(msg)
Widely used but deprecated fatal error macros from legacy code.
Definition: AppAbort.h:27
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
PropertySetType Observables
data
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ resetTargetParticleSet()

void resetTargetParticleSet ( ParticleSet P)

Definition at line 978 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, and QMCHamiltonian::H.

979 {
980  for (int i = 0; i < H.size(); i++)
981  H[i]->resetTargetParticleSet(P);
982  for (int i = 0; i < auxH.size(); i++)
984 }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians
void resetTargetParticleSet(ParticleSet &P)

◆ setNonLocalMoves() [1/2]

void setNonLocalMoves ( xmlNodePtr  cur)

set non local moves options

Parameters
curthe xml input

Definition at line 996 of file QMCHamiltonian.cpp.

References QMCHamiltonian::nlpp_ptr, and NonLocalECPotential::setNonLocalMoves().

Referenced by QMCUpdateBase::put(), and DMCBatched::setNonLocalMoveHandler().

997 {
998  if (nlpp_ptr != nullptr)
1000 }
NonLocalECPotential * nlpp_ptr
pointer to NonLocalECP
void setNonLocalMoves(xmlNodePtr cur)
set non local moves options

◆ setNonLocalMoves() [2/2]

void setNonLocalMoves ( const std::string &  non_local_move_option,
const double  tau,
const double  alpha,
const double  gamma 
)

Definition at line 1002 of file QMCHamiltonian.cpp.

References QMCHamiltonian::nlpp_ptr, and NonLocalECPotential::setNonLocalMoves().

1006 {
1007  if (nlpp_ptr != nullptr)
1008  nlpp_ptr->setNonLocalMoves(non_local_move_option, tau, alpha, gamma);
1009 }
NonLocalECPotential * nlpp_ptr
pointer to NonLocalECP
void setNonLocalMoves(xmlNodePtr cur)
set non local moves options

◆ setPrimary()

void setPrimary ( bool  primary)
inline

set PRIMARY bit of all the components

Definition at line 224 of file QMCHamiltonian.h.

References QMCHamiltonian::H, and OperatorBase::PRIMARY.

Referenced by QMCDriverFactory::createQMCDriver(), WaveFunctionTester::runDerivNLPPTest(), and WaveFunctionTester::runDerivTest().

225  {
226  for (int i = 0; i < H.size(); i++)
227  H[i]->getUpdateMode().set(OperatorBase::PRIMARY, primary);
228  }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ setProperty()

void setProperty ( IT  first)
inline

Definition at line 203 of file QMCHamiltonian.h.

References copy(), QMCHamiltonian::myIndex, and QMCHamiltonian::Observables.

204  {
205  // LocalEnergy=first[WP::LOCALENERGY];
206  // KineticEnergy=LocalEnergy-first[LOCALPOTENTIAL];
207  copy(first + myIndex, first + myIndex + Observables.size(), Observables.begin());
208  }
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639
int myIndex
starting index
PropertySetType Observables
data

◆ setRandomGenerator()

void setRandomGenerator ( RandomBase< FullPrecRealType > *  rng)

Definition at line 986 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, QMCHamiltonian::H, QMCHamiltonian::nlpp_ptr, and NonLocalECPotential::setRandomGenerator().

Referenced by QMCCostFunctionBatched::checkConfigurations(), QMCCostFunctionBatched::correlatedSampling(), WaveFunctionTester::run(), WaveFunctionTester::runDerivCloneTest(), Crowd::setRNGForHamiltonian(), qmcplusplus::TEST_CASE(), and qmcplusplus::test_hcpBe_rotation().

987 {
988  for (int i = 0; i < H.size(); i++)
989  H[i]->setRandomGenerator(rng);
990  for (int i = 0; i < auxH.size(); i++)
991  auxH[i]->setRandomGenerator(rng);
992  if (nlpp_ptr)
994 }
NonLocalECPotential * nlpp_ptr
pointer to NonLocalECP
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
void setRandomGenerator(RandomBase< FullPrecRealType > *rng) override
set the internal RNG pointer as the given pointer
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians
void setRandomGenerator(RandomBase< FullPrecRealType > *rng)

◆ size()

int size ( void  ) const
inline

return the number of Hamiltonians

Definition at line 86 of file QMCHamiltonian.h.

References QMCHamiltonian::H.

Referenced by qmcplusplus::TEST_CASE(), and QMCMain::validateXML().

86 { return H.size(); }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ total_size()

int total_size ( ) const
inline

return the total number of Hamiltonians (physical + aux)

Definition at line 89 of file QMCHamiltonian.h.

References QMCHamiltonian::auxH, and QMCHamiltonian::H.

Referenced by qmcplusplus::TEST_CASE().

89 { return H.size() + auxH.size(); }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians

◆ updateComponent()

void updateComponent ( OperatorBase op,
QMCHamiltonian ham,
ParticleSet pset 
)
static

accumulate local energy and update Observables and PropertyList

Definition at line 573 of file QMCHamiltonian.cpp.

References OperatorBase::getName(), OperatorBase::getValue(), qmcplusplus::ham, qmcplusplus::isnan(), QMCHamiltonian::LocalEnergy, QMCHamiltonian::myIndex, QMCHamiltonian::Observables, qmcplusplus::pset, OperatorBase::setObservables(), and OperatorBase::setParticlePropertyList().

Referenced by QMCHamiltonian::evaluate(), QMCHamiltonian::evaluateDeterministic(), QMCHamiltonian::evaluateWithToperator(), QMCHamiltonian::mw_evaluate(), QMCHamiltonian::mw_evaluateValueAndDerivatives(), and QMCHamiltonian::mw_evaluateWithToperator().

574 {
575  // It's much better to be able to see where this is coming from. It is caught just fine.
576  if (qmcplusplus::isnan(op.getValue()))
577  {
578  std::ostringstream msg;
579  msg << "QMCHamiltonian::updateComponent component " << op.getName() << " returns NaN." << std::endl;
580  pset.print(msg);
581  throw std::runtime_error(msg.str());
582  }
583  // The following is a ridiculous breach of encapsulation.
584  ham.LocalEnergy += op.getValue();
585  op.setObservables(ham.Observables);
586  op.setParticlePropertyList(pset.PropertyList, ham.myIndex);
587 }
int myIndex
starting index
FullPrecRealType LocalEnergy
Current Local Energy.
PropertySetType Observables
data
bool isnan(float a)
return true if the value is NaN.
Definition: math.cpp:18

◆ updateKinetic()

void updateKinetic ( QMCHamiltonian ham,
ParticleSet pset 
)
static

extract kinetic and potential energies.

Definition at line 589 of file QMCHamiltonian.cpp.

References QMCHamiltonian::H, qmcplusplus::ham, QMCHamiltonian::KineticEnergy, QMCHamiltonian::LocalEnergy, and qmcplusplus::pset.

Referenced by QMCHamiltonian::evaluate(), QMCHamiltonian::evaluateDeterministic(), QMCHamiltonian::evaluateWithToperator(), QMCHamiltonian::mw_evaluate(), QMCHamiltonian::mw_evaluateValueAndDerivatives(), and QMCHamiltonian::mw_evaluateWithToperator().

590 {
591  ham.KineticEnergy = ham.H[0]->getValue();
592  pset.PropertyList[WP::LOCALENERGY] = ham.LocalEnergy;
593  pset.PropertyList[WP::LOCALPOTENTIAL] = ham.LocalEnergy - ham.KineticEnergy;
594 }
FullPrecRealType KineticEnergy
Current Kinetic Energy.
FullPrecRealType LocalEnergy
Current Local Energy.
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians

◆ updateSource()

void updateSource ( ParticleSet s)

remove a named Hamiltonian from the list

Parameters
anamethe name of the Hamiltonian
Returns
true, if the request hamiltonian exists and is removed.

Definition at line 192 of file QMCHamiltonian.cpp.

References QMCHamiltonian::auxH, QMCHamiltonian::H, and qmcplusplus::Units::time::s.

Referenced by QMCMain::executeCMCSection().

193 {
194  for (int i = 0; i < H.size(); ++i)
195  H[i]->updateSource(s);
196  for (int i = 0; i < auxH.size(); ++i)
197  auxH[i]->updateSource(s);
198 }
std::vector< std::unique_ptr< OperatorBase > > H
vector of Hamiltonians
std::vector< std::unique_ptr< OperatorBase > > auxH
vector of Hamiltonians
void updateSource(ParticleSet &s)
remove a named Hamiltonian from the list

Friends And Related Function Documentation

◆ HamiltonianFactory

friend class HamiltonianFactory
friend

Definition at line 51 of file QMCHamiltonian.h.

Member Data Documentation

◆ age_sample

◆ auxH

◆ available_quantities_

constexpr std::array<std::string_view, 8> available_quantities_
staticprivate
Initial value:
{"weight", "LocalEnergy", "LocalPotential",
"Vq", "Vc", "Vqq",
"Vqc", "Vcc"}

Definition at line 431 of file QMCHamiltonian.h.

◆ eval_ion_derivs_fast_timer_

NewTimer& eval_ion_derivs_fast_timer_
private

Total timer for H ion deriv evaluation;.

Definition at line 460 of file QMCHamiltonian.h.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

◆ eval_vals_derivs_timer_

NewTimer& eval_vals_derivs_timer_
private

Total timer for H evaluation.

Definition at line 458 of file QMCHamiltonian.h.

Referenced by QMCHamiltonian::evaluateValueAndDerivatives().

◆ gen_sample

Array<TraceInt, 1>* gen_sample
private

◆ H

◆ ham_timer_

NewTimer& ham_timer_
private

Total timer for H evaluation.

Definition at line 456 of file QMCHamiltonian.h.

Referenced by QMCHamiltonian::evaluate(), QMCHamiltonian::evaluateDeterministic(), and QMCHamiltonian::evaluateWithToperator().

◆ id_sample

Array<TraceInt, 1>* id_sample
private

◆ KineticEnergy

◆ l2_ptr

◆ LocalEnergy

◆ mult_sample

◆ mw_res_handle_

◆ my_timers_

std::vector<std::reference_wrapper<NewTimer> > my_timers_
private

◆ myIndex

◆ myName

const std::string myName
private

getName is in the way

Definition at line 446 of file QMCHamiltonian.h.

Referenced by QMCHamiltonian::getName(), and QMCHamiltonian::makeClone().

◆ NewLocalEnergy

FullPrecRealType NewLocalEnergy
private

Current Local Energy for the proposed move.

Definition at line 444 of file QMCHamiltonian.h.

◆ nlpp_ptr

◆ numCollectables

int numCollectables
private

◆ Observables

◆ operator_types

std::map<std::string, std::string> operator_types
private

types of component operators

Definition at line 464 of file QMCHamiltonian.h.

Referenced by QMCHamiltonian::addOperatorType(), and QMCHamiltonian::getOperatorType().

◆ pid_sample

Array<TraceInt, 1>* pid_sample
private

◆ position_sample

◆ request

◆ step_sample

Array<TraceInt, 1>* step_sample
private

◆ streaming_position

bool streaming_position
private

◆ weight_sample


The documentation for this class was generated from the following files: