DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
AuxDensityMatrix.h
Go to the documentation of this file.
1//
2// Created by Arghadwip Paul.
3//
4
5#ifndef DFTFE_AUXDM_AUXDENSITYMATRIX_H
6#define DFTFE_AUXDM_AUXDENSITYMATRIX_H
7
8#include <vector>
9#include <utility>
10#include <map>
11#include <string>
12#include <unordered_map>
13#include <mpi.h>
14#include <dftUtils.h>
15
16namespace dftfe
17{
30
37
38 template <dftfe::utils::MemorySpace memorySpace>
40 {
41 public:
42 /**
43 * @brief compute local descriptors of the aux basis electron-density
44 * representation at the supplied range of Quadrature index range
45 */
46 virtual void
48 const std::pair<dftfe::uInt, dftfe::uInt> &quadIndexRange,
49 std::unordered_map<DensityDescriptorDataAttributes, std::vector<double>>
50 &densityData) = 0;
51
52
53
54 virtual void
56 const std::pair<dftfe::uInt, dftfe::uInt> &quadIndexRange,
57 std::unordered_map<WfcDescriptorDataAttributes, std::vector<double>>
58 &wfcData) = 0;
59
60 /**
61 * @brief Compute aux basis overlap matrix batchwise contribution from
62 * supplied set of quadrature points and their associated weights
63 */
64 virtual void
65 evalOverlapMatrixStart(const std::vector<double> &quadpts,
66 const std::vector<double> &quadWt) = 0;
67
68 /**
69 * @brief for MPI accumulation
70 */
71 virtual void
72 evalOverlapMatrixEnd(const MPI_Comm &mpiComm) = 0;
73
74 // FIXME: to be extended to memoryspace
75 /**
76 * @brief Projects the KS density matrix to aux basis (L2 projection) batch wise
77 */
78 virtual void
80 const std::unordered_map<std::string, std::vector<dataTypes::number>>
81 &projectionInputsDataType,
82 const std::unordered_map<std::string, std::vector<double>>
83 &projectionInputsReal,
84 const dftfe::Int iSpin) = 0;
85
86 /**
87 * @brief for MPI accumulation
88 */
89 virtual void
90 projectDensityMatrixEnd(const MPI_Comm &mpiComm) = 0;
91
92
93 /**
94 * @brief Projects the quadrature density to aux basis (L2 projection) batch wise
95 */
96 virtual void
98 const std::unordered_map<std::string, std::vector<double>>
99 &projectionInputs) = 0;
100
101 /**
102 * @brief for MPI accumulation
103 */
104 virtual void
105 projectDensityEnd(const MPI_Comm &mpiComm) = 0;
106 };
107} // namespace dftfe
108
109#endif // DFTFE_AUXDM_AUXDENSITYMATRIX_H
Definition AuxDensityMatrix.h:40
virtual void evalOverlapMatrixEnd(const MPI_Comm &mpiComm)=0
for MPI accumulation
virtual void evalOverlapMatrixStart(const std::vector< double > &quadpts, const std::vector< double > &quadWt)=0
Compute aux basis overlap matrix batchwise contribution from supplied set of quadrature points and th...
virtual void applyLocalOperations(const std::pair< dftfe::uInt, dftfe::uInt > &quadIndexRange, std::unordered_map< DensityDescriptorDataAttributes, std::vector< double > > &densityData)=0
compute local descriptors of the aux basis electron-density representation at the supplied range of Q...
virtual void applyLocalOperations(const std::pair< dftfe::uInt, dftfe::uInt > &quadIndexRange, std::unordered_map< WfcDescriptorDataAttributes, std::vector< double > > &wfcData)=0
virtual void projectDensityStart(const std::unordered_map< std::string, std::vector< double > > &projectionInputs)=0
Projects the quadrature density to aux basis (L2 projection) batch wise.
virtual void projectDensityMatrixEnd(const MPI_Comm &mpiComm)=0
for MPI accumulation
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)=0
Projects the KS density matrix to aux basis (L2 projection) batch wise.
virtual void projectDensityEnd(const MPI_Comm &mpiComm)=0
for MPI accumulation
Definition pseudoPotentialToDftfeConverter.cc:34
DensityDescriptorDataAttributes
Definition AuxDensityMatrix.h:19
@ gradValuesSpinUp
Definition AuxDensityMatrix.h:23
@ laplacianSpinUp
Definition AuxDensityMatrix.h:27
@ gradValuesSpinDown
Definition AuxDensityMatrix.h:24
@ valuesSpinUp
Definition AuxDensityMatrix.h:21
@ laplacianSpinDown
Definition AuxDensityMatrix.h:28
@ hessianSpinUp
Definition AuxDensityMatrix.h:25
@ hessianSpinDown
Definition AuxDensityMatrix.h:26
@ valuesTotal
Definition AuxDensityMatrix.h:20
@ valuesSpinDown
Definition AuxDensityMatrix.h:22
std::int32_t Int
Definition TypeConfig.h:11
WfcDescriptorDataAttributes
Definition AuxDensityMatrix.h:32
@ tauSpinUp
Definition AuxDensityMatrix.h:34
@ tauSpinDown
Definition AuxDensityMatrix.h:35
@ tauTotal
Definition AuxDensityMatrix.h:33