DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
AuxDensityMatrixAtomicBasis.h
Go to the documentation of this file.
1//
2// Created by Arghadwip Paul, Sambit Das
3//
4
5#ifndef DFTFE_AUXDM_AUXDENSITYMATRIXATOMICBASIS_H
6#define DFTFE_AUXDM_AUXDENSITYMATRIXATOMICBASIS_H
7
8#include "AuxDensityMatrix.h"
9#include "AtomicBasis.h"
10#include "AtomicBasisData.h"
11#include <vector>
12#include <utility>
13#include <map>
14#include <algorithm>
15
16
17namespace dftfe
18{
19 template <dftfe::utils::MemorySpace memorySpace>
21 {
22 public:
23 void
25 const AtomicBasis::BasisType basisType,
26 const std::vector<std::pair<std::string, std::vector<double>>>
27 & atomCoords,
28 const std::unordered_map<std::string, std::string> &atomBasisFileNames,
29 const int nSpin,
30 const int maxDerOrder);
31
32 void
34 const std::vector<double> &quadpts,
35 std::unordered_map<DensityDescriptorDataAttributes, std::vector<double>>
36 &densityData) override;
37
38 void
39 evalOverlapMatrixStart(const std::vector<double> &quadpts,
40 const std::vector<double> &quadWt) override;
41
42 void
43 evalOverlapMatrixEnd(const MPI_Comm &mpiComm) override;
44
45 /**
46 *
47 * @param projectionInputs is a map from string to inputs needed
48 * for projection.
49 * eg - projectionInputsReal["quadpts"],
50 * projectionInputsReal["quadWt"],
51 * projectionInputsDataType["psiFunc"],
52 * projectionInputsReal["fValues"]
53 *
54 * psiFunc The SCF wave function or eigen function in FE Basis.
55 * psiFunc(quad_index, wfc_index),
56 * quad_index is fastest.
57 * fValues are the occupancies.
58 *
59 * @param iSpin indicates up (iSpin = 0) or down (iSpin = 0) spin.
60 *
61 */
62 virtual void
64 const std::unordered_map<std::string, std::vector<dataTypes::number>>
65 &projectionInputsDataType,
66 const std::unordered_map<std::string, std::vector<double>>
67 & projectionInputsReal,
68 const int iSpin) override;
69
70 void
71 projectDensityMatrixEnd(const MPI_Comm &mpiComm) override;
72
73 void
75 const std::unordered_map<std::string, std::vector<double>>
76 &projectionInputs) override;
77
78 void
79 projectDensityEnd(const MPI_Comm &mpiComm) override;
80
81
82 private:
85 std::unique_ptr<AtomicBasis> d_atomicBasisPtr;
87
90
91 std::vector<double> d_DM;
92 std::vector<double> d_SMatrix;
93 std::vector<double> d_SMatrixInv;
94 std::vector<double> d_basisWFCInnerProducts;
95 std::vector<double> d_fValues;
96
97 int d_nWFC;
99
100
101 void
103 std::vector<double> &
105 };
106} // namespace dftfe
107#endif // DFTFE_AUXDM_AUXDENSITYMATRIXATOMICBASIS_H
Definition AtomicBasisData.h:15
BasisType
Definition AtomicBasis.h:34
Definition AuxDensityMatrixAtomicBasis.h:21
int d_nSpin
Definition AuxDensityMatrixAtomicBasis.h:84
std::vector< double > & getOverlapMatrixInv()
int d_maxDerOrder
Definition AuxDensityMatrixAtomicBasis.h:89
int d_nBasis
Definition AuxDensityMatrixAtomicBasis.h:88
void reinit(const AtomicBasis::BasisType basisType, const std::vector< std::pair< std::string, std::vector< double > > > &atomCoords, const std::unordered_map< std::string, std::string > &atomBasisFileNames, const int nSpin, const int maxDerOrder)
int d_nWFC
Definition AuxDensityMatrixAtomicBasis.h:97
int d_iSpin
Definition AuxDensityMatrixAtomicBasis.h:98
void projectDensityMatrixEnd(const MPI_Comm &mpiComm) override
for MPI accumulation
AtomicBasisData d_atomicBasisData
Definition AuxDensityMatrixAtomicBasis.h:86
void projectDensityEnd(const MPI_Comm &mpiComm) override
for MPI accumulation
std::vector< double > d_SMatrix
Definition AuxDensityMatrixAtomicBasis.h:92
std::unique_ptr< AtomicBasis > d_atomicBasisPtr
Definition AuxDensityMatrixAtomicBasis.h:85
void projectDensityStart(const std::unordered_map< std::string, std::vector< double > > &projectionInputs) override
Projects the quadrature density to aux basis (L2 projection) batch wise.
void applyLocalOperations(const std::vector< double > &quadpts, std::unordered_map< DensityDescriptorDataAttributes, std::vector< double > > &densityData) override
compute local descriptors of the aux basis electron-density representation at the supplied set of poi...
std::vector< double > d_basisWFCInnerProducts
Definition AuxDensityMatrixAtomicBasis.h:94
std::vector< double > d_SMatrixInv
Definition AuxDensityMatrixAtomicBasis.h:93
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_fValues
Definition AuxDensityMatrixAtomicBasis.h:95
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) override
std::vector< double > d_DM
Definition AuxDensityMatrixAtomicBasis.h:91
void evalOverlapMatrixEnd(const MPI_Comm &mpiComm) override
for MPI accumulation
int d_nQuad
Definition AuxDensityMatrixAtomicBasis.h:83
Definition AuxDensityMatrix.h:33
Definition pseudoPotentialToDftfeConverter.cc:34
DensityDescriptorDataAttributes
Definition AuxDensityMatrix.h:19