QMCPACK
LatticeGaussianProduct Class Reference

A composite Orbital. More...

+ Inheritance diagram for LatticeGaussianProduct:
+ Collaboration diagram for LatticeGaussianProduct:

Public Member Functions

 LatticeGaussianProduct (ParticleSet &centers, ParticleSet &ptcls)
 
 ~LatticeGaussianProduct () override
 
std::string getClassName () const override
 return class name More...
 
LogValue evaluateLog (const ParticleSet &P, ParticleSet::ParticleGradient &G, ParticleSet::ParticleLaplacian &L) override
 
PsiValue ratio (ParticleSet &P, int iat) override
 evaluate the ratio $exp(U(iat)-U_0(iat))$ More...
 
void acceptMove (ParticleSet &P, int iat, bool safe_to_delay=false) override
 a move for iat-th particle is accepted. More...
 
void restore (int iat) override
 If a move for iat-th particle is rejected, restore to the content. More...
 
void registerData (ParticleSet &P, WFBufferType &buf) override
 equivalent to evalaute with additional data management More...
 
LogValue updateBuffer (ParticleSet &P, WFBufferType &buf, bool fromscratch) override
 For particle-by-particle move. More...
 
void copyFromBuffer (ParticleSet &P, WFBufferType &buf) override
 copy the current data from a buffer More...
 
GradType evalGrad (ParticleSet &P, int iat) override
 return the current gradient for the iat-th particle More...
 
PsiValue ratioGrad (ParticleSet &P, int iat, GradType &grad_iat) override
 evaluate the ratio of the new to old WaveFunctionComponent value and the new gradient More...
 
std::unique_ptr< WaveFunctionComponentmakeClone (ParticleSet &tqp) const override
 make clone More...
 
void evaluateDerivatives (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, Vector< ValueType > &dhpsioverpsi) override
 
void evaluateLogAndStore (const ParticleSet &P, ParticleSet::ParticleGradient &dG, ParticleSet::ParticleLaplacian &dL)
 
- Public Member Functions inherited from WaveFunctionComponent
const LogValueget_log_value () const
 
 WaveFunctionComponent (const std::string &obj_name="")
 default constructor More...
 
virtual ~WaveFunctionComponent ()
 default destructor More...
 
virtual void checkSanity () const
 Validate the internal consistency of the object. More...
 
const std::string & getName () const
 return object name More...
 
PsiValue getValue () const
 assembles the full value More...
 
virtual bool isFermionic () const
 true, if this component is fermionic More...
 
virtual bool isMultiDet () const
 true, if this component is multi-determinant More...
 
virtual void checkOutVariables (const opt_variables_type &active)
 check out variational optimizable variables More...
 
virtual void registerTWFFastDerivWrapper (const ParticleSet &P, TWFFastDerivWrapper &twf) const
 Register the component with the TWFFastDerivWrapper wrapper. More...
 
virtual void mw_evaluateLog (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, const RefVector< ParticleSet::ParticleGradient > &G_list, const RefVector< ParticleSet::ParticleLaplacian > &L_list) const
 evaluate from scratch the same type WaveFunctionComponent of multiple walkers More...
 
virtual void recompute (const ParticleSet &P)
 recompute the value of the WaveFunctionComponents which require critical accuracy. More...
 
virtual void mw_recompute (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, const std::vector< bool > &recompute) const
 
virtual void evaluateHessian (ParticleSet &P, HessVector &grad_grad_psi_all)
 
virtual void prepareGroup (ParticleSet &P, int ig)
 Prepare internal data for updating WFC correspond to a particle group It should be called before moving particles of a given group. More...
 
virtual void mw_prepareGroup (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int ig) const
 
virtual GradType evalGradWithSpin (ParticleSet &P, int iat, ComplexType &spingrad)
 return the current spin gradient for the iat-th particle Default implementation assumes that WaveFunctionComponent does not explicitly depend on Spin. More...
 
template<CoordsType CT>
void mw_evalGrad (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, const int iat, TWFGrads< CT > &grads_now) const
 compute the current gradients for the iat-th particle of multiple walkers More...
 
virtual void mw_evalGrad (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int iat, std::vector< GradType > &grad_now) const
 compute the current gradients for the iat-th particle of multiple walkers More...
 
virtual GradType evalGradSource (ParticleSet &P, ParticleSet &source, int iat)
 return the logarithmic gradient for the iat-th particle of the source particleset More...
 
virtual GradType evalGradSource (ParticleSet &P, ParticleSet &source, int iat, TinyVector< ParticleSet::ParticleGradient, OHMMS_DIM > &grad_grad, TinyVector< ParticleSet::ParticleLaplacian, OHMMS_DIM > &lapl_grad)
 Adds the gradient w.r.t. More...
 
virtual PsiValue ratioGradWithSpin (ParticleSet &P, int iat, GradType &grad_iat, ComplexType &spingrad_iat)
 evaluate the ratio of the new to old WaveFunctionComponent value and the new spin gradient Default implementation assumes that WaveFunctionComponent does not explicitly depend on Spin. More...
 
template<CoordsType CT>
void mw_ratioGrad (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int iat, std::vector< PsiValue > &ratios, TWFGrads< CT > &grad_new) const
 
virtual void mw_ratioGrad (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int iat, std::vector< PsiValue > &ratios, std::vector< GradType > &grad_new) const
 compute the ratio of the new to old WaveFunctionComponent value and the new gradient of multiple walkers More...
 
virtual void mw_accept_rejectMove (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int iat, const std::vector< bool > &isAccepted, bool safe_to_delay=false) const
 moves of the iat-th particle on some walkers in a batch is accepted. More...
 
virtual void completeUpdates ()
 complete all the delayed or asynchronous operations before leaving the p-by-p move region. More...
 
virtual void mw_completeUpdates (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list) const
 complete all the delayed or asynchronous operations for all the walkers in a batch before leaving the p-by-p move region. More...
 
virtual void mw_calcRatio (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int iat, std::vector< PsiValue > &ratios) const
 compute the ratio of the new to old WaveFunctionComponent value of multiple walkers More...
 
virtual LogValue evaluateGL (const ParticleSet &P, ParticleSet::ParticleGradient &G, ParticleSet::ParticleLaplacian &L, bool fromscratch)
 compute gradients and laplacian of the TWF with respect to each particle. More...
 
virtual void mw_evaluateGL (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, const RefVector< ParticleSet::ParticleGradient > &G_list, const RefVector< ParticleSet::ParticleLaplacian > &L_list, bool fromscratch) const
 evaluate gradients and laplacian of the same type WaveFunctionComponent of multiple walkers More...
 
virtual void createResource (ResourceCollection &collection) const
 initialize a shared resource and hand it to a collection More...
 
virtual void acquireResource (ResourceCollection &collection, const RefVectorWithLeader< WaveFunctionComponent > &wfc_list) const
 acquire a shared resource from a collection More...
 
virtual void releaseResource (ResourceCollection &collection, const RefVectorWithLeader< WaveFunctionComponent > &wfc_list) const
 return a shared resource to a collection More...
 
virtual RealType KECorrection ()
 Return the Chiesa kinetic energy correction. More...
 
virtual bool isOptimizable () const
 if true, this contains optimizable components More...
 
virtual void extractOptimizableObjectRefs (UniqueOptObjRefs &opt_obj_refs)
 extract underlying OptimizableObject references More...
 
virtual void evaluateDerivatives (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, Vector< ValueType > &dhpsioverpsi)=0
 Compute the derivatives of both the log of the wavefunction and kinetic energy with respect to optimizable parameters. More...
 
virtual void evaluateDerivativesWF (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi)
 Compute the derivatives of the log of the wavefunction with respect to optimizable parameters. More...
 
virtual void evaluateGradDerivatives (const ParticleSet::ParticleGradient &G_in, std::vector< ValueType > &dgradlogpsi)
 Calculates the derivatives of $ \nabla \textnormal{log} \psi_f $ with respect to the optimizable parameters, and the dot product of this is then performed with the passed-in G_in gradient vector. More...
 
virtual void finalizeOptimization ()
 
virtual void evaluateRatiosAlltoOne (ParticleSet &P, std::vector< ValueType > &ratios)
 evaluate the ratios of one virtual move with respect to all the particles More...
 
virtual void evaluateRatios (const VirtualParticleSet &VP, std::vector< ValueType > &ratios)
 evaluate ratios to evaluate the non-local PP More...
 
virtual void evaluateSpinorRatios (const VirtualParticleSet &VP, const std::pair< ValueVector, ValueVector > &spinor_multiplier, std::vector< ValueType > &ratios)
 Used by SOECPComponent for faster SOC evaluation. More...
 
virtual void mw_evaluateRatios (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< const VirtualParticleSet > &vp_list, std::vector< std::vector< ValueType >> &ratios) const
 evaluate ratios to evaluate the non-local PP multiple walkers More...
 
virtual void evaluateDerivRatios (const VirtualParticleSet &VP, const opt_variables_type &optvars, std::vector< ValueType > &ratios, Matrix< ValueType > &dratios)
 evaluate ratios to evaluate the non-local PP More...
 
virtual void mw_evalGradWithSpin (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int iat, std::vector< GradType > &grad_now, std::vector< ComplexType > &spingrad_now) const
 compute the current gradients and spin gradients for the iat-th particle of multiple walkers More...
 
virtual void mw_ratioGradWithSpin (const RefVectorWithLeader< WaveFunctionComponent > &wfc_list, const RefVectorWithLeader< ParticleSet > &p_list, int iat, std::vector< PsiValue > &ratios, std::vector< GradType > &grad_new, std::vector< ComplexType > &spingrad_new) const
 compute the ratio of the new to old WaveFunctionComponent value and the new gradient/spingradient of multiple walkers More...
 

Public Attributes

std::vector< RealTypeParticleAlpha
 
std::vector< int > ParticleCenter
 
- Public Attributes inherited from WaveFunctionComponent
int UpdateMode
 current update mode More...
 
opt_variables_type myVars
 list of variables this WaveFunctionComponent handles More...
 
size_t Bytes_in_WFBuffer
 Bytes in WFBuffer. More...
 

Private Attributes

ParticleAttrib< RealTypeU
 
ParticleAttrib< RealTyped2U
 
ParticleAttrib< PosTypedU
 
RealTypeFirstAddressOfdU
 
RealTypeLastAddressOfdU
 
int myTableID
 table index More...
 
ParticleSetCenterRef
 orbital centers More...
 
int NumTargetPtcls
 
int NumCenters
 
RealType curVal
 
RealType curLap
 
PosType curGrad
 

Additional Inherited Members

- Public Types inherited from WaveFunctionComponent
enum  {
  ORB_PBYP_RATIO, ORB_PBYP_ALL, ORB_PBYP_PARTIAL, ORB_WALKER,
  ORB_ALLWALKER
}
 enum for a update mode More...
 
using Walker_t = ParticleSet::Walker_t
 
using WFBufferType = Walker_t::WFBuffer_t
 
using BufferType = Walker_t::Buffer_t
 
using RealMatrix_t = OrbitalSetTraits< RealType >::ValueMatrix
 
using ValueVector = OrbitalSetTraits< ValueType >::ValueVector
 
using ValueMatrix = OrbitalSetTraits< ValueType >::ValueMatrix
 
using GradMatrix = OrbitalSetTraits< ValueType >::GradMatrix
 
using HessType = OrbitalSetTraits< ValueType >::HessType
 
using HessVector = OrbitalSetTraits< ValueType >::HessVector
 
using LogValue = std::complex< QTFull::RealType >
 
using PsiValue = QTFull::ValueType
 
- Public Types inherited from QMCTraits
enum  { DIM = OHMMS_DIM, DIM_VGL = OHMMS_DIM + 2 }
 
using QTBase = QMCTypes< OHMMS_PRECISION, DIM >
 
using QTFull = QMCTypes< OHMMS_PRECISION_FULL, DIM >
 
using RealType = QTBase::RealType
 
using ComplexType = QTBase::ComplexType
 
using ValueType = QTBase::ValueType
 
using PosType = QTBase::PosType
 
using GradType = QTBase::GradType
 
using TensorType = QTBase::TensorType
 
using IndexType = OHMMS_INDEXTYPE
 define other types More...
 
using FullPrecRealType = QTFull::RealType
 
using FullPrecValueType = QTFull::ValueType
 
using PropertySetType = RecordNamedProperty< FullPrecRealType >
 define PropertyList_t More...
 
using PtclGrpIndexes = std::vector< std::pair< int, int > >
 
- Protected Attributes inherited from WaveFunctionComponent
const std::string my_name_
 Name of the object It is required to be different for objects of the same derived type like multiple J1. More...
 
LogValue log_value_
 Current $\log\phi $. More...
 

Detailed Description

A composite Orbital.

Definition at line 27 of file LatticeGaussianProduct.h.

Constructor & Destructor Documentation

◆ LatticeGaussianProduct()

LatticeGaussianProduct ( ParticleSet centers,
ParticleSet ptcls 
)

Definition at line 24 of file LatticeGaussianProduct.cpp.

References ParticleSet::addTable(), LatticeGaussianProduct::CenterRef, LatticeGaussianProduct::d2U, LatticeGaussianProduct::dU, LatticeGaussianProduct::FirstAddressOfdU, ParticleSet::getTotalNum(), LatticeGaussianProduct::LastAddressOfdU, LatticeGaussianProduct::myTableID, LatticeGaussianProduct::NumCenters, LatticeGaussianProduct::NumTargetPtcls, OHMMS_DIM, Vector< T, Alloc >::resize(), and LatticeGaussianProduct::U.

24  : CenterRef(centers)
25 {
26  NumTargetPtcls = ptcls.getTotalNum();
27  NumCenters = centers.getTotalNum();
28  myTableID = ptcls.addTable(CenterRef);
30  dU.resize(NumTargetPtcls);
32  FirstAddressOfdU = &(dU[0][0]);
34 }
void resize(size_type n, Type_t val=Type_t())
Resize the container.
Definition: OhmmsVector.h:166
ParticleSet & CenterRef
orbital centers
#define OHMMS_DIM
Definition: config.h:64

◆ ~LatticeGaussianProduct()

~LatticeGaussianProduct ( )
overridedefault

Member Function Documentation

◆ acceptMove()

void acceptMove ( ParticleSet P,
int  iat,
bool  safe_to_delay = false 
)
overridevirtual

a move for iat-th particle is accepted.

Update the current content.

Parameters
Ptarget ParticleSet
iatindex of the particle whose new position was proposed
safe_to_delayif true, delayed accept is safe.

Implements WaveFunctionComponent.

Definition at line 127 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::curGrad, LatticeGaussianProduct::curLap, LatticeGaussianProduct::curVal, LatticeGaussianProduct::d2U, LatticeGaussianProduct::dU, and LatticeGaussianProduct::U.

◆ copyFromBuffer()

void copyFromBuffer ( ParticleSet P,
WFBufferType buf 
)
overridevirtual

copy the current data from a buffer

Parameters
Pthe ParticleSet to operate on
bufPooledData which stores the data for each walker

copyFromBuffer uses the data stored by registerData

Implements WaveFunctionComponent.

Definition at line 191 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::d2U, Vector< T, Alloc >::first_address(), LatticeGaussianProduct::FirstAddressOfdU, PooledMemory< T_scalar, Alloc >::get(), Vector< T, Alloc >::last_address(), LatticeGaussianProduct::LastAddressOfdU, and LatticeGaussianProduct::U.

◆ evalGrad()

GradType evalGrad ( ParticleSet P,
int  iat 
)
overridevirtual

return the current gradient for the iat-th particle

Parameters
Pquantum particle set
iatparticle index
Returns
the gradient of the iat-th particle

Reimplemented from WaveFunctionComponent.

Definition at line 97 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::curGrad, ParticleSet::getDistTableAB(), LatticeGaussianProduct::myTableID, LatticeGaussianProduct::ParticleAlpha, and LatticeGaussianProduct::ParticleCenter.

98 {
99  const auto& d_table = P.getDistTableAB(myTableID);
100  int icent = ParticleCenter[iat];
101  if (icent == -1)
102  return GradType();
103  RealType a = ParticleAlpha[iat];
104  PosType newdisp = -1.0 * d_table.getTempDispls()[icent];
105  curGrad = -2.0 * a * newdisp;
106  return curGrad;
107 }
QTBase::GradType GradType
Definition: Configuration.h:62
QMCTraits::PosType PosType
QMCTraits::RealType RealType

◆ evaluateDerivatives()

void evaluateDerivatives ( ParticleSet P,
const opt_variables_type optvars,
Vector< ValueType > &  dlogpsi,
Vector< ValueType > &  dhpsioverpsi 
)
inlineoverride

Definition at line 73 of file LatticeGaussianProduct.h.

77  {}

◆ evaluateLog()

Parameters
Pinput configuration containing N particles
Ga vector containing N gradients
La vector containing N laplacians
Returns
The wavefunction value $exp(-J({\bf R}))$

Upon exit, the gradient $G[i]={\bf \nabla}_i J({\bf R})$ and the laplacian $L[i]=\nabla^2_i J({\bf R})$ are accumulated. While evaluating the value of the Jastrow for a set of particles add the gradient and laplacian contribution of the Jastrow to G(radient) and L(aplacian) for local energy calculations such that

\[ G[i]+={\bf \nabla}_i J({\bf R}) \]

and

\[ L[i]+=\nabla^2_i J({\bf R}). \]

Implements WaveFunctionComponent.

Definition at line 52 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::d2U, LatticeGaussianProduct::dU, ParticleSet::getDistTableAB(), WaveFunctionComponent::log_value_, LatticeGaussianProduct::myTableID, LatticeGaussianProduct::NumTargetPtcls, LatticeGaussianProduct::ParticleAlpha, and LatticeGaussianProduct::U.

55 {
56  const auto& d_table = P.getDistTableAB(myTableID);
57  int icent = 0;
58  log_value_ = 0.0;
59  RealType dist = 0.0;
60  PosType disp = 0.0;
61  for (int iat = 0; iat < NumTargetPtcls; iat++)
62  {
63  U[iat] = 0.0;
64  dU[iat] = 0.0;
65  d2U[iat] = 0.0;
66  RealType a = ParticleAlpha[iat];
67  if (a > 0.0)
68  {
69  dist = d_table.getDistRow(iat)[icent];
70  disp = -1.0 * d_table.getDisplRow(iat)[icent];
71  log_value_ -= a * dist * dist;
72  U[iat] += a * dist * dist;
73  G[iat] -= 2.0 * a * disp;
74  L[iat] -= 6.0 * a;
75  icent++;
76  }
77  }
78  return log_value_;
79 }
QMCTraits::PosType PosType
QMCTraits::RealType RealType

◆ evaluateLogAndStore()

void evaluateLogAndStore ( const ParticleSet P,
ParticleSet::ParticleGradient dG,
ParticleSet::ParticleLaplacian dL 
)

Definition at line 134 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::d2U, LatticeGaussianProduct::dU, ParticleSet::getDistTableAB(), WaveFunctionComponent::log_value_, LatticeGaussianProduct::myTableID, LatticeGaussianProduct::NumTargetPtcls, LatticeGaussianProduct::ParticleAlpha, and LatticeGaussianProduct::U.

Referenced by LatticeGaussianProduct::registerData(), and LatticeGaussianProduct::updateBuffer().

137 {
138  const auto& d_table = P.getDistTableAB(myTableID);
139  RealType dist = 0.0;
140  PosType disp = 0.0;
141  int icent = 0;
142  log_value_ = 0.0;
143  U = 0.0;
144  dU = 0.0;
145  d2U = 0.0;
146  for (int iat = 0; iat < NumTargetPtcls; iat++)
147  {
148  RealType a = ParticleAlpha[iat];
149  if (a > 0.0)
150  {
151  dist = d_table.getDistRow(iat)[icent];
152  disp = -1.0 * d_table.getDisplRow(iat)[icent];
153  log_value_ -= a * dist * dist;
154  U[iat] += a * dist * dist;
155  dU[iat] -= 2.0 * a * disp;
156  d2U[iat] -= 6.0 * a;
157  dG[iat] -= 2.0 * a * disp;
158  dL[iat] -= 6.0 * a;
159  icent++;
160  }
161  }
162 }
QMCTraits::PosType PosType
QMCTraits::RealType RealType

◆ getClassName()

std::string getClassName ( ) const
inlineoverridevirtual

return class name

Implements WaveFunctionComponent.

Definition at line 49 of file LatticeGaussianProduct.h.

49 { return "LatticeGaussianProduct"; }

◆ makeClone()

std::unique_ptr< WaveFunctionComponent > makeClone ( ParticleSet tqp) const
overridevirtual

make clone

Parameters
tqptarget Quantum ParticleSet
deepcopyif true, make a decopy

If not true, return a proxy class

Reimplemented from WaveFunctionComponent.

Definition at line 198 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::CenterRef, LatticeGaussianProduct::ParticleAlpha, and LatticeGaussianProduct::ParticleCenter.

199 {
200  auto j1copy = std::make_unique<LatticeGaussianProduct>(CenterRef, tqp);
201  j1copy->ParticleAlpha = ParticleAlpha;
202  j1copy->ParticleCenter = ParticleCenter;
203  return j1copy;
204 }
ParticleSet & CenterRef
orbital centers

◆ ratio()

PsiValue ratio ( ParticleSet P,
int  iat 
)
overridevirtual

evaluate the ratio $exp(U(iat)-U_0(iat))$

Parameters
Pactive particle set
iatparticle that has been moved.

Implements WaveFunctionComponent.

Definition at line 85 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::curVal, qmcplusplus::exp(), ParticleSet::getDistTableAB(), LatticeGaussianProduct::myTableID, LatticeGaussianProduct::ParticleAlpha, LatticeGaussianProduct::ParticleCenter, and LatticeGaussianProduct::U.

86 {
87  const auto& d_table = P.getDistTableAB(myTableID);
88  int icent = ParticleCenter[iat];
89  if (icent == -1)
90  return 1.0;
91  RealType newdist = d_table.getTempDists()[icent];
92  curVal = ParticleAlpha[iat] * (newdist * newdist);
93  return std::exp(static_cast<PsiValue>(U[iat] - curVal));
94 }
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
QMCTraits::RealType RealType

◆ ratioGrad()

PsiValue ratioGrad ( ParticleSet P,
int  iat,
GradType grad_iat 
)
overridevirtual

evaluate the ratio of the new to old WaveFunctionComponent value and the new gradient

Parameters
Pthe active ParticleSet
iatthe index of a particle
grad_iatGradient for the active particle

Reimplemented from WaveFunctionComponent.

Definition at line 110 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::curGrad, LatticeGaussianProduct::curVal, qmcplusplus::exp(), ParticleSet::getDistTableAB(), LatticeGaussianProduct::myTableID, LatticeGaussianProduct::ParticleAlpha, LatticeGaussianProduct::ParticleCenter, and LatticeGaussianProduct::U.

111 {
112  const auto& d_table = P.getDistTableAB(myTableID);
113  int icent = ParticleCenter[iat];
114  if (icent == -1)
115  return 1.0;
116  RealType a = ParticleAlpha[iat];
117  RealType newdist = d_table.getTempDists()[icent];
118  PosType newdisp = -1.0 * d_table.getTempDispls()[icent];
119  curVal = a * newdist * newdist;
120  curGrad = -2.0 * a * newdisp;
121  grad_iat += curGrad;
122  return std::exp(static_cast<PsiValue>(U[iat] - curVal));
123 }
QMCTraits::PosType PosType
MakeReturn< UnaryNode< FnExp, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t exp(const Vector< T1, C1 > &l)
QMCTraits::RealType RealType

◆ registerData()

void registerData ( ParticleSet P,
WFBufferType buf 
)
overridevirtual

equivalent to evalaute with additional data management

Implements WaveFunctionComponent.

Definition at line 165 of file LatticeGaussianProduct.cpp.

References PooledMemory< T_scalar, Alloc >::add(), Vector< T, Alloc >::begin(), LatticeGaussianProduct::d2U, Vector< T, Alloc >::end(), LatticeGaussianProduct::evaluateLogAndStore(), LatticeGaussianProduct::FirstAddressOfdU, ParticleSet::G, ParticleSet::L, LatticeGaussianProduct::LastAddressOfdU, and LatticeGaussianProduct::U.

166 {
167  evaluateLogAndStore(P, P.G, P.L);
168  // Add U, d2U and dU. Keep the order!!!
169  buf.add(U.begin(), U.end());
170  buf.add(d2U.begin(), d2U.end());
172 }
void evaluateLogAndStore(const ParticleSet &P, ParticleSet::ParticleGradient &dG, ParticleSet::ParticleLaplacian &dL)

◆ restore()

void restore ( int  iat)
overridevirtual

If a move for iat-th particle is rejected, restore to the content.

Parameters
iatindex of the particle whose new position was proposed

Ye: hopefully we can gradually move away from restore

Implements WaveFunctionComponent.

Definition at line 125 of file LatticeGaussianProduct.cpp.

125 {}

◆ updateBuffer()

LatticeGaussianProduct::LogValue updateBuffer ( ParticleSet P,
WFBufferType buf,
bool  fromscratch 
)
overridevirtual

For particle-by-particle move.

Put the objects of this class in the walker buffer or forward the memory cursor.

Parameters
Pparticle set
bufAnonymous storage
fromscratchrequest recomputing the precision critical pieces of wavefunction from scratch
Returns
log value of the wavefunction.

Implements WaveFunctionComponent.

Definition at line 174 of file LatticeGaussianProduct.cpp.

References LatticeGaussianProduct::d2U, LatticeGaussianProduct::evaluateLogAndStore(), Vector< T, Alloc >::first_address(), LatticeGaussianProduct::FirstAddressOfdU, ParticleSet::G, ParticleSet::L, Vector< T, Alloc >::last_address(), LatticeGaussianProduct::LastAddressOfdU, WaveFunctionComponent::log_value_, PooledMemory< T_scalar, Alloc >::put(), and LatticeGaussianProduct::U.

177 {
178  evaluateLogAndStore(P, P.G, P.L);
179  buf.put(U.first_address(), U.last_address());
180  buf.put(d2U.first_address(), d2U.last_address());
182  return log_value_;
183 }
void evaluateLogAndStore(const ParticleSet &P, ParticleSet::ParticleGradient &dG, ParticleSet::ParticleLaplacian &dL)
pointer last_address()
Definition: OhmmsVector.h:258
pointer first_address()
Definition: OhmmsVector.h:255

Member Data Documentation

◆ CenterRef

ParticleSet& CenterRef
private

◆ curGrad

◆ curLap

RealType curLap
private

Definition at line 38 of file LatticeGaussianProduct.h.

Referenced by LatticeGaussianProduct::acceptMove().

◆ curVal

◆ d2U

◆ dU

◆ FirstAddressOfdU

◆ LastAddressOfdU

◆ myTableID

◆ NumCenters

int NumCenters
private

◆ NumTargetPtcls

◆ ParticleAlpha

◆ ParticleCenter

◆ U


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