20#ifndef DFTFE_ONCVCLASS_H
21#define DFTFE_ONCVCLASS_H
26 template <
typename ValueType, dftfe::utils::MemorySpace memorySpace>
31 const std::string &scratchFolderName,
32 const std::set<dftfe::uInt> &atomTypes,
33 const bool floatingNuclearCharges,
35 const std::map<dftfe::uInt, dftfe::uInt> &atomAttributes,
36 const bool reproducibleOutput,
39 const bool memOptMode);
60 FEBasisOperations<ValueType, double, dftfe::utils::MemorySpace::HOST>>
61 basisOperationsHostPtr,
62#
if defined(DFTFE_WITH_DEVICE)
67 basisOperationsDevicePtr,
72#
if defined(DFTFE_WITH_DEVICE)
83 const std::vector<std::vector<double>> &atomLocations,
85 const bool singlePrecNonLocalOperator,
86 const bool computeSphericalFnTimesXNonLocalOperator =
true);
106 const std::vector<std::vector<double>> &atomLocations,
107 const std::vector<dftfe::Int> &imageIds,
108 const std::vector<std::vector<double>> &periodicCoords,
109 const std::vector<double> &kPointWeights,
110 const std::vector<double> &kPointCoordinates,
111 const bool updateNonlocalSparsity);
116 const std::vector<std::vector<double>> &atomLocations,
117 const std::vector<dftfe::Int> &imageIds,
118 const std::vector<std::vector<double>> &periodicCoords,
119 const std::vector<double> &kPointWeights,
120 const std::vector<double> &kPointCoordinates,
121 const bool updateNonlocalSparsity,
122 const std::map<
dftfe::uInt, std::vector<dftfe::Int>> &sparsityPattern,
123 const std::vector<std::vector<dealii::CellId>>
124 &elementIdsInAtomCompactSupport,
125 const std::vector<std::vector<dftfe::uInt>>
126 &elementIndexesInAtomCompactSupport,
127 const std::vector<dftfe::uInt> &atomIdsInCurrentProcess,
140 std::vector<double> &Val);
151 std::vector<double> &Val);
183 const std::shared_ptr<
210 const std::vector<dftfe::Int> &imageIds,
211 const std::vector<std::vector<double>> &periodicCoords,
212 std::vector<dftfe::uInt> &imageIdsTemp,
213 std::vector<double> &imageCoordsTemp);
230#if defined(DFTFE_WITH_DEVICE)
233 d_BLASWrapperDevicePtr;
236 std::map<dftfe::uInt, std::vector<double>>
246 std::vector<std::shared_ptr<AtomCenteredSphericalFunctionBase>>
248 std::shared_ptr<AtomCenteredSphericalFunctionContainer>
250 std::map<std::pair<dftfe::uInt, dftfe::uInt>,
251 std::shared_ptr<AtomCenteredSphericalFunctionBase>>
271 FEBasisOperations<ValueType, double, dftfe::utils::MemorySpace::HOST>>
273#if defined(DFTFE_WITH_DEVICE)
276 FEBasisOperations<ValueType, double, dftfe::utils::MemorySpace::DEVICE>>
277 d_BasisOperatorDevicePtr;
294 std::shared_ptr<AtomicCenteredNonLocalOperator<ValueType, memorySpace>>
303 std::vector<std::shared_ptr<AtomCenteredSphericalFunctionBase>>
306 std::map<dftfe::uInt, std::shared_ptr<AtomCenteredSphericalFunctionBase>>>
309 std::map<dftfe::uInt, std::shared_ptr<AtomCenteredSphericalFunctionBase>>>
312 std::map<dftfe::uInt, std::shared_ptr<AtomCenteredSphericalFunctionBase>>>
Definition AtomicCenteredNonLocalOperator.h:58
Definition FEBasisOperations.h:85
Definition excManager.h:28
Definition BLASWrapper.h:35
std::map< dftfe::uInt, dftfe::uInt > d_atomTypeAtributes
Definition oncvClass.h:317
std::vector< std::map< dftfe::uInt, std::shared_ptr< AtomCenteredSphericalFunctionBase > > > d_atomicCoreDensityVector
Definition oncvClass.h:313
double getRadialLocalPseudo(dftfe::uInt Znum, double rad)
std::vector< std::map< dftfe::uInt, std::shared_ptr< AtomCenteredSphericalFunctionBase > > > d_atomicValenceDensityVector
Definition oncvClass.h:310
void createAtomCenteredSphericalFunctionsForDensities()
Creating Density splines for all atomTypes.
const dftfe::utils::MemoryStorage< ValueType, memorySpace > & getCouplingMatrix(CouplingType couplingtype=CouplingType::HamiltonianEntries)
dftfe::uInt d_nuclearChargeQuadratureIdElectro
Definition oncvClass.h:264
dftfe::uInt d_densityQuadratureId
Definition oncvClass.h:262
double getRadialCoreDensity(dftfe::uInt Znum, double rad)
dftfe::Int d_verbosity
Definition oncvClass.h:283
void computeNonlocalPseudoPotentialConstants()
oncvClass(const MPI_Comm &mpi_comm_parent, const std::string &scratchFolderName, const std::set< dftfe::uInt > &atomTypes, const bool floatingNuclearCharges, const dftfe::uInt nOMPThreads, const std::map< dftfe::uInt, dftfe::uInt > &atomAttributes, const bool reproducibleOutput, const dftfe::Int verbosity, const bool useDevice, const bool memOptMode)
const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueType, memorySpace > > getNonLocalOperator()
void setImageCoordinates(const std::vector< std::vector< double > > &atomLocations, const std::vector< dftfe::Int > &imageIds, const std::vector< std::vector< double > > &periodicCoords, std::vector< dftfe::uInt > &imageIdsTemp, std::vector< double > &imageCoordsTemp)
Converts the periodic image data structure to relevant form for the container class.
std::vector< std::vector< double > > d_nonLocalPseudoPotentialConstants
Definition oncvClass.h:235
void getRadialCoreDensity(dftfe::uInt Znum, double rad, std::vector< double > &Val)
std::shared_ptr< excManager< memorySpace > > d_excManagerPtr
Definition oncvClass.h:268
dftfe::uInt d_localContributionQuadratureId
Definition oncvClass.h:263
std::set< dftfe::uInt > d_atomTypes
Definition oncvClass.h:285
std::map< dftfe::uInt, std::vector< dftfe::uInt > > d_atomTypesList
Definition oncvClass.h:286
std::shared_ptr< AtomicCenteredNonLocalOperator< ValueType, memorySpace > > d_nonLocalOperator
Definition oncvClass.h:295
void createAtomCenteredSphericalFunctionsForLocalPotential()
bool d_HamiltonianCouplingMatrixSinglePrecEntriesUpdated
Definition oncvClass.h:245
dftfe::uInt getTotalNumberOfSphericalFunctionsForAtomId(dftfe::uInt atomId)
std::map< dftfe::uInt, std::vector< double > > d_atomicNonLocalPseudoPotentialConstants
Definition oncvClass.h:237
std::vector< std::vector< double > > d_atomLocations
Definition oncvClass.h:284
std::vector< std::vector< double > > d_imagePositions
Definition oncvClass.h:289
dftfe::utils::MemoryStorage< ValueType, memorySpace > d_couplingMatrixEntries
Definition oncvClass.h:238
const MPI_Comm d_mpiCommParent
Definition oncvClass.h:255
bool d_singlePrecNonLocalOperator
Definition oncvClass.h:282
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > d_couplingMatrixEntriesSinglePrec
Definition oncvClass.h:242
void initialiseNonLocalContribution(const std::vector< std::vector< double > > &atomLocations, const std::vector< dftfe::Int > &imageIds, const std::vector< std::vector< double > > &periodicCoords, const std::vector< double > &kPointWeights, const std::vector< double > &kPointCoordinates, const bool updateNonlocalSparsity)
Initialises all the data members with addresses/values to/of dftClass.
std::shared_ptr< AtomCenteredSphericalFunctionContainer > d_atomicProjectorFnsContainer
Definition oncvClass.h:249
std::vector< std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_atomicProjectorFnsVector
Definition oncvClass.h:304
double getRadialValenceDensity(dftfe::uInt Znum, double rad)
bool d_HamiltonianCouplingMatrixEntriesUpdated
Definition oncvClass.h:244
bool d_floatingNuclearCharges
Definition oncvClass.h:281
const std::shared_ptr< AtomicCenteredNonLocalOperator< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > > getNonLocalOperatorSinglePrec()
std::shared_ptr< AtomicCenteredNonLocalOperator< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > > d_nonLocalOperatorSinglePrec
Definition oncvClass.h:300
void getRadialValenceDensity(dftfe::uInt Znum, double rad, std::vector< double > &Val)
bool coreNuclearDensityPresent(dftfe::uInt Znum)
std::vector< std::map< dftfe::uInt, std::shared_ptr< AtomCenteredSphericalFunctionBase > > > d_atomicLocalPotVector
Definition oncvClass.h:307
dftfe::uInt d_numEigenValues
Definition oncvClass.h:290
std::vector< std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_atomicWaveFnsVector
Definition oncvClass.h:247
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > d_BLASWrapperHostPtr
Definition oncvClass.h:229
void initLocalPotential()
Initialises local potential.
dealii::ConditionalOStream pcout
Definition oncvClass.h:259
void initialise(std::shared_ptr< dftfe::basis::FEBasisOperations< ValueType, double, dftfe::utils::MemorySpace::HOST > > basisOperationsHostPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperPtrHost, dftfe::uInt densityQuadratureId, dftfe::uInt localContributionQuadratureId, dftfe::uInt sparsityPatternQuadratureId, dftfe::uInt nlpspQuadratureId, dftfe::uInt densityQuadratureIdElectro, std::shared_ptr< excManager< memorySpace > > excFunctionalPtr, const std::vector< std::vector< double > > &atomLocations, dftfe::uInt numEigenValues, const bool singlePrecNonLocalOperator, const bool computeSphericalFnTimesXNonLocalOperator=true)
Initialises all the data members with addresses/values to/of dftClass.
void initialiseNonLocalContribution(const std::vector< std::vector< double > > &atomLocations, const std::vector< dftfe::Int > &imageIds, const std::vector< std::vector< double > > &periodicCoords, const std::vector< double > &kPointWeights, const std::vector< double > &kPointCoordinates, const bool updateNonlocalSparsity, const std::map< dftfe::uInt, std::vector< dftfe::Int > > &sparsityPattern, const std::vector< std::vector< dealii::CellId > > &elementIdsInAtomCompactSupport, const std::vector< std::vector< dftfe::uInt > > &elementIndexesInAtomCompactSupport, const std::vector< dftfe::uInt > &atomIdsInCurrentProcess, dftfe::uInt numberElements)
dftfe::uInt d_densityQuadratureIdElectro
Definition oncvClass.h:265
void createAtomCenteredSphericalFunctionsForProjectors()
bool d_reproducible_output
Definition oncvClass.h:314
dftfe::uInt d_nlpspQuadratureId
Definition oncvClass.h:267
bool d_memoryOptMode
Definition oncvClass.h:261
bool d_useDevice
Definition oncvClass.h:260
dftfe::uInt getAtomIdInCurrentProcessor(dftfe::uInt iAtom)
dftfe::uInt d_nOMPThreads
Definition oncvClass.h:291
std::shared_ptr< dftfe::basis::FEBasisOperations< ValueType, double, dftfe::utils::MemorySpace::HOST > > d_BasisOperatorHostPtr
Definition oncvClass.h:272
double getRmaxValenceDensity(dftfe::uInt Znum)
std::map< std::pair< dftfe::uInt, dftfe::uInt >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_atomicProjectorFnsMap
Definition oncvClass.h:252
dftfe::uInt getTotalNumberOfAtomsInCurrentProcessor()
std::map< dftfe::uInt, bool > d_atomTypeCoreFlagMap
Definition oncvClass.h:280
double getRmaxLocalPot(dftfe::uInt Znum)
double getRmaxCoreDensity(dftfe::uInt Znum)
std::vector< dftfe::Int > d_imageIds
Definition oncvClass.h:288
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > & getCouplingMatrixSinglePrec(CouplingType couplingtype=CouplingType::HamiltonianEntries)
dftfe::uInt d_sparsityPatternQuadratureId
Definition oncvClass.h:266
std::string d_dftfeScratchFolderName
Definition oncvClass.h:287
const dftfe::uInt d_this_mpi_process
Definition oncvClass.h:256
Definition pseudopotentialBaseClass.h:60
Definition MemoryStorage.h:33
Definition FEBasisOperations.h:30
@ DEVICE
Definition MemorySpaceType.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
CouplingType
Definition pseudopotentialBaseClass.h:48
@ HamiltonianEntries
Definition pseudopotentialBaseClass.h:49
std::uint32_t uInt
Definition TypeConfig.h:10
std::int32_t Int
Definition TypeConfig.h:11
T type
Definition dftfeDataTypes.h:113