QMCPACK
SplineR2R< ST > Class Template Reference

class to match ST real spline with BsplineSet::ValueType (real) SPOs More...

+ Inheritance diagram for SplineR2R< ST >:
+ Collaboration diagram for SplineR2R< ST >:

Public Types

using SplineType = typename bspline_traits< ST, 3 >::SplineType
 
using BCType = typename bspline_traits< ST, 3 >::BCType
 
using DataType = ST
 
using PointType = TinyVector< ST, 3 >
 
using SingleSplineType = UBspline_3d_d
 
using TT = typename BsplineSet::ValueType
 
using vContainer_type = Vector< ST, aligned_allocator< ST > >
 
using gContainer_type = VectorSoaContainer< ST, 3 >
 
using hContainer_type = VectorSoaContainer< ST, 6 >
 
using ghContainer_type = VectorSoaContainer< ST, 10 >
 
- Public Types inherited from SPOSet
using ValueVector = OrbitalSetTraits< ValueType >::ValueVector
 
using ValueMatrix = OrbitalSetTraits< ValueType >::ValueMatrix
 
using GradVector = OrbitalSetTraits< ValueType >::GradVector
 
using GradMatrix = OrbitalSetTraits< ValueType >::GradMatrix
 
using HessVector = OrbitalSetTraits< ValueType >::HessVector
 
using HessMatrix = OrbitalSetTraits< ValueType >::HessMatrix
 
using GGGVector = OrbitalSetTraits< ValueType >::GradHessVector
 
using GGGMatrix = OrbitalSetTraits< ValueType >::GradHessMatrix
 
using SPOMap = std::map< std::string, const std::unique_ptr< const SPOSet > >
 
using OffloadMWVGLArray = Array< ValueType, 3, OffloadPinnedAllocator< ValueType > >
 
using OffloadMWVArray = Array< ValueType, 2, OffloadPinnedAllocator< ValueType > >
 
template<typename DT >
using OffloadMatrix = Matrix< DT, OffloadPinnedAllocator< DT > >
 
- 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 Member Functions

 SplineR2R (const std::string &my_name)
 
 SplineR2R (const SplineR2R &in)
 
virtual std::string getClassName () const override
 return class name More...
 
virtual std::string getKeyword () const override
 
bool isComplex () const override
 
bool isRotationSupported () const override
 return true if this SPOSet can be wrappered by RotatedSPO More...
 
std::unique_ptr< SPOSetmakeClone () const override
 make a clone of itself every derived class must implement this to have threading working correctly. More...
 
void storeParamsBeforeRotation () override
 Store an original copy of the spline coefficients for orbital rotation. More...
 
void applyRotation (const ValueMatrix &rot_mat, bool use_stored_copy) override
 apply rotation to all the orbitals More...
 
void resizeStorage (size_t n, size_t nvals)
 
void bcast_tables (Communicate *comm)
 
void gather_tables (Communicate *comm)
 
template<typename GT , typename BCT >
void create_spline (GT &xyz_g, BCT &xyz_bc)
 
void flush_zero ()
 
void set_spline (SingleSplineType *spline_r, SingleSplineType *spline_i, int twist, int ispline, int level)
 
bool read_splines (hdf_archive &h5f)
 
bool write_splines (hdf_archive &h5f)
 
int convertPos (const PointType &r, PointType &ru)
 convert position in PrimLattice unit and return sign More...
 
void assign_v (int bc_sign, const vContainer_type &myV, ValueVector &psi, int first, int last) const
 
void evaluateValue (const ParticleSet &P, const int iat, ValueVector &psi) override
 evaluate the values of this single-particle orbital set More...
 
void evaluateDetRatios (const VirtualParticleSet &VP, ValueVector &psi, const ValueVector &psiinv, std::vector< TT > &ratios) override
 
void assign_vgl (int bc_sign, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi, int first, int last) const
 
void assign_vgl_from_l (int bc_sign, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi)
 assign_vgl_from_l can be used when myL is precomputed and myV,myG,myL in cartesian More...
 
void evaluateVGL (const ParticleSet &P, const int iat, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi) override
 evaluate the values, gradients and laplacians of this single-particle orbital set More...
 
void assign_vgh (int bc_sign, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi, int first, int last) const
 
void evaluateVGH (const ParticleSet &P, const int iat, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi) override
 evaluate the values, gradients and hessians of this single-particle orbital set More...
 
void assign_vghgh (int bc_sign, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi, GGGVector &grad_grad_grad_psi, int first=0, int last=-1) const
 
void evaluateVGHGH (const ParticleSet &P, const int iat, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi, GGGVector &grad_grad_grad_psi) override
 evaluate the values, gradients, hessians, and grad hessians of this single-particle orbital set More...
 
- Public Member Functions inherited from BsplineSet
 BsplineSet (const std::string &my_name)
 
auto & getHalfG () const
 
void init_base (int n)
 
int remap_kpoints ()
 remap kpoints to group general kpoints & special kpoints More...
 
void setOrbitalSetSize (int norbs) override
 set the OrbitalSetSize More...
 
void evaluate_notranspose (const ParticleSet &P, int first, int last, ValueMatrix &logdet, GradMatrix &dlogdet, ValueMatrix &d2logdet) override
 evaluate the values, gradients and laplacians of this single-particle orbital for [first,last) particles More...
 
void mw_evaluate_notranspose (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int first, int last, const RefVector< ValueMatrix > &logdet_list, const RefVector< GradMatrix > &dlogdet_list, const RefVector< ValueMatrix > &d2logdet_list) const override
 
void evaluate_notranspose (const ParticleSet &P, int first, int last, ValueMatrix &logdet, GradMatrix &dlogdet, HessMatrix &grad_grad_logdet) override
 evaluate the values, gradients and hessians of this single-particle orbital for [first,last) particles More...
 
void evaluate_notranspose (const ParticleSet &P, int first, int last, ValueMatrix &logdet, GradMatrix &dlogdet, HessMatrix &grad_grad_logdet, GGGMatrix &grad_grad_grad_logdet) override
 evaluate the values, gradients, hessians and third derivatives of this single-particle orbital for [first,last) particles More...
 
void evaluateGradSource (const ParticleSet &P, int first, int last, const ParticleSet &source, int iat_src, GradMatrix &gradphi) override
 evaluate the gradients of this single-particle orbital for [first,last) target particles with respect to the given source particle More...
 
void evaluateGradSource (const ParticleSet &P, int first, int last, const ParticleSet &source, int iat_src, GradMatrix &grad_phi, HessMatrix &grad_grad_phi, GradMatrix &grad_lapl_phi) override
 evaluate the gradients of values, gradients, laplacians of this single-particle orbital for [first,last) target particles with respect to the given source particle More...
 
virtual void evaluateDetRatios (const VirtualParticleSet &VP, ValueVector &psi, const ValueVector &psiinv, std::vector< ValueType > &ratios)
 evaluate determinant ratios for virtual moves, e.g., sphere move for nonlocalPP More...
 
virtual void evaluateValue (const ParticleSet &P, int iat, ValueVector &psi)=0
 evaluate the values of this single-particle orbital set More...
 
virtual void evaluateVGH (const ParticleSet &P, int iat, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi)
 evaluate the values, gradients and hessians of this single-particle orbital set More...
 
virtual void evaluateVGHGH (const ParticleSet &P, int iat, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi, GGGVector &grad_grad_grad_psi)
 evaluate the values, gradients, hessians, and grad hessians of this single-particle orbital set More...
 
virtual void evaluateVGL (const ParticleSet &P, int iat, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi)=0
 evaluate the values, gradients and laplacians of this single-particle orbital set More...
 
virtual void finalizeConstruction ()
 finalize the construction of SPOSet More...
 
virtual void mw_evaluateDetRatios (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< const VirtualParticleSet > &vp_list, const RefVector< ValueVector > &psi_list, const std::vector< const ValueType * > &invRow_ptr_list, std::vector< std::vector< ValueType >> &ratios_list) const
 evaluate determinant ratios for virtual moves, e.g., sphere move for nonlocalPP, of multiple walkers More...
 
virtual void mw_evaluateVGL (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int iat, const RefVector< ValueVector > &psi_v_list, const RefVector< GradVector > &dpsi_v_list, const RefVector< ValueVector > &d2psi_v_list) const
 evaluate the values, gradients and laplacians of this single-particle orbital sets of multiple walkers More...
 
virtual void mw_evaluateVGLandDetRatioGrads (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int iat, const std::vector< const ValueType * > &invRow_ptr_list, OffloadMWVGLArray &phi_vgl_v, std::vector< ValueType > &ratios, std::vector< GradType > &grads) const
 evaluate the values, gradients and laplacians of this single-particle orbital sets and determinant ratio and grads of multiple walkers. More...
 
virtual void acquireResource (ResourceCollection &collection, const RefVectorWithLeader< SPOSet > &spo_list) const
 acquire a shared resource from collection More...
 
virtual void createResource (ResourceCollection &collection) const
 initialize a shared resource and hand it to collection More...
 
virtual void releaseResource (ResourceCollection &collection, const RefVectorWithLeader< SPOSet > &spo_list) const
 return a shared resource to collection More...
 
- Public Member Functions inherited from SPOSet
 SPOSet (const std::string &my_name)
 constructor More...
 
virtual ~SPOSet ()=default
 destructor More...
 
int size () const
 return the size of the orbital set Ye: this needs to be replaced by getOrbitalSetSize(); More...
 
void basic_report (const std::string &pad="") const
 print basic SPOSet information More...
 
virtual void report (const std::string &pad="") const
 print SPOSet information More...
 
int getOrbitalSetSize () const
 return the size of the orbitals More...
 
virtual bool isOptimizable () const
 Query if this SPOSet is optimizable. More...
 
virtual void extractOptimizableObjectRefs (UniqueOptObjRefs &opt_obj_refs)
 extract underlying OptimizableObject references More...
 
virtual void checkOutVariables (const opt_variables_type &active)
 check out variational optimizable variables More...
 
virtual bool isOMPoffload () const
 Query if this SPOSet uses OpenMP offload. More...
 
virtual bool hasIonDerivs () const
 Query if this SPOSet has an explicit ion dependence. More...
 
virtual void checkObject () const
 check a few key parameters before putting the SPO into a determinant More...
 
virtual void evaluateDerivatives (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, Vector< ValueType > &dhpsioverpsi, const int &FirstIndex, const int &LastIndex)
 Parameter derivatives of the wavefunction and the Laplacian of the wavefunction. More...
 
virtual void evaluateDerivativesWF (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, int FirstIndex, int LastIndex)
 Parameter derivatives of the wavefunction. More...
 
virtual void evaluateDerivatives (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, Vector< ValueType > &dhpsioverpsi, const ValueType &psiCurrent, const std::vector< ValueType > &Coeff, const std::vector< size_t > &C2node_up, const std::vector< size_t > &C2node_dn, const ValueVector &detValues_up, const ValueVector &detValues_dn, const GradMatrix &grads_up, const GradMatrix &grads_dn, const ValueMatrix &lapls_up, const ValueMatrix &lapls_dn, const ValueMatrix &M_up, const ValueMatrix &M_dn, const ValueMatrix &Minv_up, const ValueMatrix &Minv_dn, const GradMatrix &B_grad, const ValueMatrix &B_lapl, const std::vector< int > &detData_up, const size_t N1, const size_t N2, const size_t NP1, const size_t NP2, const std::vector< std::vector< int >> &lookup_tbl)
 Evaluate the derivative of the optimized orbitals with respect to the parameters this is used only for MSD, to be refined for better serving both single and multi SD. More...
 
virtual void evaluateDerivativesWF (ParticleSet &P, const opt_variables_type &optvars, Vector< ValueType > &dlogpsi, const QTFull::ValueType &psiCurrent, const std::vector< ValueType > &Coeff, const std::vector< size_t > &C2node_up, const std::vector< size_t > &C2node_dn, const ValueVector &detValues_up, const ValueVector &detValues_dn, const ValueMatrix &M_up, const ValueMatrix &M_dn, const ValueMatrix &Minv_up, const ValueMatrix &Minv_dn, const std::vector< int > &detData_up, const std::vector< std::vector< int >> &lookup_tbl)
 Evaluate the derivative of the optimized orbitals with respect to the parameters this is used only for MSD, to be refined for better serving both single and multi SD. More...
 
virtual void evaluateDetRatios (const VirtualParticleSet &VP, ValueVector &psi, const ValueVector &psiinv, std::vector< ValueType > &ratios)
 evaluate determinant ratios for virtual moves, e.g., sphere move for nonlocalPP More...
 
virtual void evaluateDetSpinorRatios (const VirtualParticleSet &VP, ValueVector &psi, const std::pair< ValueVector, ValueVector > &spinor_multiplier, const ValueVector &invrow, std::vector< ValueType > &ratios)
 evaluate determinant ratios for virtual moves, specifically for Spinor SPOSets More...
 
virtual void evaluateDerivRatios (const VirtualParticleSet &VP, const opt_variables_type &optvars, ValueVector &psi, const ValueVector &psiinv, std::vector< ValueType > &ratios, Matrix< ValueType > &dratios, int FirstIndex, int LastIndex)
 Determinant ratios and parameter derivatives of the wavefunction for virtual moves. More...
 
virtual void mw_evaluateDetRatios (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< const VirtualParticleSet > &vp_list, const RefVector< ValueVector > &psi_list, const std::vector< const ValueType *> &invRow_ptr_list, std::vector< std::vector< ValueType >> &ratios_list) const
 evaluate determinant ratios for virtual moves, e.g., sphere move for nonlocalPP, of multiple walkers More...
 
virtual void evaluateVGL_spin (const ParticleSet &P, int iat, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi, ValueVector &dspin)
 evaluate the values, gradients and laplacians and spin gradient of this single-particle orbital set More...
 
virtual void mw_evaluateValue (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int iat, const RefVector< ValueVector > &psi_v_list) const
 evaluate the values this single-particle orbital sets of multiple walkers More...
 
virtual void mw_evaluateVGL (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int iat, const RefVector< ValueVector > &psi_v_list, const RefVector< GradVector > &dpsi_v_list, const RefVector< ValueVector > &d2psi_v_list) const
 evaluate the values, gradients and laplacians of this single-particle orbital sets of multiple walkers More...
 
virtual void mw_evaluateVGLWithSpin (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int iat, const RefVector< ValueVector > &psi_v_list, const RefVector< GradVector > &dpsi_v_list, const RefVector< ValueVector > &d2psi_v_list, OffloadMatrix< ComplexType > &mw_dspin) const
 evaluate the values, gradients and laplacians and spin gradient of this single-particle orbital sets of multiple walkers More...
 
virtual void mw_evaluateVGLandDetRatioGrads (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int iat, const std::vector< const ValueType *> &invRow_ptr_list, OffloadMWVGLArray &phi_vgl_v, std::vector< ValueType > &ratios, std::vector< GradType > &grads) const
 evaluate the values, gradients and laplacians of this single-particle orbital sets and determinant ratio and grads of multiple walkers. More...
 
virtual void mw_evaluateVGLandDetRatioGradsWithSpin (const RefVectorWithLeader< SPOSet > &spo_list, const RefVectorWithLeader< ParticleSet > &P_list, int iat, const std::vector< const ValueType *> &invRow_ptr_list, OffloadMWVGLArray &phi_vgl_v, std::vector< ValueType > &ratios, std::vector< GradType > &grads, std::vector< ValueType > &spingrads) const
 evaluate the values, gradients and laplacians of this single-particle orbital sets and determinant ratio and grads of multiple walkers. More...
 
virtual void evaluate_spin (const ParticleSet &P, int iat, ValueVector &psi, ValueVector &dpsi)
 evaluate the values of this single-particle orbital set More...
 
virtual void evaluateThirdDeriv (const ParticleSet &P, int first, int last, GGGMatrix &grad_grad_grad_logdet)
 evaluate the third derivatives of this single-particle orbital set More...
 
virtual void evaluate_notranspose_spin (const ParticleSet &P, int first, int last, ValueMatrix &logdet, GradMatrix &dlogdet, ValueMatrix &d2logdet, ValueMatrix &dspinlogdet)
 evaluate the values, gradients and laplacians of this single-particle orbital for [first,last) particles, including the spin gradient More...
 
virtual void evaluateGradSourceRow (const ParticleSet &P, int iel, const ParticleSet &source, int iat_src, GradVector &gradphi)
 Returns a row of d/dR_iat phi_j(r) evaluated at position r. More...
 
virtual PosType get_k (int orb)
 access the k point related to the given orbital More...
 
virtual void createResource (ResourceCollection &collection) const
 initialize a shared resource and hand it to collection More...
 
virtual void acquireResource (ResourceCollection &collection, const RefVectorWithLeader< SPOSet > &spo_list) const
 acquire a shared resource from collection More...
 
virtual void releaseResource (ResourceCollection &collection, const RefVectorWithLeader< SPOSet > &spo_list) const
 return a shared resource to collection More...
 
virtual bool transformSPOSet ()
 Used only by cusp correction in AOS LCAO. More...
 
virtual void finalizeConstruction ()
 finalize the construction of SPOSet More...
 
const std::string & getName () const
 return object name More...
 

Protected Attributes

CrystalLattice< ST, 3 > PrimLattice
 primitive cell More...
 
vContainer_type myV
 intermediate result vectors More...
 
vContainer_type myL
 
gContainer_type myG
 
hContainer_type myH
 
ghContainer_type mygH
 
- Protected Attributes inherited from BsplineSet
size_t MyIndex
 Index of this adoptor, when multiple adoptors are used for NUMA or distributed cases. More...
 
size_t first_spo
 first index of the SPOs this Spline handles More...
 
size_t last_spo
 last index of the SPOs this Spline handles More...
 
TinyVector< int, DHalfG
 sign bits at the G/2 boundaries More...
 
std::vector< bool > MakeTwoCopies
 flags to unpack sin/cos More...
 
std::vector< SPOSet::PosTypekPoints
 kpoints for each unique orbitals. More...
 
aligned_vector< int > BandIndexMap
 remap splines to orbitals More...
 
std::vector< int > offset
 band offsets used for communication More...
 
- Protected Attributes inherited from SPOSet
const std::string my_name_
 name of the object, unique identifier More...
 
IndexType OrbitalSetSize
 number of Single-particle orbitals More...
 
opt_variables_type myVars
 Optimizable variables. More...
 

Private Attributes

bool IsGamma
 
Tensor< ST, 3 > GGt
 $GGt=G^t G $, transformation for tensor in LatticeUnit to CartesianUnit, e.g. Hessian More...
 
std::shared_ptr< MultiBspline< ST > > SplineInst
 multi bspline set More...
 
std::shared_ptr< std::vector< ST > > coef_copy_
 Copy of original splines for orbital rotation. More...
 
Matrix< TTratios_private
 thread private ratios for reduction when using nested threading, numVP x numThread More...
 

Friends

template<class BSPLINESPO >
class SplineSetReader
 
struct BsplineReader
 

Additional Inherited Members

- Static Protected Attributes inherited from BsplineSet
static const int D = DIM
 

Detailed Description

template<typename ST>
class qmcplusplus::SplineR2R< ST >

class to match ST real spline with BsplineSet::ValueType (real) SPOs

Template Parameters
STprecision of spline

Requires temporage storage and multiplication of the sign of the real part of the phase Internal storage ST type arrays are aligned and padded.

Definition at line 34 of file SplineR2R.h.

Member Typedef Documentation

◆ BCType

using BCType = typename bspline_traits<ST, 3>::BCType

Definition at line 38 of file SplineR2R.h.

◆ DataType

using DataType = ST

Definition at line 39 of file SplineR2R.h.

◆ gContainer_type

Definition at line 50 of file SplineR2R.h.

◆ ghContainer_type

Definition at line 52 of file SplineR2R.h.

◆ hContainer_type

Definition at line 51 of file SplineR2R.h.

◆ PointType

using PointType = TinyVector<ST, 3>

Definition at line 40 of file SplineR2R.h.

◆ SingleSplineType

using SingleSplineType = UBspline_3d_d

Definition at line 41 of file SplineR2R.h.

◆ SplineType

using SplineType = typename bspline_traits<ST, 3>::SplineType

Definition at line 37 of file SplineR2R.h.

◆ TT

using TT = typename BsplineSet::ValueType

Definition at line 43 of file SplineR2R.h.

◆ vContainer_type

Definition at line 49 of file SplineR2R.h.

Constructor & Destructor Documentation

◆ SplineR2R() [1/2]

SplineR2R ( const std::string &  my_name)
inline

Definition at line 79 of file SplineR2R.h.

79 : BsplineSet(my_name) {}
BsplineSet(const std::string &my_name)
Definition: BsplineSet.h:59

◆ SplineR2R() [2/2]

SplineR2R ( const SplineR2R< ST > &  in)
default

Member Function Documentation

◆ applyRotation()

void applyRotation ( const ValueMatrix rot_mat,
bool  use_stored_copy 
)
overridevirtual

apply rotation to all the orbitals

Reimplemented from SPOSet.

Definition at line 106 of file SplineR2R.cpp.

References qmcplusplus::syclBLAS::copy_n(), and BLAS::gemm().

107 {
108  // SplineInst is a MultiBspline. See src/spline2/MultiBspline.hpp
109  const auto spline_ptr = SplineInst->getSplinePtr();
110  assert(spline_ptr != nullptr);
111  const auto spl_coefs = spline_ptr->coefs;
112  const auto Nsplines = spline_ptr->num_splines; // May include padding
113  const auto coefs_tot_size = spline_ptr->coefs_size;
114  const auto BasisSetSize = coefs_tot_size / Nsplines;
115  const auto TrueNOrbs = rot_mat.size1(); // == Nsplines - padding
116  assert(OrbitalSetSize == rot_mat.rows());
117  assert(OrbitalSetSize == rot_mat.cols());
118 
119  if (!use_stored_copy)
120  {
121  assert(coef_copy_ != nullptr);
122  std::copy_n(spl_coefs, coefs_tot_size, coef_copy_->begin());
123  }
124 
125 
126  if constexpr (std::is_same_v<ST, RealType>)
127  {
128  //Here, ST should be equal to ValueType, which will be double for R2R. Using BLAS to make things faster
129  BLAS::gemm('N', 'N', OrbitalSetSize, BasisSetSize, OrbitalSetSize, ST(1.0), rot_mat.data(), OrbitalSetSize,
130  coef_copy_->data(), Nsplines, ST(0.0), spl_coefs, Nsplines);
131  }
132  else
133  {
134  //Here, ST is float but ValueType is double for R2R. Due to issues with type conversions, just doing naive matrix multiplication in this case to not lose precision on rot_mat
135  for (IndexType i = 0; i < BasisSetSize; i++)
136  for (IndexType j = 0; j < OrbitalSetSize; j++)
137  {
138  const auto cur_elem = Nsplines * i + j;
139  FullPrecValueType newval{0.};
140  for (IndexType k = 0; k < OrbitalSetSize; k++)
141  {
142  const auto index = i * Nsplines + k;
143  newval += (*coef_copy_)[index] * rot_mat[k][j];
144  }
145  spl_coefs[cur_elem] = newval;
146  }
147  }
148 }
std::shared_ptr< std::vector< ST > > coef_copy_
Copy of original splines for orbital rotation.
Definition: SplineR2R.h:62
QTFull::ValueType FullPrecValueType
Definition: Configuration.h:67
IndexType OrbitalSetSize
number of Single-particle orbitals
Definition: SPOSet.h:566
OHMMS_INDEXTYPE IndexType
define other types
Definition: Configuration.h:65
sycl::event copy_n(sycl::queue &aq, const T1 *restrict VA, size_t array_size, T2 *restrict VC, const std::vector< sycl::event > &events)
Definition: syclBLAS.cpp:548
static void gemm(char Atrans, char Btrans, int M, int N, int K, double alpha, const double *A, int lda, const double *restrict B, int ldb, double beta, double *restrict C, int ldc)
Definition: BLAS.hpp:235
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ assign_v()

void assign_v ( int  bc_sign,
const vContainer_type myV,
ValueVector psi,
int  first,
int  last 
) const
inline

Definition at line 152 of file SplineR2R.cpp.

References omptarget::min().

154 {
155  // protect last against kPoints.size() and psi.size()
156  size_t last_real = std::min(kPoints.size(), psi.size());
157  last = last > last_real ? last_real : last;
158 
159  const ST signed_one = (bc_sign & 1) ? -1 : 1;
160 #pragma omp simd
161  for (size_t j = first; j < last; ++j)
162  psi[first_spo + j] = signed_one * myV[j];
163 }
T min(T a, T b)
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
size_t first_spo
first index of the SPOs this Spline handles
Definition: BsplineSet.h:41
std::vector< SPOSet::PosType > kPoints
kpoints for each unique orbitals.
Definition: BsplineSet.h:52

◆ assign_vgh()

void assign_vgh ( int  bc_sign,
ValueVector psi,
GradVector dpsi,
HessVector grad_grad_psi,
int  first,
int  last 
) const

Definition at line 310 of file SplineR2R.cpp.

References omptarget::min(), and qmcplusplus::v_m_v().

316 {
317  // protect last against kPoints.size() and psi.size()
318  const size_t last_real = std::min(kPoints.size(), psi.size());
319  last = last > last_real ? last_real : last;
320 
321  const ST signed_one = (bc_sign & 1) ? -1 : 1;
322  const ST g00 = PrimLattice.G(0), g01 = PrimLattice.G(1), g02 = PrimLattice.G(2), g10 = PrimLattice.G(3),
323  g11 = PrimLattice.G(4), g12 = PrimLattice.G(5), g20 = PrimLattice.G(6), g21 = PrimLattice.G(7),
324  g22 = PrimLattice.G(8);
325 
326  const ST* restrict g0 = myG.data(0);
327  const ST* restrict g1 = myG.data(1);
328  const ST* restrict g2 = myG.data(2);
329  const ST* restrict h00 = myH.data(0);
330  const ST* restrict h01 = myH.data(1);
331  const ST* restrict h02 = myH.data(2);
332  const ST* restrict h11 = myH.data(3);
333  const ST* restrict h12 = myH.data(4);
334  const ST* restrict h22 = myH.data(5);
335 
336 #pragma omp simd
337  for (size_t j = first; j < last; ++j)
338  {
339  //dot(PrimLattice.G,myG[j])
340  const ST dX_r = g00 * g0[j] + g01 * g1[j] + g02 * g2[j];
341  const ST dY_r = g10 * g0[j] + g11 * g1[j] + g12 * g2[j];
342  const ST dZ_r = g20 * g0[j] + g21 * g1[j] + g22 * g2[j];
343 
344  const size_t psiIndex = j + first_spo;
345  psi[psiIndex] = signed_one * myV[j];
346  dpsi[psiIndex][0] = signed_one * dX_r;
347  dpsi[psiIndex][1] = signed_one * dY_r;
348  dpsi[psiIndex][2] = signed_one * dZ_r;
349 
350  const ST h_xx_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g00, g01, g02, g00, g01, g02);
351  const ST h_xy_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g00, g01, g02, g10, g11, g12);
352  const ST h_xz_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g00, g01, g02, g20, g21, g22);
353  const ST h_yx_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g10, g11, g12, g00, g01, g02);
354  const ST h_yy_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g10, g11, g12, g10, g11, g12);
355  const ST h_yz_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g10, g11, g12, g20, g21, g22);
356  const ST h_zx_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g20, g21, g22, g00, g01, g02);
357  const ST h_zy_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g20, g21, g22, g10, g11, g12);
358  const ST h_zz_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g20, g21, g22, g20, g21, g22);
359 
360  grad_grad_psi[psiIndex][0] = signed_one * h_xx_r;
361  grad_grad_psi[psiIndex][1] = signed_one * h_xy_r;
362  grad_grad_psi[psiIndex][2] = signed_one * h_xz_r;
363  grad_grad_psi[psiIndex][3] = signed_one * h_yx_r;
364  grad_grad_psi[psiIndex][4] = signed_one * h_yy_r;
365  grad_grad_psi[psiIndex][5] = signed_one * h_yz_r;
366  grad_grad_psi[psiIndex][6] = signed_one * h_zx_r;
367  grad_grad_psi[psiIndex][7] = signed_one * h_zy_r;
368  grad_grad_psi[psiIndex][8] = signed_one * h_zz_r;
369  }
370 }
CrystalLattice< ST, 3 > PrimLattice
primitive cell
Definition: SplineR2R.h:70
hContainer_type myH
Definition: SplineR2R.h:75
T v_m_v(T h00, T h01, T h02, T h11, T h12, T h22, T g1x, T g1y, T g1z, T g2x, T g2y, T g2z)
compute vector[3]^T x matrix[3][3] x vector[3]
Tensor_t G
Reciprocal unit vectors. G(j,i) i=vector and j=x,y,z.
T min(T a, T b)
gContainer_type myG
Definition: SplineR2R.h:74
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
size_t first_spo
first index of the SPOs this Spline handles
Definition: BsplineSet.h:41
std::vector< SPOSet::PosType > kPoints
kpoints for each unique orbitals.
Definition: BsplineSet.h:52

◆ assign_vghgh()

void assign_vghgh ( int  bc_sign,
ValueVector psi,
GradVector dpsi,
HessVector grad_grad_psi,
GGGVector grad_grad_grad_psi,
int  first = 0,
int  last = -1 
) const

Definition at line 394 of file SplineR2R.cpp.

References omptarget::min(), qmcplusplus::t3_contract(), and qmcplusplus::v_m_v().

401 {
402  // protect last against kPoints.size() and psi.size()
403  const size_t last_real = std::min(kPoints.size(), psi.size());
404  last = last < 0 ? last_real : (last > last_real ? last_real : last);
405 
406  const ST signed_one = (bc_sign & 1) ? -1 : 1;
407  const ST g00 = PrimLattice.G(0), g01 = PrimLattice.G(1), g02 = PrimLattice.G(2), g10 = PrimLattice.G(3),
408  g11 = PrimLattice.G(4), g12 = PrimLattice.G(5), g20 = PrimLattice.G(6), g21 = PrimLattice.G(7),
409  g22 = PrimLattice.G(8);
410 
411  const ST* restrict g0 = myG.data(0);
412  const ST* restrict g1 = myG.data(1);
413  const ST* restrict g2 = myG.data(2);
414  const ST* restrict h00 = myH.data(0);
415  const ST* restrict h01 = myH.data(1);
416  const ST* restrict h02 = myH.data(2);
417  const ST* restrict h11 = myH.data(3);
418  const ST* restrict h12 = myH.data(4);
419  const ST* restrict h22 = myH.data(5);
420 
421  const ST* restrict gh000 = mygH.data(0);
422  const ST* restrict gh001 = mygH.data(1);
423  const ST* restrict gh002 = mygH.data(2);
424  const ST* restrict gh011 = mygH.data(3);
425  const ST* restrict gh012 = mygH.data(4);
426  const ST* restrict gh022 = mygH.data(5);
427  const ST* restrict gh111 = mygH.data(6);
428  const ST* restrict gh112 = mygH.data(7);
429  const ST* restrict gh122 = mygH.data(8);
430  const ST* restrict gh222 = mygH.data(9);
431 
432  //SIMD doesn't work quite right yet. Comment out until further debugging.
433  //#pragma omp simd
434  for (size_t j = first; j < last; ++j)
435  {
436  const ST val_r = myV[j];
437 
438 
439  //dot(PrimLattice.G,myG[j])
440  const ST dX_r = g00 * g0[j] + g01 * g1[j] + g02 * g2[j];
441  const ST dY_r = g10 * g0[j] + g11 * g1[j] + g12 * g2[j];
442  const ST dZ_r = g20 * g0[j] + g21 * g1[j] + g22 * g2[j];
443 
444  const size_t psiIndex = j + first_spo;
445  psi[psiIndex] = signed_one * val_r;
446  dpsi[psiIndex][0] = signed_one * dX_r;
447  dpsi[psiIndex][1] = signed_one * dY_r;
448  dpsi[psiIndex][2] = signed_one * dZ_r;
449 
450  //intermediates for computation of hessian. \partial_i \partial_j phi in cartesian coordinates.
451  const ST f_xx_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g00, g01, g02, g00, g01, g02);
452  const ST f_xy_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g00, g01, g02, g10, g11, g12);
453  const ST f_xz_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g00, g01, g02, g20, g21, g22);
454  const ST f_yy_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g10, g11, g12, g10, g11, g12);
455  const ST f_yz_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g10, g11, g12, g20, g21, g22);
456  const ST f_zz_r = v_m_v(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], g20, g21, g22, g20, g21, g22);
457 
458  /* const ST h_xx_r=f_xx_r;
459  const ST h_xy_r=f_xy_r+(kX*dY_i+kY*dX_i)-kX*kY*val_r;
460  const ST h_xz_r=f_xz_r+(kX*dZ_i+kZ*dX_i)-kX*kZ*val_r;
461  const ST h_yy_r=f_yy_r+2*kY*dY_i-kY*kY*val_r;
462  const ST h_yz_r=f_yz_r+(kY*dZ_i+kZ*dY_i)-kY*kZ*val_r;
463  const ST h_zz_r=f_zz_r+2*kZ*dZ_i-kZ*kZ*val_r; */
464 
465  grad_grad_psi[psiIndex][0] = f_xx_r * signed_one;
466  grad_grad_psi[psiIndex][1] = f_xy_r * signed_one;
467  grad_grad_psi[psiIndex][2] = f_xz_r * signed_one;
468  grad_grad_psi[psiIndex][4] = f_yy_r * signed_one;
469  grad_grad_psi[psiIndex][5] = f_yz_r * signed_one;
470  grad_grad_psi[psiIndex][8] = f_zz_r * signed_one;
471 
472  //symmetry:
473  grad_grad_psi[psiIndex][3] = grad_grad_psi[psiIndex][1];
474  grad_grad_psi[psiIndex][6] = grad_grad_psi[psiIndex][2];
475  grad_grad_psi[psiIndex][7] = grad_grad_psi[psiIndex][5];
476  //These are the real and imaginary components of the third SPO derivative. _xxx denotes
477  // third derivative w.r.t. x, _xyz, a derivative with resepect to x,y, and z, and so on.
478 
479  const ST f3_xxx_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
480  gh122[j], gh222[j], g00, g01, g02, g00, g01, g02, g00, g01, g02);
481  const ST f3_xxy_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
482  gh122[j], gh222[j], g00, g01, g02, g00, g01, g02, g10, g11, g12);
483  const ST f3_xxz_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
484  gh122[j], gh222[j], g00, g01, g02, g00, g01, g02, g20, g21, g22);
485  const ST f3_xyy_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
486  gh122[j], gh222[j], g00, g01, g02, g10, g11, g12, g10, g11, g12);
487  const ST f3_xyz_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
488  gh122[j], gh222[j], g00, g01, g02, g10, g11, g12, g20, g21, g22);
489  const ST f3_xzz_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
490  gh122[j], gh222[j], g00, g01, g02, g20, g21, g22, g20, g21, g22);
491  const ST f3_yyy_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
492  gh122[j], gh222[j], g10, g11, g12, g10, g11, g12, g10, g11, g12);
493  const ST f3_yyz_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
494  gh122[j], gh222[j], g10, g11, g12, g10, g11, g12, g20, g21, g22);
495  const ST f3_yzz_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
496  gh122[j], gh222[j], g10, g11, g12, g20, g21, g22, g20, g21, g22);
497  const ST f3_zzz_r = t3_contract(gh000[j], gh001[j], gh002[j], gh011[j], gh012[j], gh022[j], gh111[j], gh112[j],
498  gh122[j], gh222[j], g20, g21, g22, g20, g21, g22, g20, g21, g22);
499 
500  //Here is where we build up the components of the physical hessian gradient, namely, d^3/dx^3(e^{-ik*r}\phi(r)
501  /* const ST gh_xxx_r= f3_xxx_r + 3*kX*f_xx_i - 3*kX*kX*dX_r - kX*kX*kX*val_i;
502  const ST gh_xxy_r= f3_xxy_r +(kY*f_xx_i+2*kX*f_xy_i) - (kX*kX*dY_r+2*kX*kY*dX_r)-kX*kX*kY*val_i;
503  const ST gh_xxz_r= f3_xxz_r +(kZ*f_xx_i+2*kX*f_xz_i) - (kX*kX*dZ_r+2*kX*kZ*dX_r)-kX*kX*kZ*val_i;
504  const ST gh_xyy_r= f3_xyy_r +(2*kY*f_xy_i+kX*f_yy_i) - (2*kX*kY*dY_r+kY*kY*dX_r)-kX*kY*kY*val_i;
505  const ST gh_xyz_r= f3_xyz_r +(kX*f_yz_i+kY*f_xz_i+kZ*f_xy_i)-(kX*kY*dZ_r+kY*kZ*dX_r+kZ*kX*dY_r) - kX*kY*kZ*val_i;
506  const ST gh_xzz_r= f3_xzz_r +(2*kZ*f_xz_i+kX*f_zz_i) - (2*kX*kZ*dZ_r+kZ*kZ*dX_r)-kX*kZ*kZ*val_i;
507  const ST gh_yyy_r= f3_yyy_r + 3*kY*f_yy_i - 3*kY*kY*dY_r - kY*kY*kY*val_i;
508  const ST gh_yyz_r= f3_yyz_r +(kZ*f_yy_i+2*kY*f_yz_i) - (kY*kY*dZ_r+2*kY*kZ*dY_r)-kY*kY*kZ*val_i;
509  const ST gh_yzz_r= f3_yzz_r +(2*kZ*f_yz_i+kY*f_zz_i) - (2*kY*kZ*dZ_r+kZ*kZ*dY_r)-kY*kZ*kZ*val_i;
510  const ST gh_zzz_r= f3_zzz_r + 3*kZ*f_zz_i - 3*kZ*kZ*dZ_r - kZ*kZ*kZ*val_i;*/
511  //[x][xx] //These are the unique entries
512  grad_grad_grad_psi[psiIndex][0][0] = signed_one * f3_xxx_r;
513  grad_grad_grad_psi[psiIndex][0][1] = signed_one * f3_xxy_r;
514  grad_grad_grad_psi[psiIndex][0][2] = signed_one * f3_xxz_r;
515  grad_grad_grad_psi[psiIndex][0][4] = signed_one * f3_xyy_r;
516  grad_grad_grad_psi[psiIndex][0][5] = signed_one * f3_xyz_r;
517  grad_grad_grad_psi[psiIndex][0][8] = signed_one * f3_xzz_r;
518 
519  //filling in the symmetric terms. Filling out the xij terms
520  grad_grad_grad_psi[psiIndex][0][3] = grad_grad_grad_psi[psiIndex][0][1];
521  grad_grad_grad_psi[psiIndex][0][6] = grad_grad_grad_psi[psiIndex][0][2];
522  grad_grad_grad_psi[psiIndex][0][7] = grad_grad_grad_psi[psiIndex][0][5];
523 
524  //Now for everything that's a permutation of the above:
525  grad_grad_grad_psi[psiIndex][1][0] = grad_grad_grad_psi[psiIndex][0][1];
526  grad_grad_grad_psi[psiIndex][1][1] = grad_grad_grad_psi[psiIndex][0][4];
527  grad_grad_grad_psi[psiIndex][1][2] = grad_grad_grad_psi[psiIndex][0][5];
528  grad_grad_grad_psi[psiIndex][1][3] = grad_grad_grad_psi[psiIndex][0][4];
529  grad_grad_grad_psi[psiIndex][1][6] = grad_grad_grad_psi[psiIndex][0][5];
530 
531  grad_grad_grad_psi[psiIndex][2][0] = grad_grad_grad_psi[psiIndex][0][2];
532  grad_grad_grad_psi[psiIndex][2][1] = grad_grad_grad_psi[psiIndex][0][5];
533  grad_grad_grad_psi[psiIndex][2][2] = grad_grad_grad_psi[psiIndex][0][8];
534  grad_grad_grad_psi[psiIndex][2][3] = grad_grad_grad_psi[psiIndex][0][5];
535  grad_grad_grad_psi[psiIndex][2][6] = grad_grad_grad_psi[psiIndex][0][8];
536 
537  grad_grad_grad_psi[psiIndex][1][4] = signed_one * f3_yyy_r;
538  grad_grad_grad_psi[psiIndex][1][5] = signed_one * f3_yyz_r;
539  grad_grad_grad_psi[psiIndex][1][8] = signed_one * f3_yzz_r;
540 
541  grad_grad_grad_psi[psiIndex][1][7] = grad_grad_grad_psi[psiIndex][1][5];
542  grad_grad_grad_psi[psiIndex][2][4] = grad_grad_grad_psi[psiIndex][1][5];
543  grad_grad_grad_psi[psiIndex][2][5] = grad_grad_grad_psi[psiIndex][1][8];
544  grad_grad_grad_psi[psiIndex][2][7] = grad_grad_grad_psi[psiIndex][1][8];
545 
546  grad_grad_grad_psi[psiIndex][2][8] = signed_one * f3_zzz_r;
547  }
548 }
CrystalLattice< ST, 3 > PrimLattice
primitive cell
Definition: SplineR2R.h:70
hContainer_type myH
Definition: SplineR2R.h:75
T v_m_v(T h00, T h01, T h02, T h11, T h12, T h22, T g1x, T g1y, T g1z, T g2x, T g2y, T g2z)
compute vector[3]^T x matrix[3][3] x vector[3]
T t3_contract(T h000, T h001, T h002, T h011, T h012, T h022, T h111, T h112, T h122, T h222, T g1x, T g1y, T g1z, T g2x, T g2y, T g2z, T g3x, T g3y, T g3z)
Coordinate transform for a 3rd rank symmetric tensor representing coordinate derivatives (hence t3_co...
Tensor_t G
Reciprocal unit vectors. G(j,i) i=vector and j=x,y,z.
T min(T a, T b)
ghContainer_type mygH
Definition: SplineR2R.h:76
gContainer_type myG
Definition: SplineR2R.h:74
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
size_t first_spo
first index of the SPOs this Spline handles
Definition: BsplineSet.h:41
std::vector< SPOSet::PosType > kPoints
kpoints for each unique orbitals.
Definition: BsplineSet.h:52

◆ assign_vgl()

void assign_vgl ( int  bc_sign,
ValueVector psi,
GradVector dpsi,
ValueVector d2psi,
int  first,
int  last 
) const
inline

Definition at line 226 of file SplineR2R.cpp.

References omptarget::min(), and qmcplusplus::SymTrace().

232 {
233  // protect last against kPoints.size() and psi.size()
234  size_t last_real = std::min(kPoints.size(), psi.size());
235  last = last > last_real ? last_real : last;
236 
237  const ST signed_one = (bc_sign & 1) ? -1 : 1;
238  const ST g00 = PrimLattice.G(0), g01 = PrimLattice.G(1), g02 = PrimLattice.G(2), g10 = PrimLattice.G(3),
239  g11 = PrimLattice.G(4), g12 = PrimLattice.G(5), g20 = PrimLattice.G(6), g21 = PrimLattice.G(7),
240  g22 = PrimLattice.G(8);
241  const ST symGG[6] = {GGt[0], GGt[1] + GGt[3], GGt[2] + GGt[6], GGt[4], GGt[5] + GGt[7], GGt[8]};
242 
243  const ST* restrict g0 = myG.data(0);
244  const ST* restrict g1 = myG.data(1);
245  const ST* restrict g2 = myG.data(2);
246  const ST* restrict h00 = myH.data(0);
247  const ST* restrict h01 = myH.data(1);
248  const ST* restrict h02 = myH.data(2);
249  const ST* restrict h11 = myH.data(3);
250  const ST* restrict h12 = myH.data(4);
251  const ST* restrict h22 = myH.data(5);
252 
253 #pragma omp simd
254  for (size_t j = first; j < last; ++j)
255  {
256  const size_t psiIndex = first_spo + j;
257  psi[psiIndex] = signed_one * myV[j];
258  dpsi[psiIndex][0] = signed_one * (g00 * g0[j] + g01 * g1[j] + g02 * g2[j]);
259  dpsi[psiIndex][1] = signed_one * (g10 * g0[j] + g11 * g1[j] + g12 * g2[j]);
260  dpsi[psiIndex][2] = signed_one * (g20 * g0[j] + g21 * g1[j] + g22 * g2[j]);
261  d2psi[psiIndex] = signed_one * SymTrace(h00[j], h01[j], h02[j], h11[j], h12[j], h22[j], symGG);
262  }
263 }
CrystalLattice< ST, 3 > PrimLattice
primitive cell
Definition: SplineR2R.h:70
T SymTrace(T h00, T h01, T h02, T h11, T h12, T h22, const T gg[6])
compute Trace(H*G)
Tensor< ST, 3 > GGt
, transformation for tensor in LatticeUnit to CartesianUnit, e.g. Hessian
Definition: SplineR2R.h:57
hContainer_type myH
Definition: SplineR2R.h:75
Tensor_t G
Reciprocal unit vectors. G(j,i) i=vector and j=x,y,z.
T min(T a, T b)
gContainer_type myG
Definition: SplineR2R.h:74
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
size_t first_spo
first index of the SPOs this Spline handles
Definition: BsplineSet.h:41
std::vector< SPOSet::PosType > kPoints
kpoints for each unique orbitals.
Definition: BsplineSet.h:52

◆ assign_vgl_from_l()

void assign_vgl_from_l ( int  bc_sign,
ValueVector psi,
GradVector dpsi,
ValueVector d2psi 
)
inline

assign_vgl_from_l can be used when myL is precomputed and myV,myG,myL in cartesian

Definition at line 268 of file SplineR2R.cpp.

269 {
270  const ST signed_one = (bc_sign & 1) ? -1 : 1;
271  const ST* restrict g0 = myG.data(0);
272  const ST* restrict g1 = myG.data(1);
273  const ST* restrict g2 = myG.data(2);
274 
275  const size_t last_real = last_spo > psi.size() ? psi.size() : last_spo;
276 #pragma omp simd
277  for (int psiIndex = first_spo; psiIndex < last_real; ++psiIndex)
278  {
279  const size_t j = psiIndex - first_spo;
280  psi[psiIndex] = signed_one * myV[j];
281  dpsi[psiIndex][0] = signed_one * g0[j];
282  dpsi[psiIndex][1] = signed_one * g1[j];
283  dpsi[psiIndex][2] = signed_one * g2[j];
284  d2psi[psiIndex] = signed_one * myL[j];
285  }
286 }
vContainer_type myL
Definition: SplineR2R.h:73
gContainer_type myG
Definition: SplineR2R.h:74
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
size_t first_spo
first index of the SPOs this Spline handles
Definition: BsplineSet.h:41
size_t last_spo
last index of the SPOs this Spline handles
Definition: BsplineSet.h:43

◆ bcast_tables()

void bcast_tables ( Communicate comm)
inline

Definition at line 120 of file SplineR2R.h.

References qmcplusplus::comm, and SplineR2R< ST >::SplineInst.

120 { chunked_bcast(comm, SplineInst->getSplinePtr()); }
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ convertPos()

int convertPos ( const PointType r,
PointType ru 
)
inline

convert position in PrimLattice unit and return sign

Definition at line 153 of file SplineR2R.h.

References BsplineSet::D, qmcplusplus::floor(), BsplineSet::HalfG, SplineR2R< ST >::PrimLattice, and CrystalLattice< T, D >::toUnit().

154  {
155  ru = PrimLattice.toUnit(r);
156  int bc_sign = 0;
157  for (int i = 0; i < D; i++)
158  if (-std::numeric_limits<ST>::epsilon() < ru[i] && ru[i] < 0)
159  ru[i] = ST(0.0);
160  else
161  {
162  ST img = std::floor(ru[i]);
163  ru[i] -= img;
164  bc_sign += HalfG[i] * (int)img;
165  }
166  return bc_sign;
167  }
CrystalLattice< ST, 3 > PrimLattice
primitive cell
Definition: SplineR2R.h:70
static const int D
Definition: BsplineSet.h:37
SingleParticlePos toUnit(const TinyVector< T1, D > &r) const
Convert a cartesian vector to a unit vector.
TinyVector< int, D > HalfG
sign bits at the G/2 boundaries
Definition: BsplineSet.h:45
MakeReturn< UnaryNode< FnFloor, typename CreateLeaf< Vector< T1, C1 > >::Leaf_t > >::Expression_t floor(const Vector< T1, C1 > &l)

◆ create_spline()

void create_spline ( GT &  xyz_g,
BCT &  xyz_bc 
)
inline

Definition at line 134 of file SplineR2R.h.

References qmcplusplus::app_log(), qmcplusplus::dot(), CrystalLattice< T, D >::G, SplineR2R< ST >::GGt, SplineR2R< ST >::myV, SplineR2R< ST >::PrimLattice, Vector< T, Alloc >::size(), SplineR2R< ST >::SplineInst, and qmcplusplus::transpose().

135  {
137  SplineInst = std::make_shared<MultiBspline<ST>>();
138  SplineInst->create(xyz_g, xyz_bc, myV.size());
139 
140  app_log() << "MEMORY " << SplineInst->sizeInByte() / (1 << 20) << " MB allocated "
141  << "for the coefficients in 3D spline orbital representation" << std::endl;
142  }
CrystalLattice< ST, 3 > PrimLattice
primitive cell
Definition: SplineR2R.h:70
Tensor< ST, 3 > GGt
, transformation for tensor in LatticeUnit to CartesianUnit, e.g. Hessian
Definition: SplineR2R.h:57
std::ostream & app_log()
Definition: OutputManager.h:65
Tensor_t G
Reciprocal unit vectors. G(j,i) i=vector and j=x,y,z.
size_type size() const
return the current size
Definition: OhmmsVector.h:162
AntiSymTensor< T, D > transpose(const AntiSymTensor< T, D > &rhs)
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
Tensor< typename BinaryReturn< T1, T2, OpMultiply >::Type_t, D > dot(const AntiSymTensor< T1, D > &lhs, const AntiSymTensor< T2, D > &rhs)
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ evaluateDetRatios()

void evaluateDetRatios ( const VirtualParticleSet VP,
ValueVector psi,
const ValueVector psiinv,
std::vector< TT > &  ratios 
)
override

Definition at line 183 of file SplineR2R.cpp.

References ParticleSet::activeR(), qmcplusplus::C2C::assign_v(), qmcplusplus::simd::dot(), FairDivideAligned(), VirtualParticleSet::getTotalNum(), omp_get_num_threads(), and omp_get_thread_num().

187 {
188  const bool need_resize = ratios_private.rows() < VP.getTotalNum();
189 
190 #pragma omp parallel
191  {
192  int tid = omp_get_thread_num();
193  // initialize thread private ratios
194  if (need_resize)
195  {
196  if (tid == 0) // just like #pragma omp master, but one fewer call to the runtime
197  ratios_private.resize(VP.getTotalNum(), omp_get_num_threads());
198 #pragma omp barrier
199  }
200  int first, last;
201  FairDivideAligned(psi.size(), getAlignment<ST>(), omp_get_num_threads(), tid, first, last);
202  const int last_real = kPoints.size() < last ? kPoints.size() : last;
203 
204  for (int iat = 0; iat < VP.getTotalNum(); ++iat)
205  {
206  const PointType& r = VP.activeR(iat);
207  PointType ru;
208  int bc_sign = convertPos(r, ru);
209 
210  spline2::evaluate3d(SplineInst->getSplinePtr(), ru, myV, first, last);
211  assign_v(bc_sign, myV, psi, first, last_real);
212  ratios_private[iat][tid] = simd::dot(psi.data() + first, psiinv.data() + first, last_real - first);
213  }
214  }
215 
216  // do the reduction manually
217  for (int iat = 0; iat < VP.getTotalNum(); ++iat)
218  {
219  ratios[iat] = TT(0);
220  for (int tid = 0; tid < ratios_private.cols(); tid++)
221  ratios[iat] += ratios_private[iat][tid];
222  }
223 }
int convertPos(const PointType &r, PointType &ru)
convert position in PrimLattice unit and return sign
Definition: SplineR2R.h:153
Matrix< TT > ratios_private
thread private ratios for reduction when using nested threading, numVP x numThread ...
Definition: SplineR2R.h:65
T dot(const T *restrict a, const T *restrict b, int n, TRES res=TRES())
dot product
TinyVector< ST, 3 > PointType
Definition: SplineR2R.h:40
void resize(size_type n, size_type m)
Resize the container.
Definition: OhmmsMatrix.h:99
size_type cols() const
Definition: OhmmsMatrix.h:78
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
void FairDivideAligned(const int ntot, const int base, const int npart, const int me, int &first, int &last)
Partition ntot over npart and the size of each partition is a multiple of base size.
Definition: FairDivide.h:96
typename BsplineSet::ValueType TT
Definition: SplineR2R.h:43
omp_int_t omp_get_num_threads()
Definition: OpenMP.h:27
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
size_type rows() const
Definition: OhmmsMatrix.h:77
std::vector< SPOSet::PosType > kPoints
kpoints for each unique orbitals.
Definition: BsplineSet.h:52
void assign_v(int bc_sign, const vContainer_type &myV, ValueVector &psi, int first, int last) const
Definition: SplineR2R.cpp:152
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ evaluateValue()

void evaluateValue ( const ParticleSet P,
const int  iat,
ValueVector psi 
)
overridevirtual

evaluate the values of this single-particle orbital set

Parameters
Pcurrent ParticleSet
iatactive particle
psivalues of the SPO

Implements SPOSet.

Definition at line 166 of file SplineR2R.cpp.

References ParticleSet::activeR(), qmcplusplus::C2C::assign_v(), FairDivideAligned(), omp_get_num_threads(), and omp_get_thread_num().

167 {
168  const PointType& r = P.activeR(iat);
169  PointType ru;
170  int bc_sign = convertPos(r, ru);
171 
172 #pragma omp parallel
173  {
174  int first, last;
175  FairDivideAligned(psi.size(), getAlignment<ST>(), omp_get_num_threads(), omp_get_thread_num(), first, last);
176 
177  spline2::evaluate3d(SplineInst->getSplinePtr(), ru, myV, first, last);
178  assign_v(bc_sign, myV, psi, first, last);
179  }
180 }
int convertPos(const PointType &r, PointType &ru)
convert position in PrimLattice unit and return sign
Definition: SplineR2R.h:153
TinyVector< ST, 3 > PointType
Definition: SplineR2R.h:40
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
void FairDivideAligned(const int ntot, const int base, const int npart, const int me, int &first, int &last)
Partition ntot over npart and the size of each partition is a multiple of base size.
Definition: FairDivide.h:96
omp_int_t omp_get_num_threads()
Definition: OpenMP.h:27
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
void assign_v(int bc_sign, const vContainer_type &myV, ValueVector &psi, int first, int last) const
Definition: SplineR2R.cpp:152
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ evaluateVGH()

void evaluateVGH ( const ParticleSet P,
const int  iat,
ValueVector psi,
GradVector dpsi,
HessVector grad_grad_psi 
)
overridevirtual

evaluate the values, gradients and hessians of this single-particle orbital set

Parameters
Pcurrent ParticleSet
iatactive particle
psivalues of the SPO
dpsigradients of the SPO
grad_grad_psihessians of the SPO

Reimplemented from SPOSet.

Definition at line 373 of file SplineR2R.cpp.

References ParticleSet::activeR(), FairDivideAligned(), omp_get_num_threads(), and omp_get_thread_num().

378 {
379  const PointType& r = P.activeR(iat);
380  PointType ru;
381  int bc_sign = convertPos(r, ru);
382 
383 #pragma omp parallel
384  {
385  int first, last;
386  FairDivideAligned(psi.size(), getAlignment<ST>(), omp_get_num_threads(), omp_get_thread_num(), first, last);
387 
388  spline2::evaluate3d_vgh(SplineInst->getSplinePtr(), ru, myV, myG, myH, first, last);
389  assign_vgh(bc_sign, psi, dpsi, grad_grad_psi, first, last);
390  }
391 }
int convertPos(const PointType &r, PointType &ru)
convert position in PrimLattice unit and return sign
Definition: SplineR2R.h:153
TinyVector< ST, 3 > PointType
Definition: SplineR2R.h:40
hContainer_type myH
Definition: SplineR2R.h:75
void assign_vgh(int bc_sign, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi, int first, int last) const
Definition: SplineR2R.cpp:310
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
void FairDivideAligned(const int ntot, const int base, const int npart, const int me, int &first, int &last)
Partition ntot over npart and the size of each partition is a multiple of base size.
Definition: FairDivide.h:96
gContainer_type myG
Definition: SplineR2R.h:74
omp_int_t omp_get_num_threads()
Definition: OpenMP.h:27
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ evaluateVGHGH()

void evaluateVGHGH ( const ParticleSet P,
const int  iat,
ValueVector psi,
GradVector dpsi,
HessVector grad_grad_psi,
GGGVector grad_grad_grad_psi 
)
overridevirtual

evaluate the values, gradients, hessians, and grad hessians of this single-particle orbital set

Parameters
Pcurrent ParticleSet
iatactive particle
psivalues of the SPO
dpsigradients of the SPO
grad_grad_psihessians of the SPO
grad_grad_grad_psigrad hessians of the SPO

Reimplemented from SPOSet.

Definition at line 551 of file SplineR2R.cpp.

References ParticleSet::activeR(), FairDivideAligned(), omp_get_num_threads(), and omp_get_thread_num().

557 {
558  const PointType& r = P.activeR(iat);
559  PointType ru;
560  int bc_sign = convertPos(r, ru);
561 
562 #pragma omp parallel
563  {
564  int first, last;
565  FairDivideAligned(psi.size(), getAlignment<ST>(), omp_get_num_threads(), omp_get_thread_num(), first, last);
566 
567  spline2::evaluate3d_vghgh(SplineInst->getSplinePtr(), ru, myV, myG, myH, mygH, first, last);
568  assign_vghgh(bc_sign, psi, dpsi, grad_grad_psi, grad_grad_grad_psi, first, last);
569  }
570 }
int convertPos(const PointType &r, PointType &ru)
convert position in PrimLattice unit and return sign
Definition: SplineR2R.h:153
TinyVector< ST, 3 > PointType
Definition: SplineR2R.h:40
hContainer_type myH
Definition: SplineR2R.h:75
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
void FairDivideAligned(const int ntot, const int base, const int npart, const int me, int &first, int &last)
Partition ntot over npart and the size of each partition is a multiple of base size.
Definition: FairDivide.h:96
ghContainer_type mygH
Definition: SplineR2R.h:76
void assign_vghgh(int bc_sign, ValueVector &psi, GradVector &dpsi, HessVector &grad_grad_psi, GGGVector &grad_grad_grad_psi, int first=0, int last=-1) const
Definition: SplineR2R.cpp:394
gContainer_type myG
Definition: SplineR2R.h:74
omp_int_t omp_get_num_threads()
Definition: OpenMP.h:27
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ evaluateVGL()

void evaluateVGL ( const ParticleSet P,
const int  iat,
ValueVector psi,
GradVector dpsi,
ValueVector d2psi 
)
overridevirtual

evaluate the values, gradients and laplacians of this single-particle orbital set

Parameters
Pcurrent ParticleSet
iatactive particle
psivalues of the SPO
dpsigradients of the SPO
d2psilaplacians of the SPO

Implements SPOSet.

Definition at line 289 of file SplineR2R.cpp.

References ParticleSet::activeR(), qmcplusplus::C2C::assign_vgl(), FairDivideAligned(), omp_get_num_threads(), and omp_get_thread_num().

294 {
295  const PointType& r = P.activeR(iat);
296  PointType ru;
297  int bc_sign = convertPos(r, ru);
298 
299 #pragma omp parallel
300  {
301  int first, last;
302  FairDivideAligned(psi.size(), getAlignment<ST>(), omp_get_num_threads(), omp_get_thread_num(), first, last);
303 
304  spline2::evaluate3d_vgh(SplineInst->getSplinePtr(), ru, myV, myG, myH, first, last);
305  assign_vgl(bc_sign, psi, dpsi, d2psi, first, last);
306  }
307 }
void assign_vgl(int bc_sign, ValueVector &psi, GradVector &dpsi, ValueVector &d2psi, int first, int last) const
Definition: SplineR2R.cpp:226
int convertPos(const PointType &r, PointType &ru)
convert position in PrimLattice unit and return sign
Definition: SplineR2R.h:153
TinyVector< ST, 3 > PointType
Definition: SplineR2R.h:40
hContainer_type myH
Definition: SplineR2R.h:75
omp_int_t omp_get_thread_num()
Definition: OpenMP.h:25
void FairDivideAligned(const int ntot, const int base, const int npart, const int me, int &first, int &last)
Partition ntot over npart and the size of each partition is a multiple of base size.
Definition: FairDivide.h:96
gContainer_type myG
Definition: SplineR2R.h:74
omp_int_t omp_get_num_threads()
Definition: OpenMP.h:27
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ flush_zero()

void flush_zero ( )
inline

Definition at line 144 of file SplineR2R.h.

References SplineR2R< ST >::SplineInst.

144 { SplineInst->flush_zero(); }
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ gather_tables()

void gather_tables ( Communicate comm)
inline

Definition at line 122 of file SplineR2R.h.

References qmcplusplus::comm, FairDivideLow(), gatherv(), BsplineSet::kPoints, BsplineSet::offset, Communicate::size(), and SplineR2R< ST >::SplineInst.

123  {
124  if (comm->size() == 1)
125  return;
126  const int Nbands = kPoints.size();
127  const int Nbandgroups = comm->size();
128  offset.resize(Nbandgroups + 1, 0);
129  FairDivideLow(Nbands, Nbandgroups, offset);
130  gatherv(comm, SplineInst->getSplinePtr(), SplineInst->getSplinePtr()->z_stride, offset);
131  }
int size() const
return the number of tasks
Definition: Communicate.h:118
void FairDivideLow(int ntot, int npart, IV &adist)
partition ntot elements among npart
Definition: FairDivide.h:114
std::vector< int > offset
band offsets used for communication
Definition: BsplineSet.h:56
void gatherv(T *sb, T *rb, int n, IT &counts, IT &displ, int dest)
std::vector< SPOSet::PosType > kPoints
kpoints for each unique orbitals.
Definition: BsplineSet.h:52
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ getClassName()

virtual std::string getClassName ( ) const
inlineoverridevirtual

return class name

Implements SPOSet.

Definition at line 82 of file SplineR2R.h.

82 { return "SplineR2R"; }

◆ getKeyword()

virtual std::string getKeyword ( ) const
inlineoverridevirtual

Implements BsplineSet.

Definition at line 83 of file SplineR2R.h.

83 { return "SplineR2R"; }

◆ isComplex()

bool isComplex ( ) const
inlineoverridevirtual

Implements BsplineSet.

Definition at line 84 of file SplineR2R.h.

84 { return false; };

◆ isRotationSupported()

bool isRotationSupported ( ) const
inlineoverridevirtual

return true if this SPOSet can be wrappered by RotatedSPO

Reimplemented from SPOSet.

Definition at line 85 of file SplineR2R.h.

85 { return true; }

◆ makeClone()

std::unique_ptr<SPOSet> makeClone ( ) const
inlineoverridevirtual

make a clone of itself every derived class must implement this to have threading working correctly.

Implements BsplineSet.

Definition at line 87 of file SplineR2R.h.

87 { return std::make_unique<SplineR2R>(*this); }

◆ read_splines()

bool read_splines ( hdf_archive h5f)

Definition at line 39 of file SplineR2R.cpp.

References hdf_archive::readEntry().

40 {
41  std::ostringstream o;
42  o << "spline_" << MyIndex;
43  einspline_engine<SplineType> bigtable(SplineInst->getSplinePtr());
44  return h5f.readEntry(bigtable, o.str().c_str()); //"spline_0");
45 }
size_t MyIndex
Index of this adoptor, when multiple adoptors are used for NUMA or distributed cases.
Definition: BsplineSet.h:39
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ resizeStorage()

void resizeStorage ( size_t  n,
size_t  nvals 
)
inline

Definition at line 107 of file SplineR2R.h.

References BsplineSet::HalfG, BsplineSet::init_base(), SplineR2R< ST >::IsGamma, SplineR2R< ST >::myG, SplineR2R< ST >::mygH, SplineR2R< ST >::myH, SplineR2R< ST >::myL, SplineR2R< ST >::myV, qmcplusplus::n, VectorSoaContainer< T, D, Alloc >::resize(), and Vector< T, Alloc >::resize().

108  {
109  init_base(n);
110  const size_t npad = getAlignedSize<ST>(n);
111  myV.resize(npad);
112  myG.resize(npad);
113  myL.resize(npad);
114  myH.resize(npad);
115  mygH.resize(npad);
116 
117  IsGamma = ((HalfG[0] == 0) && (HalfG[1] == 0) && (HalfG[2] == 0));
118  }
void resize(size_type n, Type_t val=Type_t())
Resize the container.
Definition: OhmmsVector.h:166
hContainer_type myH
Definition: SplineR2R.h:75
ghContainer_type mygH
Definition: SplineR2R.h:76
vContainer_type myL
Definition: SplineR2R.h:73
gContainer_type myG
Definition: SplineR2R.h:74
vContainer_type myV
intermediate result vectors
Definition: SplineR2R.h:72
void init_base(int n)
Definition: BsplineSet.h:66
TinyVector< int, D > HalfG
sign bits at the G/2 boundaries
Definition: BsplineSet.h:45
void resize(size_type n)
resize myData

◆ set_spline()

void set_spline ( SingleSplineType spline_r,
SingleSplineType spline_i,
int  twist,
int  ispline,
int  level 
)
inline

Definition at line 29 of file SplineR2R.cpp.

34 {
35  SplineInst->copy_spline(spline_r, ispline);
36 }
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ storeParamsBeforeRotation()

void storeParamsBeforeRotation ( )
overridevirtual

Store an original copy of the spline coefficients for orbital rotation.

Reimplemented from SPOSet.

Definition at line 57 of file SplineR2R.cpp.

References qmcplusplus::syclBLAS::copy_n().

58 {
59  const auto spline_ptr = SplineInst->getSplinePtr();
60  const auto coefs_tot_size = spline_ptr->coefs_size;
61  coef_copy_ = std::make_shared<std::vector<ST>>(coefs_tot_size);
62 
63  std::copy_n(spline_ptr->coefs, coefs_tot_size, coef_copy_->begin());
64 }
std::shared_ptr< std::vector< ST > > coef_copy_
Copy of original splines for orbital rotation.
Definition: SplineR2R.h:62
sycl::event copy_n(sycl::queue &aq, const T1 *restrict VA, size_t array_size, T2 *restrict VC, const std::vector< sycl::event > &events)
Definition: syclBLAS.cpp:548
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

◆ write_splines()

bool write_splines ( hdf_archive h5f)

Definition at line 48 of file SplineR2R.cpp.

References hdf_archive::writeEntry().

49 {
50  std::ostringstream o;
51  o << "spline_" << MyIndex;
52  einspline_engine<SplineType> bigtable(SplineInst->getSplinePtr());
53  return h5f.writeEntry(bigtable, o.str().c_str()); //"spline_0");
54 }
size_t MyIndex
Index of this adoptor, when multiple adoptors are used for NUMA or distributed cases.
Definition: BsplineSet.h:39
std::shared_ptr< MultiBspline< ST > > SplineInst
multi bspline set
Definition: SplineR2R.h:59

Friends And Related Function Documentation

◆ BsplineReader

friend struct BsplineReader
friend

Definition at line 216 of file SplineR2R.h.

◆ SplineSetReader

friend class SplineSetReader
friend

Definition at line 215 of file SplineR2R.h.

Member Data Documentation

◆ coef_copy_

std::shared_ptr<std::vector<ST> > coef_copy_
private

Copy of original splines for orbital rotation.

Definition at line 62 of file SplineR2R.h.

◆ GGt

Tensor<ST, 3> GGt
private

$GGt=G^t G $, transformation for tensor in LatticeUnit to CartesianUnit, e.g. Hessian

Definition at line 57 of file SplineR2R.h.

Referenced by SplineR2R< ST >::create_spline().

◆ IsGamma

bool IsGamma
private

Definition at line 55 of file SplineR2R.h.

Referenced by SplineR2R< ST >::resizeStorage().

◆ myG

gContainer_type myG
protected

Definition at line 74 of file SplineR2R.h.

Referenced by SplineR2R< ST >::resizeStorage().

◆ mygH

ghContainer_type mygH
protected

Definition at line 76 of file SplineR2R.h.

Referenced by SplineR2R< ST >::resizeStorage().

◆ myH

hContainer_type myH
protected

Definition at line 75 of file SplineR2R.h.

Referenced by SplineR2R< ST >::resizeStorage().

◆ myL

vContainer_type myL
protected

Definition at line 73 of file SplineR2R.h.

Referenced by SplineR2R< ST >::resizeStorage().

◆ myV

vContainer_type myV
protected

intermediate result vectors

Definition at line 72 of file SplineR2R.h.

Referenced by SplineR2R< ST >::create_spline(), and SplineR2R< ST >::resizeStorage().

◆ PrimLattice

CrystalLattice<ST, 3> PrimLattice
protected

primitive cell

Definition at line 70 of file SplineR2R.h.

Referenced by SplineR2R< ST >::convertPos(), and SplineR2R< ST >::create_spline().

◆ ratios_private

Matrix<TT> ratios_private
private

thread private ratios for reduction when using nested threading, numVP x numThread

Definition at line 65 of file SplineR2R.h.

◆ SplineInst

std::shared_ptr<MultiBspline<ST> > SplineInst
private

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