DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
AuxDensityMatrixFE.h
Go to the documentation of this file.
1//
2// Created by Sambit Das.
3//
4
5#ifndef DFTFE_AUXDM_AUXDENSITYMATRIXFE_H
6#define DFTFE_AUXDM_AUXDENSITYMATRIXFE_H
7
8#include <vector>
9#include <utility>
10#include <AuxDensityMatrix.h>
11
12namespace dftfe
13{
14 template <dftfe::utils::MemorySpace memorySpace>
15 class AuxDensityMatrixFE : public AuxDensityMatrix<memorySpace>
16 {
17 public:
18 // FIXME: to be implemented
19
20 void
23 &eigenVectorsFlattenedMemSpace,
24 const std::vector<std::vector<double>> &fractionalOccupancies);
25
26
27 void
29 const std::pair<dftfe::uInt, dftfe::uInt> &quadIndexRange,
30 std::unordered_map<DensityDescriptorDataAttributes, std::vector<double>>
31 &densityData) override;
32
33
34
35 void
37 const std::pair<dftfe::uInt, dftfe::uInt> &quadIndexRange,
38 std::unordered_map<WfcDescriptorDataAttributes, std::vector<double>>
39 &wfcData) override;
40
41 void
42 evalOverlapMatrixStart(const std::vector<double> &quadpts,
43 const std::vector<double> &quadWt) override;
44
45 void
46 evalOverlapMatrixEnd(const MPI_Comm &mpiComm) override;
47
48 virtual void
50 const std::unordered_map<std::string, std::vector<dataTypes::number>>
51 &projectionInputsDataType,
52 const std::unordered_map<std::string, std::vector<double>>
53 &projectionInputsReal,
54 const dftfe::Int iSpin) override;
55
56 void
57 projectDensityMatrixEnd(const MPI_Comm &mpiComm) override;
58
59 /**
60 * @brief Projects the quadrature density to aux basis (L2 projection).
61 * This is actually a copy call. All the local partition quadrature points
62 * must to be passed to this function in one go
63 *
64 * @param projectionInputs is a map from string to inputs needed
65 * for projection.
66 * projectionInputs["quadpts"],
67 * projectionInputs["quadWt"],
68 * projectionInputs["densityFunc"]
69 * projectionInputs["gradDensityFunc"]
70 *
71 * densityFunc The density Values at quad points
72 * densityFunc(spin_index, quad_index),
73 * quad_index is fastest.
74 *
75 * gradDensityFunc The density Values at quad points
76 * gradDensityFunc(spin_index, quad_index,dim_index),
77 * dim_index is fastest.
78 *
79 */
80 void
82 const std::unordered_map<std::string, std::vector<double>>
83 &projectionInputs) override;
84
85 void
86 projectDensityEnd(const MPI_Comm &mpiComm) override;
87
88 const std::vector<std::vector<double>> *
90
93
94 private:
97
98 const std::vector<std::vector<double>> *d_fractionalOccupancies;
99
100 std::vector<double> d_densityValsTotalAllQuads;
101 std::vector<double> d_densityValsSpinUpAllQuads;
102 std::vector<double> d_densityValsSpinDownAllQuads;
105 std::vector<double> d_tauValsTotalAllQuads;
106 std::vector<double> d_tauValsSpinUpAllQuads;
107 std::vector<double> d_tauValsSpinDownAllQuads;
108
109 std::vector<double> d_quadPointsAll;
110 std::vector<double> d_quadWeightsAll;
111 };
112} // namespace dftfe
113
114#endif // DFTFE_AUXDM_AUXDENSITYMATRIXFE_H
Definition AuxDensityMatrixFE.h:16
void projectDensityEnd(const MPI_Comm &mpiComm) override
for MPI accumulation
void evalOverlapMatrixStart(const std::vector< double > &quadpts, const std::vector< double > &quadWt) override
Compute aux basis overlap matrix batchwise contribution from supplied set of quadrature points and th...
std::vector< double > d_densityValsSpinDownAllQuads
Definition AuxDensityMatrixFE.h:102
std::vector< double > d_quadPointsAll
Definition AuxDensityMatrixFE.h:109
std::vector< double > d_tauValsSpinUpAllQuads
Definition AuxDensityMatrixFE.h:106
std::vector< double > d_densityValsSpinUpAllQuads
Definition AuxDensityMatrixFE.h:101
const std::vector< std::vector< double > > * getDensityMatrixComponents_occupancies() const
void applyLocalOperations(const std::pair< dftfe::uInt, dftfe::uInt > &quadIndexRange, std::unordered_map< DensityDescriptorDataAttributes, std::vector< double > > &densityData) override
compute local descriptors of the aux basis electron-density representation at the supplied range of Q...
std::vector< double > d_tauValsSpinDownAllQuads
Definition AuxDensityMatrixFE.h:107
virtual void projectDensityMatrixStart(const std::unordered_map< std::string, std::vector< dataTypes::number > > &projectionInputsDataType, const std::unordered_map< std::string, std::vector< double > > &projectionInputsReal, const dftfe::Int iSpin) override
Projects the KS density matrix to aux basis (L2 projection) batch wise.
void applyLocalOperations(const std::pair< dftfe::uInt, dftfe::uInt > &quadIndexRange, std::unordered_map< WfcDescriptorDataAttributes, std::vector< double > > &wfcData) override
std::vector< double > d_gradDensityValsSpinUpAllQuads
Definition AuxDensityMatrixFE.h:103
const dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > * getDensityMatrixComponents_wavefunctions() const
std::vector< double > d_densityValsTotalAllQuads
Definition AuxDensityMatrixFE.h:100
void evalOverlapMatrixEnd(const MPI_Comm &mpiComm) override
for MPI accumulation
const std::vector< std::vector< double > > * d_fractionalOccupancies
Definition AuxDensityMatrixFE.h:98
void setDensityMatrixComponents(const dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > &eigenVectorsFlattenedMemSpace, const std::vector< std::vector< double > > &fractionalOccupancies)
void projectDensityStart(const std::unordered_map< std::string, std::vector< double > > &projectionInputs) override
Projects the quadrature density to aux basis (L2 projection). This is actually a copy call....
std::vector< double > d_quadWeightsAll
Definition AuxDensityMatrixFE.h:110
const dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > * d_eigenVectorsFlattenedMemSpacePtr
Definition AuxDensityMatrixFE.h:96
std::vector< double > d_gradDensityValsSpinDownAllQuads
Definition AuxDensityMatrixFE.h:104
std::vector< double > d_tauValsTotalAllQuads
Definition AuxDensityMatrixFE.h:105
void projectDensityMatrixEnd(const MPI_Comm &mpiComm) override
for MPI accumulation
Definition AuxDensityMatrix.h:40
Definition MemoryStorage.h:33
Definition pseudoPotentialToDftfeConverter.cc:34
DensityDescriptorDataAttributes
Definition AuxDensityMatrix.h:19
std::int32_t Int
Definition TypeConfig.h:11
WfcDescriptorDataAttributes
Definition AuxDensityMatrix.h:32