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
31 template <dftfe::utils::MemorySpace memorySpace>
33 {
34 public:
35 /**
36 * @brief compute local descriptors of the aux basis electron-density
37 * representation at the supplied set of points using
38 */
39 virtual void
41 const std::vector<double> &Points,
42 std::unordered_map<DensityDescriptorDataAttributes, std::vector<double>>
43 &densityData) = 0;
44
45 /**
46 * @brief Compute aux basis overlap matrix batchwise contribution from
47 * supplied set of quadrature points and their associated weights
48 */
49 virtual void
50 evalOverlapMatrixStart(const std::vector<double> &quadpts,
51 const std::vector<double> &quadWt) = 0;
52
53 /**
54 * @brief for MPI accumulation
55 */
56 virtual void
57 evalOverlapMatrixEnd(const MPI_Comm &mpiComm) = 0;
58
59 // FIXME: to be extended to memoryspace
60 /**
61 * @brief Projects the KS density matrix to aux basis (L2 projection) batch wise
62 */
63 virtual void
65 const std::unordered_map<std::string, std::vector<dataTypes::number>>
66 &projectionInputsDataType,
67 const std::unordered_map<std::string, std::vector<double>>
68 & projectionInputsReal,
69 const int iSpin) = 0;
70
71 /**
72 * @brief for MPI accumulation
73 */
74 virtual void
75 projectDensityMatrixEnd(const MPI_Comm &mpiComm) = 0;
76
77
78 /**
79 * @brief Projects the quadrature density to aux basis (L2 projection) batch wise
80 */
81 virtual void
83 const std::unordered_map<std::string, std::vector<double>>
84 &projectionInputs) = 0;
85
86 /**
87 * @brief for MPI accumulation
88 */
89 virtual void
90 projectDensityEnd(const MPI_Comm &mpiComm) = 0;
91 };
92} // namespace dftfe
93
94#endif // DFTFE_AUXDM_AUXDENSITYMATRIX_H
Definition AuxDensityMatrix.h:33
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 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 projectDensityEnd(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 int iSpin)=0
Projects the KS density matrix to aux basis (L2 projection) batch wise.
virtual void applyLocalOperations(const std::vector< double > &Points, std::unordered_map< DensityDescriptorDataAttributes, std::vector< double > > &densityData)=0
compute local descriptors of the aux basis electron-density representation at the supplied set of poi...
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