![]() |
QMCPACK
|
TWFFastDerivWrapper is a wrapper class for TrialWavefunction that provides separate and low level access to the Jastrow and SPOSet objects. More...
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... | |
SPOSet * | getSPOSet (const IndexType sid) const |
IndexType | numOrbitals (const IndexType sid) const |
Private Attributes | |
std::vector< SPOSet * > | spos_ |
std::vector< IndexType > | groups_ |
std::vector< ValueMatrix > | psi_M_ |
std::vector< ValueMatrix > | psi_M_inv_ |
std::vector< WaveFunctionComponent * > | jastrow_list_ |
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.
using GradMatrix = SPOSet::GradMatrix |
Definition at line 34 of file TWFFastDerivWrapper.h.
using GradType = QMCTraits::GradType |
Definition at line 39 of file TWFFastDerivWrapper.h.
using GradVector = SPOSet::GradVector |
Definition at line 41 of file TWFFastDerivWrapper.h.
using HessMatrix = SPOSet::HessMatrix |
Definition at line 35 of file TWFFastDerivWrapper.h.
using IndexType = QMCTraits::IndexType |
Definition at line 36 of file TWFFastDerivWrapper.h.
using RealType = QMCTraits::RealType |
Definition at line 37 of file TWFFastDerivWrapper.h.
using ValueMatrix = SPOSet::ValueMatrix |
Definition at line 33 of file TWFFastDerivWrapper.h.
using ValueType = QMCTraits::ValueType |
Definition at line 38 of file TWFFastDerivWrapper.h.
using ValueVector = SPOSet::ValueVector |
Definition at line 40 of file TWFFastDerivWrapper.h.
|
default |
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.
[in] | P. | ParticleSet |
[in] | groupid. | ParticleSet groupid to which SPOSet corresponds. |
[in] | spo. | Pointer to SPOSet. |
Definition at line 30 of file TWFFastDerivWrapper.cpp.
References TWFFastDerivWrapper::groups_, and TWFFastDerivWrapper::spos_.
|
inline |
Definition at line 58 of file TWFFastDerivWrapper.h.
References TWFFastDerivWrapper::jastrow_list_.
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.
[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. |
Definition at line 248 of file TWFFastDerivWrapper.cpp.
References B().
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().
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')).)
[in] | P | electron particle set. |
[in] | iel | electron being moved. |
[in,out] | grad | grad_iel(J(r')). So iel electron gradient at new position. |
Definition at line 87 of file TWFFastDerivWrapper.cpp.
References qmcplusplus::convertToReal(), and TWFFastDerivWrapper::jastrow_list_.
Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution().
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.
[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. |
Definition at line 215 of file TWFFastDerivWrapper.cpp.
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().
TWFFastDerivWrapper::GradType evaluateJastrowGradSource | ( | ParticleSet & | P, |
ParticleSet & | source, | ||
const int | iat | ||
) | const |
Return ionic gradient of J(r).
[in] | P | electron particle set. |
[in] | source | ion particle set. |
[in] | iat | Ion to take derivative w.r.t. |
Definition at line 102 of file TWFFastDerivWrapper.cpp.
References TWFFastDerivWrapper::jastrow_list_.
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast(), NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution(), and BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().
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)).
[in] | P | electron particle set. |
[in] | source | ion particle set. |
[in] | iat | Ion to take derivative w.r.t. |
[in,out] | grad_grad | iat ion gradient of electron gradients of J(r). |
[in,out] | grad_lapl | iat ion gradient of electron laplacians of J(r). |
Definition at line 112 of file TWFFastDerivWrapper.cpp.
References TWFFastDerivWrapper::jastrow_list_, OHMMS_DIM, and TinyVector< T, D >::size().
TWFFastDerivWrapper::RealType evaluateJastrowRatio | ( | ParticleSet & | P, |
const int | iel | ||
) | const |
Given a proposed move of electron iel from r->r', computes exp(J')/exp(J)
[in] | P | electron particle set. |
[in] | iel | electron being moved. |
Definition at line 72 of file TWFFastDerivWrapper.cpp.
References qmcplusplus::convertToReal(), and TWFFastDerivWrapper::jastrow_list_.
Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixContribution().
TWFFastDerivWrapper::RealType evaluateJastrowVGL | ( | const ParticleSet & | P, |
ParticleSet::ParticleGradient & | G, | ||
ParticleSet::ParticleLaplacian & | L | ||
) | const |
Evaluates value, gradient, and laplacian of the total jastrow.
So of J in exp(J).
[in] | Particle | set. |
[in,out] | Electron | gradients. |
[in,out] | Electron | laplacians. |
Definition at line 57 of file TWFFastDerivWrapper.cpp.
References TWFFastDerivWrapper::jastrow_list_.
Referenced by BareKineticEnergy::evaluateOneBodyOpMatrix(), and BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().
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.
[in] | P | particle set. |
[in,out] | mvec | Slater matrix M_ij=phi_j(r_i) for each species group. |
[in,out] | gmat | electron gradient of slater matrix [G_ij]_a = d/dr_a,i phi_j(r_i). a=x,y,z |
[in,out] | lmat | electron laplacian of slater matrix [L_ij] = ^2_i phi_j(r_i). |
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().
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.
[in] | A | non-rectangular matrices of SPOSet derived quantities, evaluated on all orbitals and all particles. |
[in,out] | Aslice | square matrices consistent with a ground state occupation. |
Definition at line 309 of file TWFFastDerivWrapper.cpp.
References qmcplusplus::Units::distance::A.
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().
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.
[in] | P | particle set. |
[in] | source | ion particle set. |
[in] | iat | ion ID w.r.t. which to take derivative. |
[in,out] | dmvec | Slater 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] | dlmat | Slater 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. |
Definition at line 177 of file TWFFastDerivWrapper.cpp.
References ParticleSet::first(), TWFFastDerivWrapper::groups_, ParticleSet::last(), OHMMS_DIM, and TWFFastDerivWrapper::spos_.
Referenced by BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().
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.
[in] | P | particle set. |
[in] | source | ion 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().
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.
[in] | P | particle set. |
[in,out] | mmat | Output vector of slater matrices. Each vector entry corresponds to a different particle group. |
Definition at line 39 of file TWFFastDerivWrapper.cpp.
References ParticleSet::first(), TWFFastDerivWrapper::groups_, ParticleSet::last(), and TWFFastDerivWrapper::spos_.
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().
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.
[in] | P | particle set. |
[in] | iel | particle ID. |
[in,out] | val | Vector of phi_i(r_iel) for all i=0,Norbs. |
Definition at line 324 of file TWFFastDerivWrapper.cpp.
References ParticleSet::getGroupID(), TWFFastDerivWrapper::getTWFGroupIndex(), and TWFFastDerivWrapper::spos_.
Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixContribution().
Definition at line 84 of file TWFFastDerivWrapper.h.
References TWFFastDerivWrapper::spos_.
Referenced by NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution().
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.
[in] | gid. | ParticleSet group ID to look up. |
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().
void invertMatrices | ( | const std::vector< ValueMatrix > & | M, |
std::vector< ValueMatrix > & | Minv | ||
) |
Helper function that inverts all slater matrices in our species list.
[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. |
Definition at line 237 of file TWFFastDerivWrapper.cpp.
References qmcplusplus::invert_matrix().
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().
|
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().
Definition at line 85 of file TWFFastDerivWrapper.h.
References TWFFastDerivWrapper::spos_.
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast(), BareKineticEnergy::evaluateOneBodyOpMatrix(), NonLocalECPComponent::evaluateOneBodyOpMatrixContribution(), NonLocalECPComponent::evaluateOneBodyOpMatrixdRContribution(), and BareKineticEnergy::evaluateOneBodyOpMatrixForceDeriv().
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.
Definition at line 286 of file TWFFastDerivWrapper.cpp.
References qmcplusplus::Units::distance::A, and B().
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().
void wipeMatrices | ( | std::vector< ValueMatrix > & | A | ) |
Goes through a list of matrices and zeros them out.
Does this in place.
[in,out] | A. | The list of matrices to be zeroed out. After call, A is all zeros. |
Definition at line 278 of file TWFFastDerivWrapper.cpp.
References qmcplusplus::Units::distance::A.
Referenced by QMCHamiltonian::evaluateIonDerivsDeterministicFast().
|
private |
|
private |
|
private |
Definition at line 272 of file TWFFastDerivWrapper.h.
|
private |
Definition at line 273 of file TWFFastDerivWrapper.h.
|
private |
Definition at line 270 of file TWFFastDerivWrapper.h.
Referenced by TWFFastDerivWrapper::addGroup(), TWFFastDerivWrapper::getEGradELaplM(), TWFFastDerivWrapper::getIonGradIonGradELaplM(), TWFFastDerivWrapper::getIonGradM(), TWFFastDerivWrapper::getM(), TWFFastDerivWrapper::getRowM(), TWFFastDerivWrapper::getSPOSet(), TWFFastDerivWrapper::numGroups(), and TWFFastDerivWrapper::numOrbitals().