QMCPACK
TWFFastDerivWrapper Class Reference

TWFFastDerivWrapper is a wrapper class for TrialWavefunction that provides separate and low level access to the Jastrow and SPOSet objects. More...

+ Collaboration diagram for TWFFastDerivWrapper:

Public Types

using ValueMatrix = SPOSet::ValueMatrix
 
using GradMatrix = SPOSet::GradMatrix
 
using HessMatrix = SPOSet::HessMatrix
 
using IndexType = QMCTraits::IndexType
 
using RealType = QMCTraits::RealType
 
using ValueType = QMCTraits::ValueType
 
using GradType = QMCTraits::GradType
 
using ValueVector = SPOSet::ValueVector
 
using GradVector = SPOSet::GradVector
 

Public Member Functions

 TWFFastDerivWrapper ()=default
 
void addGroup (const ParticleSet &P, const IndexType groupid, SPOSet *spo)
 Add a particle group. More...
 
void addJastrow (WaveFunctionComponent *j)
 
IndexType getTWFGroupIndex (const IndexType gid) const
 Takes particle set groupID and returns the TWF internal index for it. More...
 
void getM (const ParticleSet &P, std::vector< ValueMatrix > &mmat) const
 Returns the non-rectangular slater matrices M_ij=phi_j(r_i) (N_particle x N_orbs) for each species group. More...
 
IndexType getRowM (const ParticleSet &P, const IndexType iel, ValueVector &val) const
 Returns value of all orbitals (relevant to given species/group) at a particular particle coordinate. More...
 
void getEGradELaplM (const ParticleSet &P, std::vector< ValueMatrix > &mvec, std::vector< GradMatrix > &gmat, std::vector< ValueMatrix > &lmat) const
 Returns value, gradient, and laplacian matrices for all orbitals and all particles, species by species. More...
 
void getIonGradM (const ParticleSet &P, const ParticleSet &source, const int iat, std::vector< std::vector< ValueMatrix >> &dmvec) const
 Returns x,y,z components of ion gradient of slater matrices. More...
 
void getIonGradIonGradELaplM (const ParticleSet &P, const ParticleSet &source, int iat, std::vector< std::vector< ValueMatrix >> &dmvec, std::vector< std::vector< GradMatrix >> &dgmat, std::vector< std::vector< ValueMatrix >> &dlmat) const
 Returns x,y,z components of ion gradient of slater matrices and their laplacians. More...
 
RealType evaluateJastrowVGL (const ParticleSet &P, ParticleSet::ParticleGradient &G, ParticleSet::ParticleLaplacian &L) const
 Evaluates value, gradient, and laplacian of the total jastrow. More...
 
RealType evaluateJastrowRatio (ParticleSet &P, const int iel) const
 Given a proposed move of electron iel from r->r', computes exp(J')/exp(J) More...
 
RealType calcJastrowRatioGrad (ParticleSet &P, const int iel, GradType &grad) const
 Given a proposed move of electron iel from r->r', computes exp(J(r'))/exp(J(r)) and grad_iel(J(r')).) More...
 
GradType evaluateJastrowGradSource (ParticleSet &P, ParticleSet &source, const int iat) const
 Return ionic gradient of J(r). More...
 
GradType evaluateJastrowGradSource (ParticleSet &P, ParticleSet &source, const int iat, TinyVector< ParticleSet::ParticleGradient, OHMMS_DIM > &grad_grad, TinyVector< ParticleSet::ParticleLaplacian, OHMMS_DIM > &lapl_grad) const
 Return ionic gradients of J(r), grad_iel(J(r)), and lapl_iel(J(r)). More...
 
void getGSMatrices (const std::vector< ValueMatrix > &A, std::vector< ValueMatrix > &Aslice) const
 Takes sub matrices of full SPOSet quantities (evaluated on all particles and all orbitals), consistent with ground state occupations. More...
 
ValueType computeGSDerivative (const std::vector< ValueMatrix > &Minv, const std::vector< ValueMatrix > &X, const std::vector< ValueMatrix > &dM, const std::vector< ValueMatrix > &dB) const
 Calculates derivative of observable via Tr[M^{-1} dB - X * dM ]. More...
 
void invertMatrices (const std::vector< ValueMatrix > &M, std::vector< ValueMatrix > &Minv)
 Helper function that inverts all slater matrices in our species list. More...
 
void buildX (const std::vector< ValueMatrix > &Minv, const std::vector< ValueMatrix > &B, std::vector< ValueMatrix > &X)
 Build the auxiliary X=M^-1 B M^-1 matrix. More...
 
void wipeMatrices (std::vector< ValueMatrix > &A)
 Goes through a list of matrices and zeros them out. More...
 
ValueType trAB (const std::vector< ValueMatrix > &A, const std::vector< ValueMatrix > &B)
 Returns Tr(A*B). More...
 
IndexType numGroups () const
 These are query functions to the internal state of various groups and SPOSet info. More...
 
SPOSetgetSPOSet (const IndexType sid) const
 
IndexType numOrbitals (const IndexType sid) const
 

Private Attributes

std::vector< SPOSet * > spos_
 
std::vector< IndexTypegroups_
 
std::vector< ValueMatrixpsi_M_
 
std::vector< ValueMatrixpsi_M_inv_
 
std::vector< WaveFunctionComponent * > jastrow_list_
 

Detailed Description

TWFFastDerivWrapper is a wrapper class for TrialWavefunction that provides separate and low level access to the Jastrow and SPOSet objects.

This is so that observables can be recast in matrix form and their derivatives taken efficiently. Currently this is hard coded for ground state slater jastrow wave functions, but generalization to arbitrary occupations and multideterminants are straightforward and will come online as the code is tested and validated.

Please see : J. Chem. Phys. 144, 194105 (2016) https://doi.org/10.1063/1.4948778 for implementation details and formalism.

Definition at line 30 of file TWFFastDerivWrapper.h.

Member Typedef Documentation

◆ GradMatrix

Definition at line 34 of file TWFFastDerivWrapper.h.

◆ GradType

Definition at line 39 of file TWFFastDerivWrapper.h.

◆ GradVector

Definition at line 41 of file TWFFastDerivWrapper.h.

◆ HessMatrix

Definition at line 35 of file TWFFastDerivWrapper.h.

◆ IndexType

Definition at line 36 of file TWFFastDerivWrapper.h.

◆ RealType

Definition at line 37 of file TWFFastDerivWrapper.h.

◆ ValueMatrix

Definition at line 33 of file TWFFastDerivWrapper.h.

◆ ValueType

Definition at line 38 of file TWFFastDerivWrapper.h.

◆ ValueVector

Definition at line 40 of file TWFFastDerivWrapper.h.

Constructor & Destructor Documentation

◆ TWFFastDerivWrapper()

TWFFastDerivWrapper ( )
default

Member Function Documentation

◆ addGroup()

void addGroup ( const ParticleSet P,
const IndexType  groupid,
SPOSet spo 
)

Add a particle group.

Here, a "group" corresponds to a subset of particles which are antisymmetric with respect to eachother. TWFFastDerivWrapper ensures that there is a binding between the groupid in ParticleSet and the SPOSet associated with that particle group. This function stores the ParticleSet groupid and SPOSet in a vector for lookup and communication with QMCPACK conventions, but is agnostic to the order of group registration or evaluation.

Parameters
[in]P.ParticleSet
[in]groupid.ParticleSet groupid to which SPOSet corresponds.
[in]spo.Pointer to SPOSet.
Returns
void.

Definition at line 30 of file TWFFastDerivWrapper.cpp.

References TWFFastDerivWrapper::groups_, and TWFFastDerivWrapper::spos_.

31 {
32  if (std::find(groups_.begin(), groups_.end(), gid) == groups_.end())
33  {
34  groups_.push_back(gid);
35  spos_.push_back(spo);
36  }
37 }
std::vector< IndexType > groups_

◆ addJastrow()

void addJastrow ( WaveFunctionComponent j)
inline

Definition at line 58 of file TWFFastDerivWrapper.h.

References TWFFastDerivWrapper::jastrow_list_.

58 { jastrow_list_.push_back(j); };
std::vector< WaveFunctionComponent * > jastrow_list_

◆ buildX()

void buildX ( const std::vector< ValueMatrix > &  Minv,
const std::vector< ValueMatrix > &  B,
std::vector< ValueMatrix > &  X 
)

Build the auxiliary X=M^-1 B M^-1 matrix.

Parameters
[in]Minv.List of slater matrix inverses M^-1 for a given occupation.
[in]B.Observable auxiliary matrix for a given occupation.
[in,out]X.M^-1*B*M^-1 is stored in this list of matrices.
Returns
Void.

Definition at line 248 of file TWFFastDerivWrapper.cpp.

References B().

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

251 {
252  IndexType nspecies = Minv.size();
253 
254  for (IndexType id = 0; id < nspecies; id++)
255  {
256  int ptclnum = Minv[id].rows();
257  assert(Minv[id].rows() == Minv[id].cols());
258  ValueMatrix tmpmat;
259  X[id].resize(ptclnum, ptclnum);
260  tmpmat.resize(ptclnum, ptclnum);
261  //(B*A^-1)
262  for (int i = 0; i < ptclnum; i++)
263  for (int j = 0; j < ptclnum; j++)
264  for (int k = 0; k < ptclnum; k++)
265  {
266  tmpmat[i][j] += B[id][i][k] * Minv[id][k][j];
267  }
268  //A^{-1}*B*A^{-1}
269  for (int i = 0; i < ptclnum; i++)
270  for (int j = 0; j < ptclnum; j++)
271  for (int k = 0; k < ptclnum; k++)
272  {
273  X[id][i][j] += Minv[id][i][k] * tmpmat[k][j];
274  }
275  }
276 }
double B(double x, int k, int i, const std::vector< double > &t)

◆ calcJastrowRatioGrad()

TWFFastDerivWrapper::RealType calcJastrowRatioGrad ( ParticleSet P,
const int  iel,
GradType grad 
) const

Given a proposed move of electron iel from r->r', computes exp(J(r'))/exp(J(r)) and grad_iel(J(r')).)

Parameters
[in]Pelectron particle set.
[in]ielelectron being moved.
[in,out]gradgrad_iel(J(r')). So iel electron gradient at new position.
Returns
Jastrow ratio.

Definition at line 87 of file TWFFastDerivWrapper.cpp.

References qmcplusplus::convertToReal(), and TWFFastDerivWrapper::jastrow_list_.

Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution().

90 {
91  int iel_(iel);
93  for (int i = 0; i < jastrow_list_.size(); ++i)
94  {
95  r *= jastrow_list_[i]->ratioGrad(P, iel_, grad);
96  }
97  RealType ratio_return(1.0);
98  convertToReal(r, ratio_return);
99  return ratio_return;
100 }
std::vector< WaveFunctionComponent * > jastrow_list_
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

◆ computeGSDerivative()

TWFFastDerivWrapper::ValueType computeGSDerivative ( const std::vector< ValueMatrix > &  Minv,
const std::vector< ValueMatrix > &  X,
const std::vector< ValueMatrix > &  dM,
const std::vector< ValueMatrix > &  dB 
) const

Calculates derivative of observable via Tr[M^{-1} dB - X * dM ].

Consistent with ground state occupation.

Parameters
[in]Minv.inverse of slater matrices for ground state occupations.
[in]X.X=M^-1 B M^-1. B observable matrix, and is consistent with ground state occupation.
[in]dM.Target derivative of M, and is consistent with ground state occupation.
[in]dB.Target derivative of B, and is consistent with ground state occupation.
Returns
Derivative of O psi/psi = Tr[M^{-1} dB - X * dM ]

Definition at line 215 of file TWFFastDerivWrapper.cpp.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

219 {
220  IndexType nspecies = Minv.size();
221  ValueType dval = 0.0;
222  for (int id = 0; id < nspecies; id++)
223  {
224  int ptclnum = Minv[id].rows();
225  ValueType dval_id = 0.0;
226  for (int i = 0; i < ptclnum; i++)
227  for (int j = 0; j < ptclnum; j++)
228  {
229  //Tr[M^{-1} dB - X * dM ]
230  dval_id += Minv[id][i][j] * dB[id][j][i] - X[id][i][j] * dM[id][j][i];
231  }
232  dval += dval_id;
233  }
234  return dval;
235 }

◆ evaluateJastrowGradSource() [1/2]

TWFFastDerivWrapper::GradType evaluateJastrowGradSource ( ParticleSet P,
ParticleSet source,
const int  iat 
) const

Return ionic gradient of J(r).

Parameters
[in]Pelectron particle set.
[in]sourceion particle set.
[in]iatIon to take derivative w.r.t.
Returns
Gradient of J(r) w.r.t. ion iat.

Definition at line 102 of file TWFFastDerivWrapper.cpp.

References TWFFastDerivWrapper::jastrow_list_.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast(), NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution(), and BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().

105 {
106  GradType grad_iat = GradType();
107  for (int i = 0; i < jastrow_list_.size(); ++i)
108  grad_iat += jastrow_list_[i]->evalGradSource(P, source, iat);
109  return grad_iat;
110 }
std::vector< WaveFunctionComponent * > jastrow_list_

◆ evaluateJastrowGradSource() [2/2]

TWFFastDerivWrapper::GradType evaluateJastrowGradSource ( ParticleSet P,
ParticleSet source,
const int  iat,
TinyVector< ParticleSet::ParticleGradient, OHMMS_DIM > &  grad_grad,
TinyVector< ParticleSet::ParticleLaplacian, OHMMS_DIM > &  lapl_grad 
) const

Return ionic gradients of J(r), grad_iel(J(r)), and lapl_iel(J(r)).

Parameters
[in]Pelectron particle set.
[in]sourceion particle set.
[in]iatIon to take derivative w.r.t.
[in,out]grad_gradiat ion gradient of electron gradients of J(r).
[in,out]grad_lapliat ion gradient of electron laplacians of J(r).
Returns
Gradient of J(r) w.r.t. ion iat.

Definition at line 112 of file TWFFastDerivWrapper.cpp.

References TWFFastDerivWrapper::jastrow_list_, OHMMS_DIM, and TinyVector< T, D >::size().

118 {
119  GradType grad_iat = GradType();
120  for (int dim = 0; dim < OHMMS_DIM; dim++)
121  for (int i = 0; i < grad_grad[0].size(); i++)
122  {
123  grad_grad[dim][i] = GradType();
124  lapl_grad[dim][i] = 0.0;
125  }
126  for (int i = 0; i < jastrow_list_.size(); ++i)
127  grad_iat += jastrow_list_[i]->evalGradSource(P, source, iat, grad_grad, lapl_grad);
128  return grad_iat;
129 }
std::vector< WaveFunctionComponent * > jastrow_list_
#define OHMMS_DIM
Definition: config.h:64

◆ evaluateJastrowRatio()

TWFFastDerivWrapper::RealType evaluateJastrowRatio ( ParticleSet P,
const int  iel 
) const

Given a proposed move of electron iel from r->r', computes exp(J')/exp(J)

Parameters
[in]Pelectron particle set.
[in]ielelectron being moved.
Returns
Jastrow ratio.

Definition at line 72 of file TWFFastDerivWrapper.cpp.

References qmcplusplus::convertToReal(), and TWFFastDerivWrapper::jastrow_list_.

Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixContribution().

73 {
74  //legacy calls are hit and miss with const. Remove const for index.
75  int iel_(iel);
77  for (int i = 0; i < jastrow_list_.size(); ++i)
78  {
79  r *= jastrow_list_[i]->ratio(P, iel_);
80  }
81 
82  RealType ratio_return(1.0);
83  convertToReal(r, ratio_return);
84  return ratio_return;
85 }
std::vector< WaveFunctionComponent * > jastrow_list_
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

◆ evaluateJastrowVGL()

Evaluates value, gradient, and laplacian of the total jastrow.

So of J in exp(J).

Parameters
[in]Particleset.
[in,out]Electrongradients.
[in,out]Electronlaplacians.
Returns
Jastrow value

Definition at line 57 of file TWFFastDerivWrapper.cpp.

References TWFFastDerivWrapper::jastrow_list_.

Referenced by BareKineticEnergy::evaluateOneBodyOpMatrix(), and BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().

60 {
62  G = 0.0;
63  L = 0.0;
64  for (int i = 0; i < jastrow_list_.size(); ++i)
65  {
66  logpsi += jastrow_list_[i]->evaluateLog(P, G, L);
67  }
68  RealType rval = std::real(logpsi);
69  return rval;
70 }
std::vector< WaveFunctionComponent * > jastrow_list_
QMCTraits::RealType real
std::complex< QTFull::RealType > LogValue
QMCTraits::RealType RealType

◆ getEGradELaplM()

void getEGradELaplM ( const ParticleSet P,
std::vector< ValueMatrix > &  mvec,
std::vector< GradMatrix > &  gmat,
std::vector< ValueMatrix > &  lmat 
) const

Returns value, gradient, and laplacian matrices for all orbitals and all particles, species by species.

Parameters
[in]Pparticle set.
[in,out]mvecSlater matrix M_ij=phi_j(r_i) for each species group.
[in,out]gmatelectron gradient of slater matrix [G_ij]_a = d/dr_a,i phi_j(r_i). a=x,y,z
[in,out]lmatelectron laplacian of slater matrix [L_ij] = ^2_i phi_j(r_i).
Returns
Void

Definition at line 131 of file TWFFastDerivWrapper.cpp.

References ParticleSet::first(), TWFFastDerivWrapper::groups_, ParticleSet::last(), and TWFFastDerivWrapper::spos_.

Referenced by BareKineticEnergy::evaluateOneBodyOpMatrix(), and BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().

135 {
136  IndexType ngroups = mvec.size();
137  for (IndexType i = 0; i < ngroups; i++)
138  {
139  const IndexType gid = groups_[i];
140  const IndexType first = P.first(i);
141  const IndexType last = P.last(i);
142  const IndexType nptcls = last - first;
143  const IndexType norbs = spos_[i]->getOrbitalSetSize();
144  spos_[i]->evaluate_notranspose(P, first, last, mvec[i], gmat[i], lmat[i]);
145  }
146 }
std::vector< IndexType > groups_

◆ getGSMatrices()

void getGSMatrices ( const std::vector< ValueMatrix > &  A,
std::vector< ValueMatrix > &  Aslice 
) const

Takes sub matrices of full SPOSet quantities (evaluated on all particles and all orbitals), consistent with ground state occupations.

Parameters
[in]Anon-rectangular matrices of SPOSet derived quantities, evaluated on all orbitals and all particles.
[in,out]Aslicesquare matrices consistent with a ground state occupation.
Returns
Void

Definition at line 309 of file TWFFastDerivWrapper.cpp.

References qmcplusplus::Units::distance::A.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

310 {
311  IndexType nspecies = A.size();
312  Aslice.resize(nspecies);
313  for (IndexType id = 0; id < nspecies; id++)
314  {
315  IndexType ptclnum = A[id].rows();
316  Aslice[id].resize(ptclnum, ptclnum);
317  for (IndexType i = 0; i < ptclnum; i++)
318  for (IndexType j = 0; j < ptclnum; j++)
319  Aslice[id][i][j] = A[id][i][j];
320  }
321 }

◆ getIonGradIonGradELaplM()

void getIonGradIonGradELaplM ( const ParticleSet P,
const ParticleSet source,
int  iat,
std::vector< std::vector< ValueMatrix >> &  dmvec,
std::vector< std::vector< GradMatrix >> &  dgmat,
std::vector< std::vector< ValueMatrix >> &  dlmat 
) const

Returns x,y,z components of ion gradient of slater matrices and their laplacians.

Parameters
[in]Pparticle set.
[in]sourceion particle set.
[in]iation ID w.r.t. which to take derivative.
[in,out]dmvecSlater matrix d/dR_{iat,a} M_ij=d/dR_{iat,a} phi_j(r_i) for each species group.
First index is a=x,y,z.
[in,out]dlmatSlater matrix d/dR_{iat,a} L_ij=d/dR_{iat,a} ^2_i phi_j(r_i) for each species group.
First index is a=x,y,z.
Returns
Void

Definition at line 177 of file TWFFastDerivWrapper.cpp.

References ParticleSet::first(), TWFFastDerivWrapper::groups_, ParticleSet::last(), OHMMS_DIM, and TWFFastDerivWrapper::spos_.

Referenced by BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().

183 {
184  IndexType ngroups = dmvec[0].size();
185  for (IndexType i = 0; i < ngroups; i++)
186  {
187  const IndexType gid = groups_[i];
188  const IndexType first = P.first(i);
189  const IndexType last = P.last(i);
190  const IndexType nptcls = last - first;
191  const IndexType norbs = spos_[i]->getOrbitalSetSize();
192 
193  GradMatrix grad_phi;
194  HessMatrix grad_grad_phi;
195  GradMatrix grad_lapl_phi;
196 
197  grad_phi.resize(nptcls, norbs);
198  grad_grad_phi.resize(nptcls, norbs);
199  grad_lapl_phi.resize(nptcls, norbs);
200 
201  spos_[i]->evaluateGradSource(P, first, last, source, iat, grad_phi, grad_grad_phi, grad_lapl_phi);
202 
203  for (IndexType idim = 0; idim < OHMMS_DIM; idim++)
204  for (IndexType iptcl = 0; iptcl < nptcls; iptcl++)
205  for (IndexType iorb = 0; iorb < norbs; iorb++)
206  {
207  dmvec[idim][i][iptcl][iorb] += grad_phi[iptcl][iorb][idim];
208  dlmat[idim][i][iptcl][iorb] += grad_lapl_phi[iptcl][iorb][idim];
209  for (IndexType ielec = 0; ielec < OHMMS_DIM; ielec++)
210  dgmat[idim][i][iptcl][iorb][ielec] += grad_grad_phi[iptcl][iorb](idim, ielec);
211  }
212  }
213 }
#define OHMMS_DIM
Definition: config.h:64
std::vector< IndexType > groups_

◆ getIonGradM()

void getIonGradM ( const ParticleSet P,
const ParticleSet source,
const int  iat,
std::vector< std::vector< ValueMatrix >> &  dmvec 
) const

Returns x,y,z components of ion gradient of slater matrices.

Parameters
[in]Pparticle set.
[in]sourceion particle set.
[in]

Definition at line 148 of file TWFFastDerivWrapper.cpp.

References ParticleSet::first(), TWFFastDerivWrapper::groups_, ParticleSet::last(), OHMMS_DIM, and TWFFastDerivWrapper::spos_.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

152 {
153  IndexType ngroups = dmvec[0].size();
154  for (IndexType i = 0; i < ngroups; i++)
155  {
156  const IndexType gid = groups_[i];
157  const IndexType first = P.first(i);
158  const IndexType last = P.last(i);
159  const IndexType nptcls = last - first;
160  const IndexType norbs = spos_[i]->getOrbitalSetSize();
161 
162  GradMatrix grad_phi;
163 
164  grad_phi.resize(nptcls, norbs);
165 
166  spos_[i]->evaluateGradSource(P, first, last, source, iat, grad_phi);
167 
168  for (IndexType idim = 0; idim < OHMMS_DIM; idim++)
169  for (IndexType iptcl = 0; iptcl < nptcls; iptcl++)
170  for (IndexType iorb = 0; iorb < norbs; iorb++)
171  {
172  dmvec[idim][i][iptcl][iorb] += grad_phi[iptcl][iorb][idim];
173  }
174  }
175 }
#define OHMMS_DIM
Definition: config.h:64
std::vector< IndexType > groups_

◆ getM()

void getM ( const ParticleSet P,
std::vector< ValueMatrix > &  mmat 
) const

Returns the non-rectangular slater matrices M_ij=phi_j(r_i) (N_particle x N_orbs) for each species group.

Parameters
[in]Pparticle set.
[in,out]mmatOutput vector of slater matrices. Each vector entry corresponds to a different particle group.
Returns
Void

Definition at line 39 of file TWFFastDerivWrapper.cpp.

References ParticleSet::first(), TWFFastDerivWrapper::groups_, ParticleSet::last(), and TWFFastDerivWrapper::spos_.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

40 {
41  IndexType ngroups = spos_.size();
42  for (IndexType i = 0; i < ngroups; i++)
43  {
44  const IndexType gid = groups_[i];
45  const IndexType first = P.first(i);
46  const IndexType last = P.last(i);
47  const IndexType nptcls = last - first;
48  const IndexType norbs = spos_[i]->getOrbitalSetSize();
49  GradMatrix tmpgmat;
50  ValueMatrix tmplmat;
51  tmpgmat.resize(nptcls, norbs);
52  tmplmat.resize(nptcls, norbs);
53  spos_[i]->evaluate_notranspose(P, first, last, mvec[i], tmpgmat, tmplmat);
54  }
55 }
std::vector< IndexType > groups_

◆ getRowM()

TWFFastDerivWrapper::IndexType getRowM ( const ParticleSet P,
const IndexType  iel,
ValueVector val 
) const

Returns value of all orbitals (relevant to given species/group) at a particular particle coordinate.

Parameters
[in]Pparticle set.
[in]ielparticle ID.
[in,out]valVector of phi_i(r_iel) for all i=0,Norbs.
Returns
Void

Definition at line 324 of file TWFFastDerivWrapper.cpp.

References ParticleSet::getGroupID(), TWFFastDerivWrapper::getTWFGroupIndex(), and TWFFastDerivWrapper::spos_.

Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixContribution().

327 {
328  IndexType gid = P.getGroupID(iel);
329  IndexType sid = getTWFGroupIndex(gid);
330 
331  GradVector tempg;
332  ValueVector templ;
333 
334  IndexType norbs = spos_[sid]->getOrbitalSetSize();
335 
336  tempg.resize(norbs);
337  templ.resize(norbs);
338 
339  spos_[sid]->evaluateVGL(P, iel, val, tempg, templ);
340 
341  return sid;
342 }
IndexType getTWFGroupIndex(const IndexType gid) const
Takes particle set groupID and returns the TWF internal index for it.

◆ getSPOSet()

SPOSet* getSPOSet ( const IndexType  sid) const
inline

Definition at line 84 of file TWFFastDerivWrapper.h.

References TWFFastDerivWrapper::spos_.

Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution().

84 { return spos_[sid]; };

◆ getTWFGroupIndex()

TWFFastDerivWrapper::IndexType getTWFGroupIndex ( const IndexType  gid) const

Takes particle set groupID and returns the TWF internal index for it.

ParticleSet groups can be registered in whichever order. However, the internal indexing of TWFFastDerivWrapper keeps track on a first-come, first serve basis. That is, if I register particle groups 3, 1, and 0 in that order, then the internal indexing goes like 0->3, 1->1, 2->0. Thus, this function looks up where in the internal indexing scheme ParticleSet gid is located. This is necessary, since the matrix list data structures follow the internal TWF indexing.

Parameters
[in]gid.ParticleSet group ID to look up.
Returns
The corresponding internal groupID.

Definition at line 18 of file TWFFastDerivWrapper.cpp.

References TWFFastDerivWrapper::groups_.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast(), BareKineticEnergy::evaluateOneBodyOpMatrix(), NonLocalECPComponent::evaluateOneBodyOpMatrixContribution(), NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution(), BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv(), and TWFFastDerivWrapper::getRowM().

19 {
20  IndexType return_group_index(-1);
21  for (IndexType i = 0; i < groups_.size(); i++)
22  if (gid == groups_[i])
23  return_group_index = i;
24 
25  assert(return_group_index != -1);
26 
27  return return_group_index;
28 }
std::vector< IndexType > groups_

◆ invertMatrices()

void invertMatrices ( const std::vector< ValueMatrix > &  M,
std::vector< ValueMatrix > &  Minv 
)

Helper function that inverts all slater matrices in our species list.

Parameters
[in]M.List of slater matrices for each species. These are square and consistent with some occupation.
[in,out]Minv.The species by species list of inverted matrices from M.
Returns
Void.

Definition at line 237 of file TWFFastDerivWrapper.cpp.

References qmcplusplus::invert_matrix().

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

238 {
239  IndexType nspecies = M.size();
240  for (IndexType id = 0; id < nspecies; id++)
241  {
242  assert(M[id].cols() == M[id].rows());
243  Minv[id] = M[id];
244  invert_matrix(Minv[id]);
245  }
246 }
MatrixA::value_type invert_matrix(MatrixA &M, bool getdet=true)
invert a matrix

◆ numGroups()

IndexType numGroups ( ) const
inline

These are query functions to the internal state of various groups and SPOSet info.

All group indices will refer to the internal group indices. Convert from ParticleSet to TWF group first!

Source of truth for orbital sizes will be the individual SPOSets. Particle group sizes will be ParticleSet in conjunction with groupID maps.

Definition at line 83 of file TWFFastDerivWrapper.h.

References TWFFastDerivWrapper::spos_.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

83 { return spos_.size(); };

◆ numOrbitals()

◆ trAB()

TWFFastDerivWrapper::ValueType trAB ( const std::vector< ValueMatrix > &  A,
const std::vector< ValueMatrix > &  B 
)

Returns Tr(A*B).

Actually, we assume a block diagonal structure, so this is really Sum_i Tr(A_i*B_i), where i is the species index.

Parameters
[in]A.Vector of matrices A.
[in]B.Vector of matrices B.
Returns
Value of Sum_i Tr(A_i*B_i).

Definition at line 286 of file TWFFastDerivWrapper.cpp.

References qmcplusplus::Units::distance::A, and B().

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

288 {
289  IndexType nspecies = A.size();
290  assert(A.size() == B.size());
291  ValueType val = 0.0;
292  //Now to compute the kinetic energy
293  for (IndexType id = 0; id < nspecies; id++)
294  {
295  int ptclnum = A[id].rows();
296  ValueType val_id = 0.0;
297  assert(A[id].cols() == B[id].rows() && A[id].rows() == B[id].cols());
298  for (int i = 0; i < A[id].rows(); i++)
299  for (int j = 0; j < A[id].cols(); j++)
300  {
301  val_id += A[id][i][j] * B[id][j][i];
302  }
303  val += val_id;
304  }
305 
306  return val;
307 }
double B(double x, int k, int i, const std::vector< double > &t)

◆ wipeMatrices()

void wipeMatrices ( std::vector< ValueMatrix > &  A)

Goes through a list of matrices and zeros them out.

Does this in place.

Parameters
[in,out]A.The list of matrices to be zeroed out. After call, A is all zeros.
Returns
Void.

Definition at line 278 of file TWFFastDerivWrapper.cpp.

References qmcplusplus::Units::distance::A.

Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().

279 {
280  for (IndexType id = 0; id < A.size(); id++)
281  {
282  A[id] = 0.0;
283  }
284 }

Member Data Documentation

◆ groups_

◆ jastrow_list_

◆ psi_M_

std::vector<ValueMatrix> psi_M_
private

Definition at line 272 of file TWFFastDerivWrapper.h.

◆ psi_M_inv_

std::vector<ValueMatrix> psi_M_inv_
private

Definition at line 273 of file TWFFastDerivWrapper.h.

◆ spos_


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