DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
oncvClass.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17// @author Kartick Ramakrishnan, Sambit Das
18//
19
20#ifndef DFTFE_ONCVCLASS_H
21#define DFTFE_ONCVCLASS_H
22
23#include "vector"
24#include "map"
31#include <memory>
32#include <MemorySpaceType.h>
33#include <headers.h>
34#include <TypeConfig.h>
35#include <dftUtils.h>
36#include "FEBasisOperations.h"
37#include <BLASWrapper.h>
38#include <xc.h>
39#include <excManager.h>
40#ifdef _OPENMP
41# include <omp.h>
42#else
43# define omp_get_thread_num() 0
44#endif
45namespace dftfe
46{
47 template <typename ValueType, dftfe::utils::MemorySpace memorySpace>
49 {
50 public:
51 oncvClass(const MPI_Comm & mpi_comm_parent,
52 const std::string & scratchFolderName,
53 const std::set<unsigned int> &atomTypes,
54 const bool floatingNuclearCharges,
55 const unsigned int nOMPThreads,
56 const std::map<unsigned int, unsigned int> &atomAttributes,
57 const bool reproducibleOutput,
58 const int verbosity,
59 const bool useDevice,
60 const bool memOptMode);
61 /**
62 * @brief Initialises all the data members with addresses/values to/of dftClass.
63 * @param[in] densityQuadratureId quadratureId for density.
64 * @param[in] localContributionQuadratureId quadratureId for local/zero
65 * potential
66 * @param[in] nuclearChargeQuadratureIdElectro quadratureId for nuclear
67 * charges
68 * @param[in] densityQuadratureIdElectro quadratureId for density in
69 * Electrostatics mesh
70 * @param[in] excFunctionalPtr address XC functional pointer
71 * @param[in] numEigenValues number of eigenvalues
72 * @param[in] atomLocations atomic Coordinates
73 * @param[in] imageIds image IDs of periodic cell
74 * @param[in] periodicCoords coordinates of image atoms
75 */
76
77 void
79 std::shared_ptr<
81 FEBasisOperations<ValueType, double, dftfe::utils::MemorySpace::HOST>>
82 basisOperationsHostPtr,
83#if defined(DFTFE_WITH_DEVICE)
84 std::shared_ptr<
86 double,
88 basisOperationsDevicePtr,
89#endif
90 std::shared_ptr<
92 BLASWrapperPtrHost,
93#if defined(DFTFE_WITH_DEVICE)
94 std::shared_ptr<
96 BLASWrapperPtrDevice,
97#endif
98 unsigned int densityQuadratureId,
99 unsigned int localContributionQuadratureId,
100 unsigned int sparsityPatternQuadratureId,
101 unsigned int nlpspQuadratureId,
102 unsigned int densityQuadratureIdElectro,
103 std::shared_ptr<excManager<memorySpace>> excFunctionalPtr,
104 const std::vector<std::vector<double>> & atomLocations,
105 unsigned int numEigenValues,
106 const bool singlePrecNonLocalOperator,
107 const bool computeSphericalFnTimesXNonLocalOperator = true);
108
109 /**
110 * @brief Initialises all the data members with addresses/values to/of dftClass.
111 * @param[in] densityQuadratureId quadratureId for density.
112 * @param[in] localContributionQuadratureId quadratureId for local/zero
113 * potential
114 * @param[in] nuclearChargeQuadratureIdElectro quadratureId for nuclear
115 * charges
116 * @param[in] densityQuadratureIdElectro quadratureId for density in
117 * Electrostatics mesh
118 * @param[in] bQuadValuesAllAtoms address of nuclear charge field
119 * @param[in] excFunctionalPtr address XC functional pointer
120 * @param[in] numEigenValues number of eigenvalues
121 * @param[in] atomLocations atomic Coordinates
122 * @param[in] imageIds image IDs of periodic cell
123 * @param[in] periodicCoords coordinates of image atoms
124 */
125 void
127 const std::vector<std::vector<double>> &atomLocations,
128 const std::vector<int> & imageIds,
129 const std::vector<std::vector<double>> &periodicCoords,
130 const std::vector<double> & kPointWeights,
131 const std::vector<double> & kPointCoordinates,
132 const bool updateNonlocalSparsity);
133
134
135 void
137 const std::vector<std::vector<double>> & atomLocations,
138 const std::vector<int> & imageIds,
139 const std::vector<std::vector<double>> & periodicCoords,
140 const std::vector<double> & kPointWeights,
141 const std::vector<double> & kPointCoordinates,
142 const bool updateNonlocalSparsity,
143 const std::map<unsigned int, std::vector<int>> &sparsityPattern,
144 const std::vector<std::vector<dealii::CellId>>
145 &elementIdsInAtomCompactSupport,
146 const std::vector<std::vector<unsigned int>>
147 & elementIndexesInAtomCompactSupport,
148 const std::vector<unsigned int> &atomIdsInCurrentProcess,
149 unsigned int numberElements);
150
151
152 /**
153 * @brief Initialises local potential
154 */
155 void
157
158 void
159 getRadialValenceDensity(unsigned int Znum,
160 double rad,
161 std::vector<double> &Val);
162
163 double
164 getRadialValenceDensity(unsigned int Znum, double rad);
165
166 double
167 getRmaxValenceDensity(unsigned int Znum);
168
169 void
170 getRadialCoreDensity(unsigned int Znum,
171 double rad,
172 std::vector<double> &Val);
173
174 double
175 getRadialCoreDensity(unsigned int Znum, double rad);
176
177 double
178 getRmaxCoreDensity(unsigned int Znum);
179
180 double
181 getRadialLocalPseudo(unsigned int Znum, double rad);
182
183 double
184 getRmaxLocalPot(unsigned int Znum);
185
186 bool
187 coreNuclearDensityPresent(unsigned int Znum);
188 // Returns the number of Projectors for the given atomID in cooridnates List
189 unsigned int
191 // Returns the Total Number of atoms with support in the processor
192 unsigned int
194 // Returns the atomID in coordinates list for the iAtom index.
195 unsigned int
196 getAtomIdInCurrentProcessor(unsigned int iAtom);
197
198
201
202
203 const std::shared_ptr<
206
209 memorySpace> &
211
212
213 const std::shared_ptr<AtomicCenteredNonLocalOperator<
215 memorySpace>>
217
218 private:
219 /**
220 * @brief Converts the periodic image data structure to relevant form for the container class
221 * @param[in] atomLocations atomic Coordinates
222 * @param[in] imageIds image IDs of periodic cell
223 * @param[in] periodicCoords coordinates of image atoms
224 * @param[out] imageIdsTemp image IDs of periodic cell
225 * @param[out] imageCoordsTemp coordinates of image atoms
226 */
227 void
228 setImageCoordinates(const std::vector<std::vector<double>> &atomLocations,
229 const std::vector<int> & imageIds,
230 const std::vector<std::vector<double>> &periodicCoords,
231 std::vector<unsigned int> & imageIdsTemp,
232 std::vector<double> &imageCoordsTemp);
233 /**
234 * @brief Creating Density splines for all atomTypes
235 */
236 void
238
239 void
241 void
243 void
245
246 std::shared_ptr<
249#if defined(DFTFE_WITH_DEVICE)
250 std::shared_ptr<
252 d_BLASWrapperDevicePtr;
253#endif
254 std::vector<std::vector<double>> d_nonLocalPseudoPotentialConstants;
255 std::map<unsigned int, std::vector<double>>
260 memorySpace>
262
265 std::vector<std::shared_ptr<AtomCenteredSphericalFunctionBase>>
267 std::shared_ptr<AtomCenteredSphericalFunctionContainer>
269 std::map<std::pair<unsigned int, unsigned int>,
270 std::shared_ptr<AtomCenteredSphericalFunctionBase>>
272
273 // parallel communication objects
274 const MPI_Comm d_mpiCommParent;
275 const unsigned int d_this_mpi_process;
276
277 // conditional stream object
278 dealii::ConditionalOStream pcout;
287 std::shared_ptr<excManager<memorySpace>> d_excManagerPtr;
288 std::shared_ptr<
289 dftfe::basis::
290 FEBasisOperations<ValueType, double, dftfe::utils::MemorySpace::HOST>>
292#if defined(DFTFE_WITH_DEVICE)
293 std::shared_ptr<
294 dftfe::basis::
295 FEBasisOperations<ValueType, double, dftfe::utils::MemorySpace::DEVICE>>
296 d_BasisOperatorDevicePtr;
297#endif
298
299 std::map<unsigned int, bool> d_atomTypeCoreFlagMap;
303 std::vector<std::vector<double>> d_atomLocations;
304 std::set<unsigned int> d_atomTypes;
305 std::map<unsigned int, std::vector<unsigned int>> d_atomTypesList;
307 std::vector<int> d_imageIds;
308 std::vector<std::vector<double>> d_imagePositions;
309 unsigned int d_numEigenValues;
310 unsigned int d_nOMPThreads;
311
312 // Creating Object for Atom Centerd Nonlocal Operator
313 std::shared_ptr<AtomicCenteredNonLocalOperator<ValueType, memorySpace>>
315
316 std::shared_ptr<AtomicCenteredNonLocalOperator<
318 memorySpace>>
320
321
322 std::vector<std::shared_ptr<AtomCenteredSphericalFunctionBase>>
324 std::vector<std::map<unsigned int,
325 std::shared_ptr<AtomCenteredSphericalFunctionBase>>>
327 std::vector<std::map<unsigned int,
328 std::shared_ptr<AtomCenteredSphericalFunctionBase>>>
330 std::vector<std::map<unsigned int,
331 std::shared_ptr<AtomCenteredSphericalFunctionBase>>>
334 /// FIXME: eventually it should be a map of atomic number to struct-
335 /// {valence number, mesh input etc}
336 std::map<unsigned int, unsigned int> d_atomTypeAtributes;
337
338
339
340 }; // end of class
341} // end of namespace dftfe
342#endif // DFTFE_PSEUDOPOTENTIALBASE_H
Definition AtomicCenteredNonLocalOperator.h:58
Definition FEBasisOperations.h:84
Definition excManager.h:28
Definition BLASWrapper.h:35
std::map< unsigned int, unsigned int > d_atomTypeAtributes
Definition oncvClass.h:336
int d_verbosity
Definition oncvClass.h:302
const dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > & getCouplingMatrixSinglePrec()
void createAtomCenteredSphericalFunctionsForDensities()
Creating Density splines for all atomTypes.
void computeNonlocalPseudoPotentialConstants()
oncvClass(const MPI_Comm &mpi_comm_parent, const std::string &scratchFolderName, const std::set< unsigned int > &atomTypes, const bool floatingNuclearCharges, const unsigned int nOMPThreads, const std::map< unsigned int, unsigned int > &atomAttributes, const bool reproducibleOutput, const int verbosity, const bool useDevice, const bool memOptMode)
const std::shared_ptr< AtomicCenteredNonLocalOperator< ValueType, memorySpace > > getNonLocalOperator()
std::vector< std::vector< double > > d_nonLocalPseudoPotentialConstants
Definition oncvClass.h:254
std::shared_ptr< excManager< memorySpace > > d_excManagerPtr
Definition oncvClass.h:287
unsigned int getAtomIdInCurrentProcessor(unsigned int iAtom)
const unsigned int d_this_mpi_process
Definition oncvClass.h:275
std::map< unsigned int, bool > d_atomTypeCoreFlagMap
Definition oncvClass.h:299
void initialiseNonLocalContribution(const std::vector< std::vector< double > > &atomLocations, const std::vector< 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< AtomicCenteredNonLocalOperator< ValueType, memorySpace > > d_nonLocalOperator
Definition oncvClass.h:314
unsigned int getTotalNumberOfSphericalFunctionsForAtomId(unsigned int atomId)
void createAtomCenteredSphericalFunctionsForLocalPotential()
bool d_HamiltonianCouplingMatrixSinglePrecEntriesUpdated
Definition oncvClass.h:264
unsigned int d_nlpspQuadratureId
Definition oncvClass.h:286
std::vector< std::vector< double > > d_atomLocations
Definition oncvClass.h:303
void setImageCoordinates(const std::vector< std::vector< double > > &atomLocations, const std::vector< int > &imageIds, const std::vector< std::vector< double > > &periodicCoords, std::vector< unsigned int > &imageIdsTemp, std::vector< double > &imageCoordsTemp)
Converts the periodic image data structure to relevant form for the container class.
double getRmaxLocalPot(unsigned int Znum)
double getRadialValenceDensity(unsigned int Znum, double rad)
unsigned int d_localContributionQuadratureId
Definition oncvClass.h:282
std::vector< std::vector< double > > d_imagePositions
Definition oncvClass.h:308
dftfe::utils::MemoryStorage< ValueType, memorySpace > d_couplingMatrixEntries
Definition oncvClass.h:257
void getRadialValenceDensity(unsigned int Znum, double rad, std::vector< double > &Val)
const MPI_Comm d_mpiCommParent
Definition oncvClass.h:274
bool d_singlePrecNonLocalOperator
Definition oncvClass.h:301
dftfe::utils::MemoryStorage< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > d_couplingMatrixEntriesSinglePrec
Definition oncvClass.h:261
unsigned int d_nuclearChargeQuadratureIdElectro
Definition oncvClass.h:283
std::shared_ptr< AtomCenteredSphericalFunctionContainer > d_atomicProjectorFnsContainer
Definition oncvClass.h:268
double getRmaxValenceDensity(unsigned int Znum)
std::vector< std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_atomicProjectorFnsVector
Definition oncvClass.h:323
unsigned int d_sparsityPatternQuadratureId
Definition oncvClass.h:285
bool d_HamiltonianCouplingMatrixEntriesUpdated
Definition oncvClass.h:263
bool d_floatingNuclearCharges
Definition oncvClass.h:300
bool coreNuclearDensityPresent(unsigned int Znum)
std::set< unsigned int > d_atomTypes
Definition oncvClass.h:304
const std::shared_ptr< AtomicCenteredNonLocalOperator< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > > getNonLocalOperatorSinglePrec()
std::vector< std::map< unsigned int, std::shared_ptr< AtomCenteredSphericalFunctionBase > > > d_atomicCoreDensityVector
Definition oncvClass.h:332
unsigned int d_densityQuadratureId
Definition oncvClass.h:281
void initialiseNonLocalContribution(const std::vector< std::vector< double > > &atomLocations, const std::vector< 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< unsigned int, std::vector< int > > &sparsityPattern, const std::vector< std::vector< dealii::CellId > > &elementIdsInAtomCompactSupport, const std::vector< std::vector< unsigned int > > &elementIndexesInAtomCompactSupport, const std::vector< unsigned int > &atomIdsInCurrentProcess, unsigned int numberElements)
std::shared_ptr< AtomicCenteredNonLocalOperator< typename dftfe::dataTypes::singlePrecType< ValueType >::type, memorySpace > > d_nonLocalOperatorSinglePrec
Definition oncvClass.h:319
unsigned int getTotalNumberOfAtomsInCurrentProcessor()
std::vector< int > d_imageIds
Definition oncvClass.h:307
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, unsigned int densityQuadratureId, unsigned int localContributionQuadratureId, unsigned int sparsityPatternQuadratureId, unsigned int nlpspQuadratureId, unsigned int densityQuadratureIdElectro, std::shared_ptr< excManager< memorySpace > > excFunctionalPtr, const std::vector< std::vector< double > > &atomLocations, unsigned int numEigenValues, const bool singlePrecNonLocalOperator, const bool computeSphericalFnTimesXNonLocalOperator=true)
Initialises all the data members with addresses/values to/of dftClass.
const dftfe::utils::MemoryStorage< ValueType, memorySpace > & getCouplingMatrix()
std::vector< std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_atomicWaveFnsVector
Definition oncvClass.h:266
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > d_BLASWrapperHostPtr
Definition oncvClass.h:248
void initLocalPotential()
Initialises local potential.
dealii::ConditionalOStream pcout
Definition oncvClass.h:278
unsigned int d_densityQuadratureIdElectro
Definition oncvClass.h:284
std::map< unsigned int, std::vector< unsigned int > > d_atomTypesList
Definition oncvClass.h:305
double getRadialCoreDensity(unsigned int Znum, double rad)
std::map< unsigned int, std::vector< double > > d_atomicNonLocalPseudoPotentialConstants
Definition oncvClass.h:256
unsigned int d_numEigenValues
Definition oncvClass.h:309
void createAtomCenteredSphericalFunctionsForProjectors()
std::vector< std::map< unsigned int, std::shared_ptr< AtomCenteredSphericalFunctionBase > > > d_atomicValenceDensityVector
Definition oncvClass.h:329
bool d_reproducible_output
Definition oncvClass.h:333
bool d_memoryOptMode
Definition oncvClass.h:280
bool d_useDevice
Definition oncvClass.h:279
std::vector< std::map< unsigned int, std::shared_ptr< AtomCenteredSphericalFunctionBase > > > d_atomicLocalPotVector
Definition oncvClass.h:326
std::shared_ptr< dftfe::basis::FEBasisOperations< ValueType, double, dftfe::utils::MemorySpace::HOST > > d_BasisOperatorHostPtr
Definition oncvClass.h:291
double getRmaxCoreDensity(unsigned int Znum)
std::map< std::pair< unsigned int, unsigned int >, std::shared_ptr< AtomCenteredSphericalFunctionBase > > d_atomicProjectorFnsMap
Definition oncvClass.h:271
void getRadialCoreDensity(unsigned int Znum, double rad, std::vector< double > &Val)
unsigned int d_nOMPThreads
Definition oncvClass.h:310
std::string d_dftfeScratchFolderName
Definition oncvClass.h:306
double getRadialLocalPseudo(unsigned int Znum, double rad)
Definition MemoryStorage.h:33
Definition FEBasisOperations.h:30
@ DEVICE
Definition MemorySpaceType.h:36
Definition pseudoPotentialToDftfeConverter.cc:34
T type
Definition dftfeDataTypes.h:115