QMCPACK
RPAJastrow Class Reference

JastrowBuilder using RPA functor Modification of RPAJastrow. More...

+ Inheritance diagram for RPAJastrow:
+ Collaboration diagram for RPAJastrow:

Public Member Functions

 RPAJastrow (ParticleSet &target)
 
 ~RPAJastrow () override
 
std::string getClassName () const override
 return class name More...
 
bool put (xmlNodePtr cur)
 
void buildOrbital (const std::string &name, const std::string &UL, const std::string &US, const std::string &RF, RealType R, RealType K)
 
void makeShortRange ()
 
void makeLongRange ()
 
void setHandler (std::unique_ptr< HandlerType > Handler)
 
bool isOptimizable () const override
 if true, this contains optimizable components More...
 
void checkOutVariables (const opt_variables_type &o) override
 check out optimizable variables More...
 
void extractOptimizableObjectRefs (UniqueOptObjRefs &opt_obj_refs) override
 extract underlying OptimizableObject references More...
 
LogValue evaluateLog (const ParticleSet &P, ParticleSet::ParticleGradient &G, ParticleSet::ParticleLaplacian &L) override
 evaluate the value of the WaveFunctionComponent from scratch More...
 
PsiValue ratio (ParticleSet &P, int iat) override
 evaluate the ratio of the new to old WaveFunctionComponent value 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...
 
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
 For particle-by-particle move. More...
 
LogValue updateBuffer (ParticleSet &P, WFBufferType &buf, bool fromscratch=false) override
 For particle-by-particle move. More...
 
void copyFromBuffer (ParticleSet &P, WFBufferType &buf) override
 For particle-by-particle move. More...
 
std::unique_ptr< WaveFunctionComponentmakeClone (ParticleSet &tqp) const override
 this is a great deal of logic for make clone I'm wondering what is going on More...
 
void evaluateDerivatives (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, Vector< ValueType > &dhpsioverpsi) override
 Compute the derivatives of both the log of the wavefunction and kinetic energy with respect to optimizable parameters. More...
 
- 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 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 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...
 

Private Types

using HandlerType = LRHandlerBase
 
using FuncType = BsplineFunctor< RealType >
 
using GridType = LinearGrid< RealType >
 

Private Attributes

bool IgnoreSpin
 
bool DropLongRange
 
bool DropShortRange
 
RealType Rs
 
RealType Kc
 
RealType Rcut
 
std::string ID_Rs
 
std::string rpafunc
 
std::string MyName
 
std::unique_ptr< HandlerTypemyHandler
 main handler More...
 
kSpaceJastrowLongRangeRPA
 object to handle the long-range part More...
 
WaveFunctionComponentShortRangeRPA
 
FuncTypenfunc
 numerical function owned by ShortRangeRPA More...
 
GridTypemyGrid
 
ParticleSettargetPtcl
 
UPtrVector< WaveFunctionComponentPsi
 A list of WaveFunctionComponent*. More...
 

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 > >
 
- 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...
 
- 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

JastrowBuilder using RPA functor Modification of RPAJastrow.

Definition at line 32 of file RPAJastrow.h.

Member Typedef Documentation

◆ FuncType

using FuncType = BsplineFunctor<RealType>
private

Definition at line 35 of file RPAJastrow.h.

◆ GridType

using GridType = LinearGrid<RealType>
private

Definition at line 36 of file RPAJastrow.h.

◆ HandlerType

using HandlerType = LRHandlerBase
private

Definition at line 34 of file RPAJastrow.h.

Constructor & Destructor Documentation

◆ RPAJastrow()

RPAJastrow ( ParticleSet target)

Definition at line 34 of file RPAJastrow.cpp.

34 : targetPtcl(target) {}
ParticleSet & targetPtcl
Definition: RPAJastrow.h:119

◆ ~RPAJastrow()

~RPAJastrow ( )
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 246 of file RPAJastrow.cpp.

References RPAJastrow::Psi.

247 {
248  for (int i = 0; i < Psi.size(); i++)
249  Psi[i]->acceptMove(P, iat, safe_to_delay);
250 }
void acceptMove(ParticleSet &P, int iat, bool safe_to_delay=false) override
a move for iat-th particle is accepted.
Definition: RPAJastrow.cpp:246
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121

◆ buildOrbital()

void buildOrbital ( const std::string &  name,
const std::string &  UL,
const std::string &  US,
const std::string &  RF,
RealType  R,
RealType  K 
)

Definition at line 67 of file RPAJastrow.cpp.

References qmcplusplus::app_log(), RPAJastrow::DropLongRange, RPAJastrow::DropShortRange, ParticleSet::getLattice(), ParticleSet::getSimulationCell(), ParticleSet::getTotalNum(), RPAJastrow::ID_Rs, qmcplusplus::Units::energy::K, RPAJastrow::Kc, RPAJastrow::makeLongRange(), RPAJastrow::makeShortRange(), RPAJastrow::myHandler, RPAJastrow::MyName, qmcplusplus::pow(), RPAJastrow::rpafunc, RPAJastrow::Rs, and RPAJastrow::targetPtcl.

Referenced by RPAJastrow::put().

73 {
74  std::string ID_Rs = "RPA_rs";
75  MyName = name;
76  std::string useL = UL;
77  std::string useS = US;
78  rpafunc = RF;
79  Rs = R;
80  Kc = K;
81  app_log() << std::endl << " LongRangeForm is " << rpafunc << std::endl;
82  DropLongRange = (useL == "no");
83  DropShortRange = (useS == "no");
84  RealType tlen =
85  std::pow(3.0 / 4.0 / M_PI * targetPtcl.getLattice().Volume / static_cast<RealType>(targetPtcl.getTotalNum()),
86  1.0 / 3.0);
87  if (Rs < 0)
88  {
89  if (targetPtcl.getLattice().SuperCellEnum)
90  {
91  Rs = tlen;
92  }
93  else
94  {
95  std::cout << " Error finding rs. Is this an open system?!" << std::endl;
96  Rs = 100.0;
97  }
98  }
99  int indx = targetPtcl.getSimulationCell().getKLists().ksq.size() - 1;
100  double Kc_max = std::pow(targetPtcl.getSimulationCell().getKLists().ksq[indx], 0.5);
101  if (Kc < 0)
102  {
103  Kc = 2.0 * std::pow(2.25 * M_PI, 1.0 / 3.0) / tlen;
104  }
105  if (Kc > Kc_max)
106  {
107  Kc = Kc_max;
108  app_log() << " Kc set too high. Resetting to the maximum value" << std::endl;
109  }
110  app_log() << " RPAJastrowBuilder::addTwoBodyPart Rs = " << Rs << " Kc= " << Kc << std::endl;
111  if (rpafunc == "yukawa" || rpafunc == "breakup")
112  myHandler = std::make_unique<LRHandlerTemp<YukawaBreakup<RealType>, LPQHIBasis>>(targetPtcl, Kc);
113  else if (rpafunc == "rpa")
114  myHandler = std::make_unique<LRRPAHandlerTemp<RPABreakup<RealType>, LPQHIBasis>>(targetPtcl, Kc);
115  else if (rpafunc == "dyukawa")
116  myHandler = std::make_unique<LRHandlerTemp<DerivYukawaBreakup<RealType>, LPQHIBasis>>(targetPtcl, Kc);
117  else if (rpafunc == "drpa")
118  myHandler = std::make_unique<LRRPAHandlerTemp<DerivRPABreakup<RealType>, LPQHIBasis>>(targetPtcl, Kc);
119  else
120  throw std::invalid_argument("RPAJastrowBuilder::buildOrbital: Unrecognized rpa function type.\n");
121  myHandler->Breakup(targetPtcl, Rs);
122  app_log() << " Maximum K shell " << myHandler->MaxKshell << std::endl;
123  app_log() << " Number of k vectors " << myHandler->Fk.size() << std::endl;
124  if (!DropLongRange)
125  {
126  makeLongRange();
127  app_log() << " Using LongRange part" << std::endl;
128  }
129  if (!DropShortRange)
130  {
131  makeShortRange();
132  app_log() << " Using ShortRange part" << std::endl;
133  }
134 }
size_t getTotalNum() const
Definition: ParticleSet.h:493
std::ostream & app_log()
Definition: OutputManager.h:65
MakeReturn< BinaryNode< FnPow, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t, typename CreateLeaf< Vector< T2, C2 > >::Leaf_t > >::Expression_t pow(const Vector< T1, C1 > &l, const Vector< T2, C2 > &r)
const auto & getSimulationCell() const
Definition: ParticleSet.h:250
std::unique_ptr< HandlerType > myHandler
main handler
Definition: RPAJastrow.h:110
ParticleSet & targetPtcl
Definition: RPAJastrow.h:119
QMCTraits::RealType RealType
const auto & getLattice() const
Definition: ParticleSet.h:251

◆ checkOutVariables()

void checkOutVariables ( const opt_variables_type o)
overridevirtual

check out optimizable variables

Reimplemented from WaveFunctionComponent.

Definition at line 203 of file RPAJastrow.cpp.

References WaveFunctionComponent::checkOutVariables(), kSpaceJastrow::checkOutVariables(), RPAJastrow::LongRangeRPA, and RPAJastrow::ShortRangeRPA.

204 {
207 }
kSpaceJastrow * LongRangeRPA
object to handle the long-range part
Definition: RPAJastrow.h:112
WaveFunctionComponent * ShortRangeRPA
Definition: RPAJastrow.h:115
virtual void checkOutVariables(const opt_variables_type &active)
check out variational optimizable variables
void checkOutVariables(const opt_variables_type &active) override
check out variational optimizable variables

◆ copyFromBuffer()

void copyFromBuffer ( ParticleSet P,
WFBufferType buf 
)
overridevirtual

For particle-by-particle move.

Copy data or attach memory from a walker buffer to the objects of this class. The log value, P.G and P.L contribution from the objects of this class are also added.

Parameters
Pparticle set
bufAnonymous storage

Implements WaveFunctionComponent.

Definition at line 272 of file RPAJastrow.cpp.

References RPAJastrow::Psi.

273 {
274  for (int i = 0; i < Psi.size(); i++)
275  Psi[i]->copyFromBuffer(P, buf);
276 }
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121
void copyFromBuffer(ParticleSet &P, WFBufferType &buf) override
For particle-by-particle move.
Definition: RPAJastrow.cpp:272

◆ evalGrad()

RPAJastrow::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 227 of file RPAJastrow.cpp.

References RPAJastrow::Psi.

228 {
229  GradType grad(0);
230  for (int i = 0; i < Psi.size(); i++)
231  grad += Psi[i]->evalGrad(P, iat);
232  return grad;
233 }
QTBase::GradType GradType
Definition: Configuration.h:62
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121
GradType evalGrad(ParticleSet &P, int iat) override
return the current gradient for the iat-th particle
Definition: RPAJastrow.cpp:227

◆ evaluateDerivatives()

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

Compute the derivatives of both the log of the wavefunction and kinetic energy with respect to optimizable parameters.

Parameters
Pparticle set
optvarsoptimizable parameters
dlogpsiarray of derivatives of the log of the wavefunction. Add the contribution from this component.
dhpsioverpsiarray of Hamiltonian derivatives. Add the kinetic energy derivatives contribution from this component. $ -\frac{1}{2}{\partial}_\alpha \tilde L - G \cdot {\partial}_\alpha \tilde G $. $ \tilde L $ and $ \tilde G $ are from this WaveFunctionComponent. $ G $ is from TrialWaveFunction. The 1/m factor is applied in TrialWaveFunction. This is a bug when the particle set doesn't hold equal mass particles.

Implements WaveFunctionComponent.

Definition at line 90 of file RPAJastrow.h.

94  {}

◆ evaluateLog()

RPAJastrow::LogValue evaluateLog ( const ParticleSet P,
ParticleSet::ParticleGradient G,
ParticleSet::ParticleLaplacian L 
)
overridevirtual

evaluate the value of the WaveFunctionComponent from scratch

Parameters
[in]Pactive ParticleSet
[out]GGradients, $\nabla\ln\Psi$
[out]LLaplacians, $\nabla^2\ln\Psi$
Returns
the log value

Mainly for walker-by-walker move. The initial stage of particle-by-particle move also uses this. causes complete state update in WFC's

Implements WaveFunctionComponent.

Definition at line 209 of file RPAJastrow.cpp.

References WaveFunctionComponent::log_value_, and RPAJastrow::Psi.

212 {
213  log_value_ = 0.0;
214  for (int i = 0; i < Psi.size(); i++)
215  log_value_ += Psi[i]->evaluateLog(P, G, L);
216  return log_value_;
217 }
LogValue evaluateLog(const ParticleSet &P, ParticleSet::ParticleGradient &G, ParticleSet::ParticleLaplacian &L) override
evaluate the value of the WaveFunctionComponent from scratch
Definition: RPAJastrow.cpp:209
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121

◆ extractOptimizableObjectRefs()

void extractOptimizableObjectRefs ( UniqueOptObjRefs opt_obj_refs)
inlineoverridevirtual

extract underlying OptimizableObject references

Parameters
opt_obj_refsaggregated list of optimizable object references

Reimplemented from WaveFunctionComponent.

Definition at line 64 of file RPAJastrow.h.

References WaveFunctionComponent::extractOptimizableObjectRefs(), RPAJastrow::LongRangeRPA, UniqueOptObjRefs::push_back(), and RPAJastrow::ShortRangeRPA.

65  {
66  opt_obj_refs.push_back(*LongRangeRPA);
68  }
kSpaceJastrow * LongRangeRPA
object to handle the long-range part
Definition: RPAJastrow.h:112
WaveFunctionComponent * ShortRangeRPA
Definition: RPAJastrow.h:115
virtual void extractOptimizableObjectRefs(UniqueOptObjRefs &opt_obj_refs)
extract underlying OptimizableObject references

◆ getClassName()

std::string getClassName ( ) const
inlineoverridevirtual

return class name

Implements WaveFunctionComponent.

Definition at line 43 of file RPAJastrow.h.

43 { return "RPAJastrow"; }

◆ isOptimizable()

bool isOptimizable ( ) const
inlineoverridevirtual

if true, this contains optimizable components

Reimplemented from WaveFunctionComponent.

Definition at line 59 of file RPAJastrow.h.

59 { return true; }

◆ makeClone()

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

this is a great deal of logic for make clone I'm wondering what is going on

Reimplemented from WaveFunctionComponent.

Definition at line 280 of file RPAJastrow.cpp.

References RPAJastrow::DropLongRange, RPAJastrow::DropShortRange, RPAJastrow::Kc, RPAJastrow::myHandler, RPAJastrow::Rcut, and RPAJastrow::rpafunc.

281 {
282  std::unique_ptr<HandlerType> tempHandler;
283  if (rpafunc == "yukawa" || rpafunc == "breakup")
284  {
285  tempHandler =
286  std::make_unique<LRHandlerTemp<YukawaBreakup<RealType>, LPQHIBasis>>(dynamic_cast<const LRHandlerTemp<
287  YukawaBreakup<RealType>, LPQHIBasis>&>(
288  *myHandler),
289  tpq);
290  }
291  else if (rpafunc == "rpa")
292  {
293  tempHandler =
294  std::make_unique<LRRPAHandlerTemp<RPABreakup<RealType>, LPQHIBasis>>(dynamic_cast<const LRRPAHandlerTemp<
295  RPABreakup<RealType>, LPQHIBasis>&>(
296  *myHandler),
297  tpq);
298  }
299  else if (rpafunc == "dyukawa")
300  {
301  tempHandler = std::make_unique<LRHandlerTemp<
302  DerivYukawaBreakup<RealType>,
303  LPQHIBasis>>(dynamic_cast<const LRHandlerTemp<DerivYukawaBreakup<RealType>, LPQHIBasis>&>(*myHandler), tpq);
304  }
305  else if (rpafunc == "drpa")
306  {
307  tempHandler = std::make_unique<LRRPAHandlerTemp<
308  DerivRPABreakup<RealType>,
309  LPQHIBasis>>(dynamic_cast<const LRRPAHandlerTemp<DerivRPABreakup<RealType>, LPQHIBasis>&>(*myHandler), tpq);
310  }
311 
312  auto myClone = std::make_unique<RPAJastrow>(tpq);
313  myClone->Rcut = Rcut;
314  myClone->Kc = Kc;
315  myClone->setHandler(std::move(tempHandler));
316  if (!DropLongRange)
317  myClone->makeLongRange();
318  if (!DropShortRange)
319  myClone->makeShortRange();
320  return myClone;
321 }
std::unique_ptr< HandlerType > myHandler
main handler
Definition: RPAJastrow.h:110

◆ makeLongRange()

void makeLongRange ( )

Definition at line 136 of file RPAJastrow.cpp.

References ParticleSet::getLattice(), kSpaceJastrow::ISOTROPIC, RPAJastrow::Kc, RPAJastrow::LongRangeRPA, RPAJastrow::myHandler, RPAJastrow::Psi, kSpaceJastrow::setCoefficients(), and RPAJastrow::targetPtcl.

Referenced by RPAJastrow::buildOrbital().

137 {
138  // create two-body kSpaceJastrow
139  kSpaceJastrow::SymmetryType oneBodySymm, twoBodySymm;
140  bool oneBodySpin, twoBodySpin;
141  oneBodySymm = kSpaceJastrow::ISOTROPIC;
142  twoBodySymm = kSpaceJastrow::ISOTROPIC;
143  oneBodySpin = false;
144  twoBodySpin = false;
145  auto LongRangeRPA_uptr =
146  std::make_unique<kSpaceJastrow>(targetPtcl, targetPtcl, oneBodySymm, -1, "cG1", oneBodySpin, // no one-body part
147  twoBodySymm, Kc, "cG2", twoBodySpin);
148  LongRangeRPA = LongRangeRPA_uptr.get();
149  // fill in CG2 coefficients
150  std::vector<RealType> oneBodyCoefs, twoBodyCoefs;
151  twoBodyCoefs.resize(myHandler->MaxKshell);
152  // need to cancel prefactor in kSpaceJastrow
153  RealType prefactorInv = -targetPtcl.getLattice().Volume;
154  for (size_t is = 0; is < myHandler->MaxKshell; is++)
155  {
156  twoBodyCoefs[is] = prefactorInv * myHandler->Fk_symm[is];
157  }
158  LongRangeRPA->setCoefficients(oneBodyCoefs, twoBodyCoefs);
159  Psi.push_back(std::move(LongRangeRPA_uptr));
160 }
kSpaceJastrow * LongRangeRPA
object to handle the long-range part
Definition: RPAJastrow.h:112
void setCoefficients(std::vector< RealType > &oneBodyCoefs, std::vector< RealType > &twoBodyCoefs)
std::unique_ptr< HandlerType > myHandler
main handler
Definition: RPAJastrow.h:110
ParticleSet & targetPtcl
Definition: RPAJastrow.h:119
QMCTraits::RealType RealType
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121
const auto & getLattice() const
Definition: ParticleSet.h:251

◆ makeShortRange()

void makeShortRange ( )

Definition at line 162 of file RPAJastrow.cpp.

References qmcplusplus::app_log(), qmcplusplus::Units::charge::e, BsplineFunctor< REAL >::initialize(), WaveFunctionComponent::my_name_, RPAJastrow::myHandler, RPAJastrow::nfunc, RPAJastrow::Psi, RPAJastrow::Rcut, ShortRangePartAdapter< T >::setRmax(), RPAJastrow::ShortRangeRPA, and RPAJastrow::targetPtcl.

Referenced by RPAJastrow::buildOrbital().

163 {
164  // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
165  // !!! WARNING: tiny, nparam, npts should be input parameters !!!
166  // !!! fix in a future PR! !!!
167  // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
168  app_log() << " Adding Short Range part of RPA function" << std::endl;
169  //short-range uses realHandler
170  RealType tiny = 1e-6;
171  Rcut = myHandler->get_rc() - tiny;
172  //create numerical functor of type BsplineFunctor<RealType>.
173  auto nfunc_uptr = std::make_unique<FuncType>(my_name_ + "_short");
174  nfunc = nfunc_uptr.get();
175  ShortRangePartAdapter<RealType> SRA(myHandler.get());
176  SRA.setRmax(Rcut);
177  auto j2 = std::make_unique<TwoBodyJastrow<BsplineFunctor<RealType>>>("RPA", targetPtcl, false);
178  size_t nparam = 12; // number of Bspline parameters
179  size_t npts = 100; // number of 1D grid points for basis functions
180  RealType cusp = SRA.df(0);
181  RealType delta = Rcut / static_cast<double>(npts);
182  std::vector<RealType> X(npts + 1), Y(npts + 1);
183  for (size_t i = 0; i < npts; ++i)
184  {
185  X[i] = i * delta;
186  Y[i] = SRA.evaluate(X[i]);
187  }
188  X[npts] = npts * delta;
189  Y[npts] = 0.0;
190  std::string functype = "rpa";
191  std::string useit = "no";
192  nfunc->initialize(nparam, X, Y, cusp, Rcut + tiny, functype, useit);
193  for (size_t i = 0; i < npts; ++i)
194  {
195  X[i] = i * delta;
196  Y[i] = SRA.evaluate(X[i]);
197  }
198  j2->addFunc(0, 0, std::move(nfunc_uptr));
199  ShortRangeRPA = j2.get();
200  Psi.push_back(std::move(j2));
201 }
const std::string my_name_
Name of the object It is required to be different for objects of the same derived type like multiple ...
WaveFunctionComponent * ShortRangeRPA
Definition: RPAJastrow.h:115
std::ostream & app_log()
Definition: OutputManager.h:65
Real df(Real r) override
evaluate the first derivative
std::unique_ptr< HandlerType > myHandler
main handler
Definition: RPAJastrow.h:110
ParticleSet & targetPtcl
Definition: RPAJastrow.h:119
void initialize(int numPoints, std::vector< Real > &x, std::vector< Real > &y, REAL cusp, REAL rcut, std::string &id, std::string &optimize)
QMCTraits::RealType RealType
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121
FuncType * nfunc
numerical function owned by ShortRangeRPA
Definition: RPAJastrow.h:117

◆ put()

bool put ( xmlNodePtr  cur)

Definition at line 38 of file RPAJastrow.cpp.

References OhmmsAttributeSet::add(), ParameterSet::add(), qmcplusplus::app_log(), RPAJastrow::buildOrbital(), RPAJastrow::ID_Rs, RPAJastrow::Kc, RPAJastrow::MyName, ParameterSet::put(), OhmmsAttributeSet::put(), RPAJastrow::rpafunc, and RPAJastrow::Rs.

39 {
40  ReportEngine PRE("RPAJastrow", "put");
41  app_log() << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!\n";
42  app_log() << "!!! WARNING: RPAJastrow is not fully tested for production !!!\n";
43  app_log() << "!!! level calculations. Use at your own risk! !!!\n";
44  app_log() << "!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!\n";
45  //capture attribute jastrow/@name
46  MyName = "RPA_Jee";
47  std::string useL = "yes";
48  std::string useS = "yes";
49  rpafunc = "rpa";
51  a.add(MyName, "name");
52  a.add(useL, "longrange");
53  a.add(useS, "shortrange");
54  a.add(rpafunc, "function");
55  a.put(cur);
56  Rs = -1.0;
57  Kc = -1.0;
58  std::string ID_Rs = "RPA_rs";
59  ParameterSet params;
60  params.add(Rs, "rs");
61  params.add(Kc, "kc");
62  params.put(cur);
63  buildOrbital(MyName, useL, useS, rpafunc, Rs, Kc);
64  return true;
65 }
std::ostream & app_log()
Definition: OutputManager.h:65
bool put(xmlNodePtr cur)
assign attributes to the set
Definition: AttributeSet.h:55
bool put(std::istream &is) override
read from std::istream
Definition: ParameterSet.h:42
class to handle a set of parameters
Definition: ParameterSet.h:27
class to handle a set of attributes of an xmlNode
Definition: AttributeSet.h:24
void add(PDT &aparam, const std::string &aname_in, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new parameter corresponding to an xmlNode <parameter>
void buildOrbital(const std::string &name, const std::string &UL, const std::string &US, const std::string &RF, RealType R, RealType K)
Definition: RPAJastrow.cpp:67
void add(PDT &aparam, const std::string &aname, std::vector< PDT > candidate_values={}, TagStatus status=TagStatus::OPTIONAL)
add a new attribute
Definition: AttributeSet.h:42

◆ ratio()

RPAJastrow::PsiValue ratio ( ParticleSet P,
int  iat 
)
overridevirtual

evaluate the ratio of the new to old WaveFunctionComponent value

Parameters
Pthe active ParticleSet
iatthe index of a particle
Returns
$ \psi( \{ {\bf R}^{'} \} )/ \psi( \{ {\bf R}\})$

Specialized for particle-by-particle move

Implements WaveFunctionComponent.

Definition at line 219 of file RPAJastrow.cpp.

References RPAJastrow::Psi.

220 {
221  ValueType r(1.0);
222  for (int i = 0; i < Psi.size(); i++)
223  r *= Psi[i]->ratio(P, iat);
224  return static_cast<PsiValue>(r);
225 }
PsiValue ratio(ParticleSet &P, int iat) override
evaluate the ratio of the new to old WaveFunctionComponent value
Definition: RPAJastrow.cpp:219
QTBase::ValueType ValueType
Definition: Configuration.h:60
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121

◆ ratioGrad()

RPAJastrow::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 235 of file RPAJastrow.cpp.

References RPAJastrow::Psi.

236 {
237  ValueType r(1);
238  for (int i = 0; i < Psi.size(); i++)
239  {
240  r *= Psi[i]->ratioGrad(P, iat, grad_iat);
241  }
242  return static_cast<PsiValue>(r);
243 }
QTBase::ValueType ValueType
Definition: Configuration.h:60
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121

◆ registerData()

void registerData ( ParticleSet P,
WFBufferType buf 
)
overridevirtual

For particle-by-particle move.

Requests space in the buffer based on the data type sizes of the objects in this class.

Parameters
Pparticle set
bufAnonymous storage

Implements WaveFunctionComponent.

Definition at line 258 of file RPAJastrow.cpp.

References RPAJastrow::Psi.

259 {
260  for (int i = 0; i < Psi.size(); i++)
261  Psi[i]->registerData(P, buf);
262 }
void registerData(ParticleSet &P, WFBufferType &buf) override
For particle-by-particle move.
Definition: RPAJastrow.cpp:258
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121

◆ 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 252 of file RPAJastrow.cpp.

References RPAJastrow::Psi.

253 {
254  for (int i = 0; i < Psi.size(); i++)
255  Psi[i]->restore(iat);
256 }
void restore(int iat) override
If a move for iat-th particle is rejected, restore to the content.
Definition: RPAJastrow.cpp:252
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121

◆ setHandler()

void setHandler ( std::unique_ptr< HandlerType Handler)
inline

Definition at line 57 of file RPAJastrow.h.

References RPAJastrow::myHandler.

57 { myHandler = std::move(Handler); };
std::unique_ptr< HandlerType > myHandler
main handler
Definition: RPAJastrow.h:110

◆ updateBuffer()

RPAJastrow::LogValue updateBuffer ( ParticleSet P,
WFBufferType buf,
bool  fromscratch = false 
)
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 264 of file RPAJastrow.cpp.

References WaveFunctionComponent::log_value_, and RPAJastrow::Psi.

265 {
266  log_value_ = 0.0;
267  for (int i = 0; i < Psi.size(); i++)
268  log_value_ += Psi[i]->updateBuffer(P, buf, fromscratch);
269  return log_value_;
270 }
LogValue updateBuffer(ParticleSet &P, WFBufferType &buf, bool fromscratch=false) override
For particle-by-particle move.
Definition: RPAJastrow.cpp:264
UPtrVector< WaveFunctionComponent > Psi
A list of WaveFunctionComponent*.
Definition: RPAJastrow.h:121

Member Data Documentation

◆ DropLongRange

bool DropLongRange
private

Definition at line 98 of file RPAJastrow.h.

Referenced by RPAJastrow::buildOrbital(), and RPAJastrow::makeClone().

◆ DropShortRange

bool DropShortRange
private

Definition at line 99 of file RPAJastrow.h.

Referenced by RPAJastrow::buildOrbital(), and RPAJastrow::makeClone().

◆ ID_Rs

std::string ID_Rs
private

Definition at line 103 of file RPAJastrow.h.

Referenced by RPAJastrow::buildOrbital(), and RPAJastrow::put().

◆ IgnoreSpin

bool IgnoreSpin
private

Definition at line 97 of file RPAJastrow.h.

◆ Kc

◆ LongRangeRPA

kSpaceJastrow* LongRangeRPA
private

object to handle the long-range part

Definition at line 112 of file RPAJastrow.h.

Referenced by RPAJastrow::checkOutVariables(), RPAJastrow::extractOptimizableObjectRefs(), and RPAJastrow::makeLongRange().

◆ myGrid

GridType* myGrid
private

Definition at line 118 of file RPAJastrow.h.

◆ myHandler

std::unique_ptr<HandlerType> myHandler
private

◆ MyName

std::string MyName
private

Definition at line 105 of file RPAJastrow.h.

Referenced by RPAJastrow::buildOrbital(), and RPAJastrow::put().

◆ nfunc

FuncType* nfunc
private

numerical function owned by ShortRangeRPA

Definition at line 117 of file RPAJastrow.h.

Referenced by RPAJastrow::makeShortRange().

◆ Psi

◆ Rcut

RealType Rcut
private

Definition at line 102 of file RPAJastrow.h.

Referenced by RPAJastrow::makeClone(), and RPAJastrow::makeShortRange().

◆ rpafunc

std::string rpafunc
private

Definition at line 104 of file RPAJastrow.h.

Referenced by RPAJastrow::buildOrbital(), RPAJastrow::makeClone(), and RPAJastrow::put().

◆ Rs

RealType Rs
private

Definition at line 100 of file RPAJastrow.h.

Referenced by RPAJastrow::buildOrbital(), and RPAJastrow::put().

◆ ShortRangeRPA

WaveFunctionComponent* ShortRangeRPA
private

objects to handle the short-range part two-body Jastrow function

Definition at line 115 of file RPAJastrow.h.

Referenced by RPAJastrow::checkOutVariables(), RPAJastrow::extractOptimizableObjectRefs(), and RPAJastrow::makeShortRange().

◆ targetPtcl

ParticleSet& targetPtcl
private

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