-
Notifications
You must be signed in to change notification settings - Fork 34
/
KohnShamHamiltonianOperator.h
283 lines (244 loc) · 10.8 KB
/
KohnShamHamiltonianOperator.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
// ---------------------------------------------------------------------
//
// Copyright (c) 2017-2022 The Regents of the University of Michigan and DFT-FE
// authors.
//
// This file is part of the DFT-FE code.
//
// The DFT-FE code is free software; you can use it, redistribute
// it, and/or modify it under the terms of the GNU Lesser General
// Public License as published by the Free Software Foundation; either
// version 2.1 of the License, or (at your option) any later version.
// The full text of the license can be found in the file LICENSE at
// the top level of the DFT-FE distribution.
//
// ---------------------------------------------------------------------
//
#ifndef kohnShamHamiltonianOperatorClass_H_
#define kohnShamHamiltonianOperatorClass_H_
#include <constants.h>
#include <constraintMatrixInfo.h>
#include <headers.h>
#include <operator.h>
#include <BLASWrapper.h>
#include <FEBasisOperations.h>
#include <oncvClass.h>
#include <AuxDensityMatrix.h>
namespace dftfe
{
template <dftfe::utils::MemorySpace memorySpace>
class KohnShamHamiltonianOperator : public operatorDFTClass<memorySpace>
{
public:
KohnShamHamiltonianOperator(
std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
BLASWrapperPtr,
std::shared_ptr<
dftfe::basis::FEBasisOperations<dataTypes::number, double, memorySpace>>
basisOperationsPtr,
std::shared_ptr<
dftfe::basis::FEBasisOperations<dataTypes::number,
double,
dftfe::utils::MemorySpace::HOST>>
basisOperationsPtrHost,
std::shared_ptr<dftfe::oncvClass<dataTypes::number, memorySpace>>
oncvClassPtr,
std::shared_ptr<excManager> excManagerPtr,
dftParameters * dftParamsPtr,
const unsigned int densityQuadratureID,
const unsigned int lpspQuadratureID,
const unsigned int feOrderPlusOneQuadratureID,
const MPI_Comm & mpi_comm_parent,
const MPI_Comm & mpi_comm_domain);
void
init(const std::vector<double> &kPointCoordinates,
const std::vector<double> &kPointWeights);
void
resetExtPotHamFlag();
const MPI_Comm &
getMPICommunicatorDomain();
dftUtils::constraintMatrixInfo<dftfe::utils::MemorySpace::HOST> *
getOverloadedConstraintMatrixHost() const;
dftUtils::constraintMatrixInfo<memorySpace> *
getOverloadedConstraintMatrix() const
{
return &(d_basisOperationsPtr
->d_constraintInfo[d_basisOperationsPtr->d_dofHandlerID]);
}
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &
getScratchFEMultivector(const unsigned int numVectors,
const unsigned int index);
dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32, memorySpace> &
getScratchFEMultivectorSinglePrec(const unsigned int numVectors,
const unsigned int index);
/**
* @brief Computes effective potential involving exchange-correlation functionals
* @param auxDensityMatrixRepresentation core plus valence electron-density
* @param phiValues electrostatic potential arising both from electron-density and nuclear charge
*/
void
computeVEff(
std::shared_ptr<AuxDensityMatrix> auxDensityXCRepresentationPtr,
const dftfe::utils::MemoryStorage<double, dftfe::utils::MemorySpace::HOST>
& phiValues,
const unsigned int spinIndex = 0);
/**
* @brief Sets the V-eff potential
* @param vKS_quadValues the input V-KS values stored at the quadrature points
* @param spinIndex spin index
*/
void
setVEff(
const std::vector<
dftfe::utils::MemoryStorage<double, dftfe::utils::MemorySpace::HOST>>
& vKS_quadValues,
const unsigned int spinIndex);
void
computeVEffExternalPotCorr(
const std::map<dealii::CellId, std::vector<double>>
&externalPotCorrValues);
void
computeVEffPrime(
std::shared_ptr<AuxDensityMatrix> auxDensityXCRepresentationPtr,
const std::vector<
dftfe::utils::MemoryStorage<double, dftfe::utils::MemorySpace::HOST>>
&rhoPrimeValues,
const std::vector<
dftfe::utils::MemoryStorage<double, dftfe::utils::MemorySpace::HOST>>
&gradRhoPrimeValues,
const dftfe::utils::MemoryStorage<double, dftfe::utils::MemorySpace::HOST>
& phiPrimeValues,
const unsigned int spinIndex);
/**
* @brief sets the data member to appropriate kPoint and spin Index
*
* @param kPointIndex k-point Index to set
*/
void
reinitkPointSpinIndex(const unsigned int kPointIndex,
const unsigned int spinIndex);
void
reinitNumberWavefunctions(const unsigned int numWfc);
const dftfe::utils::MemoryStorage<double, memorySpace> &
getInverseSqrtMassVector();
const dftfe::utils::MemoryStorage<double, memorySpace> &
getSqrtMassVector();
void
computeCellHamiltonianMatrix(
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
void
computeCellHamiltonianMatrixExtPotContribution();
void
HX(dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &src,
const double scalarHX,
const double scalarY,
const double scalarX,
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &dst,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
void
HXCheby(
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &src,
const double scalarHX,
const double scalarY,
const double scalarX,
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &dst,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false,
const bool skip1 = false,
const bool skip2 = false,
const bool skip3 = false);
void
HXCheby(dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32,
memorySpace> &src,
const double scalarHX,
const double scalarY,
const double scalarX,
dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32,
memorySpace> &dst,
const bool onlyHPrimePartForFirstOrderDensityMatResponse,
const bool skip1,
const bool skip2,
const bool skip3);
void
HXRR(
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &src,
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &dstHX,
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> &dstMX,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false);
private:
void
setVEffExternalPotCorrToZero();
std::shared_ptr<
AtomicCenteredNonLocalOperator<dataTypes::number, memorySpace>>
d_ONCVnonLocalOperator;
std::shared_ptr<
AtomicCenteredNonLocalOperator<dataTypes::numberFP32, memorySpace>>
d_ONCVnonLocalOperatorSinglePrec;
std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace>>
d_BLASWrapperPtr;
std::shared_ptr<
dftfe::basis::FEBasisOperations<dataTypes::number, double, memorySpace>>
d_basisOperationsPtr;
std::shared_ptr<
dftfe::basis::FEBasisOperations<dataTypes::number,
double,
dftfe::utils::MemorySpace::HOST>>
d_basisOperationsPtrHost;
std::shared_ptr<dftfe::oncvClass<dataTypes::number, memorySpace>>
d_oncvClassPtr;
std::shared_ptr<excManager> d_excManagerPtr;
dftParameters * d_dftParamsPtr;
std::vector<dftfe::utils::MemoryStorage<dataTypes::number, memorySpace>>
d_cellHamiltonianMatrix;
std::vector<dftfe::utils::MemoryStorage<dataTypes::numberFP32, memorySpace>>
d_cellHamiltonianMatrixSinglePrec;
dftfe::utils::MemoryStorage<double, memorySpace>
d_cellHamiltonianMatrixExtPot;
dftfe::utils::MemoryStorage<dataTypes::number, memorySpace>
d_cellWaveFunctionMatrixSrc;
dftfe::utils::MemoryStorage<dataTypes::number, memorySpace>
d_cellWaveFunctionMatrixDst;
dftfe::utils::MemoryStorage<dataTypes::numberFP32, memorySpace>
d_cellWaveFunctionMatrixSrcSinglePrec;
dftfe::utils::MemoryStorage<dataTypes::numberFP32, memorySpace>
d_cellWaveFunctionMatrixDstSinglePrec;
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace>
d_ONCVNonLocalProjectorTimesVectorBlock;
dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32, memorySpace>
d_ONCVNonLocalProjectorTimesVectorBlockSinglePrec;
dftfe::utils::MemoryStorage<double, memorySpace> d_VeffJxW;
dftfe::utils::MemoryStorage<double, memorySpace> d_VeffExtPotJxW;
dftfe::utils::MemoryStorage<double, memorySpace>
d_invJacderExcWithSigmaTimesGradRhoJxW;
std::vector<dftfe::utils::MemoryStorage<double, memorySpace>>
d_invJacKPointTimesJxW;
// Constraints scaled with inverse sqrt diagonal Mass Matrix
std::shared_ptr<dftUtils::constraintMatrixInfo<memorySpace>>
inverseMassVectorScaledConstraintsNoneDataInfoPtr;
std::shared_ptr<dftUtils::constraintMatrixInfo<memorySpace>>
inverseSqrtMassVectorScaledConstraintsNoneDataInfoPtr;
// kPoint cartesian coordinates
std::vector<double> d_kPointCoordinates;
// k point weights
std::vector<double> d_kPointWeights;
dftfe::utils::MemoryStorage<double, memorySpace> tempHamMatrixRealBlock;
dftfe::utils::MemoryStorage<double, memorySpace> tempHamMatrixImagBlock;
const unsigned int d_densityQuadratureID;
const unsigned int d_lpspQuadratureID;
const unsigned int d_feOrderPlusOneQuadratureID;
unsigned int d_kPointIndex;
unsigned int d_spinIndex;
unsigned int d_HamiltonianIndex;
bool d_isExternalPotCorrHamiltonianComputed;
const MPI_Comm d_mpiCommParent;
const MPI_Comm d_mpiCommDomain;
const unsigned int n_mpi_processes;
const unsigned int this_mpi_process;
unsigned int d_cellsBlockSizeHamiltonianConstruction;
unsigned int d_cellsBlockSizeHX;
unsigned int d_numVectorsInternal;
dealii::ConditionalOStream pcout;
// compute-time logger
dealii::TimerOutput computing_timer;
};
} // namespace dftfe
#endif