DFT-FE 1.3.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
KohnShamDFTStandardEigenOperator.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (c) 2017-2025 The Regents of the University of Michigan and DFT-FE
4// authors.
5//
6// This file is part of the DFT-FE code.
7//
8// The DFT-FE code is free software; you can use it, redistribute
9// it, and/or modify it under the terms of the GNU Lesser General
10// Public License as published by the Free Software Foundation; either
11// version 2.1 of the License, or (at your option) any later version.
12// The full text of the license can be found in the file LICENSE at
13// the top level of the DFT-FE distribution.
14//
15// ---------------------------------------------------------------------
16//
17
18
19#ifndef KohnShamDFTStandardEigenOperatorClass_H_
20#define KohnShamDFTStandardEigenOperatorClass_H_
21#include <constants.h>
23#include <headers.h>
25#include <BLASWrapper.h>
26#include <FEBasisOperations.h>
27#include <oncvClass.h>
28#include <AuxDensityMatrix.h>
29
30#include "hubbardClass.h"
31
32namespace dftfe
33{
34 // KohnShamDFTStandardEigenOperator is derived from KohnShamDFTBaseOperator
35 // base class. This class is used where the underlying PDE is a
36 // standard eigen value problem. Currently used for the all-electron
37 // and norm conserving pseudopotential
38 template <dftfe::utils::MemorySpace memorySpace>
40 : public KohnShamDFTBaseOperator<memorySpace>
41 {
42 public:
45 BLASWrapperPtr,
46 std::shared_ptr<
48 BLASWrapperPtrHost,
49 std::shared_ptr<
51 basisOperationsPtr,
52 std::shared_ptr<
54 double,
56 basisOperationsPtrHost,
57 std::shared_ptr<
59 pseudopotentialClassPtr,
60 std::shared_ptr<excManager<memorySpace>> excManagerPtr,
61 dftParameters *dftParamsPtr,
62 const dftfe::uInt densityQuadratureID,
63 const dftfe::uInt lpspQuadratureID,
64 const dftfe::uInt feOrderPlusOneQuadratureID,
65 const MPI_Comm &mpi_comm_parent,
66 const MPI_Comm &mpi_comm_domain);
67
68
69
70 /**
71 * @brief Computing Y = scalarOX*OX + scalarX*X + scalarY*Y for a given X and Y in full precision
72 *
73 * @param src X vector
74 * @param scalarHX scalar for OX
75 * @param scalarY scalar for Y
76 * @param scalarX scalar for X
77 * @param dst Y vector
78 * @param useApproximateMatrixEntries flag to use approximate overlap matrix
79 */
80 void
83 const double scalarOX,
84 const double scalarY,
85 const double scalarX,
87 const bool useApproximateMatrixEntries = true);
88
89
90
91 /**
92 * @brief Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in full precision
93 *
94 * @param src X vector
95 * @param scalarOinvX scalar for O^{-1}X
96 * @param scalarY scalar for Y
97 * @param scalarX scalar for X
98 * @param dst Y vector
99 */
100 void
103 const double scalarOinvX,
104 const double scalarY,
105 const double scalarX,
107
108 /**
109 * @brief Computing Y = scalarOinvX*O^{-1/2}X + scalarX*X + scalarY*Y for a given X and Y in full precision
110 *
111 * @param src X vector
112 * @param scalarOinvX scalar for O^{-1}X
113 * @param scalarY scalar for Y
114 * @param scalarX scalar for X
115 * @param dst Y vector
116 */
117 void
120 const double scalarOinvX,
121 const double scalarY,
122 const double scalarX,
124
125 /**
126 * @brief Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in Reduced precision
127 *
128 * @param src X vector
129 * @param scalarOinvX scalar for O^{-1}X
130 * @param scalarY scalar for Y
131 * @param scalarX scalar for X
132 * @param dst Y vector
133 */
134 void
137 &src,
138 const double scalarOinvX,
139 const double scalarY,
140 const double scalarX,
142 &dst);
143
144 /**
145 * @brief Computing Y = scalarHX*HM^{-1}X + scalarX*X + scalarY*Y for a given X and Y in reduced precision
146 *
147 * @param src X vector
148 * @param scalarHX scalar for HX
149 * @param scalarY scalar for Y
150 * @param scalarX scalar for X
151 * @param dst Y vector
152 * @param onlyHPrimePartForFirstOrderDensityMatResponse flag to compute only HPrime part for first order density matrix response
153 * @param skip1 flag to skip extraction
154 * @param skip2 flag to skip nonLoal All Reduce
155 * @param skip3 flag to skip local HX and Assembly
156 */
157
158 void
160 memorySpace> &src,
161 const double scalarHX,
162 const double scalarY,
163 const double scalarX,
165 memorySpace> &dst,
166 const bool onlyHPrimePartForFirstOrderDensityMatResponse,
167 const bool skip1,
168 const bool skip2,
169 const bool skip3);
170
171 /**
172 * @brief Computing Y = scalarHX*M^{-1}HX + scalarX*X + scalarY*Y for a given X and Y in full precision
173 *
174 * @param src X vector
175 * @param scalarHX scalar for HX
176 * @param scalarY scalar for Y
177 * @param scalarX scalar for X
178 * @param dst Y vector
179 * @param onlyHPrimePartForFirstOrderDensityMatResponse flag to compute only HPrime part for first order density matrix response
180 * @param skip1 flag to skip extraction
181 * @param skip2 flag to skip nonLoal All Reduce
182 * @param skip3 flag to skip local HX and Assembly
183 */
184
185 void
188 const double scalarHX,
189 const double scalarY,
190 const double scalarX,
192 const bool onlyHPrimePartForFirstOrderDensityMatResponse = false,
193 const bool skip1 = false,
194 const bool skip2 = false,
195 const bool skip3 = false);
196
199
200
201 protected:
224 using KohnShamDFTBaseOperator<memorySpace>::
225 d_pseudopotentialNonLocalProjectorTimesVectorBlockSinglePrec;
229 using KohnShamDFTBaseOperator<memorySpace>::d_VeffJxW;
245 using KohnShamDFTBaseOperator<memorySpace>::d_spinIndex;
258 using KohnShamDFTBaseOperator<memorySpace>::pcout;
261 using KohnShamDFTBaseOperator<memorySpace>::d_useHubbard;
269 };
270} // namespace dftfe
271#endif
std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::numberFP32, memorySpace > > d_pseudopotentialNonLocalOperatorSinglePrec
Definition KohnShamDFTBaseOperator.h:237
dftfe::uInt d_cellsBlockSizeHamiltonianConstruction
Definition KohnShamDFTBaseOperator.h:335
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
Definition KohnShamDFTBaseOperator.h:240
std::shared_ptr< hubbard< dataTypes::number, memorySpace > > d_hubbardClassPtr
Definition KohnShamDFTBaseOperator.h:344
dftfe::utils::MemoryStorage< dataTypes::numberFP32, memorySpace > d_cellWaveFunctionMatrixSrcSinglePrec
Definition KohnShamDFTBaseOperator.h:272
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > d_basisOperationsPtrHost
Definition KohnShamDFTBaseOperator.h:251
dftfe::utils::MemoryStorage< double, memorySpace > d_VeffExtPotJxW
Definition KohnShamDFTBaseOperator.h:289
dftfe::uInt d_HamiltonianIndex
Definition KohnShamDFTBaseOperator.h:329
void reinitNumberWavefunctions(const dftfe::uInt numWfc)
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_srcNonLocalTemp
Definition KohnShamDFTBaseOperator.h:347
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_srcNonLocalTempSinglePrec
Definition KohnShamDFTBaseOperator.h:352
bool d_isExternalPotCorrHamiltonianComputed
Definition KohnShamDFTBaseOperator.h:330
dftfe::utils::MemoryStorage< double, memorySpace > tempHamMatrixRealBlock
Definition KohnShamDFTBaseOperator.h:314
std::shared_ptr< dftUtils::constraintMatrixInfo< memorySpace > > inverseMassVectorScaledConstraintsNoneDataInfoPtr
Definition KohnShamDFTBaseOperator.h:306
dftfe::utils::MemoryStorage< double, memorySpace > d_cellHamiltonianMatrixExtPot
Definition KohnShamDFTBaseOperator.h:263
void reinitkPointSpinIndex(const dftfe::uInt kPointIndex, const dftfe::uInt spinIndex)
sets the data member to appropriate kPoint and spin Index
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_tempBlockVectorOverlapInvX
Definition KohnShamDFTBaseOperator.h:282
dftfe::uInt d_nOMPThreads
Definition KohnShamDFTBaseOperator.h:338
dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > d_cellWaveFunctionMatrixSrc
Definition KohnShamDFTBaseOperator.h:267
dftfe::utils::MemoryStorage< double, memorySpace > tempHamMatrixImagBlock
Definition KohnShamDFTBaseOperator.h:315
const dftfe::uInt n_mpi_processes
Definition KohnShamDFTBaseOperator.h:333
bool d_useHubbard
Definition KohnShamDFTBaseOperator.h:345
const MPI_Comm d_mpiCommParent
Definition KohnShamDFTBaseOperator.h:331
const MPI_Comm d_mpiCommDomain
Definition KohnShamDFTBaseOperator.h:332
std::vector< dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > > d_cellHamiltonianMatrix
Definition KohnShamDFTBaseOperator.h:259
dftfe::uInt d_numVectorsInternal
Definition KohnShamDFTBaseOperator.h:337
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_pseudopotentialNonLocalProjectorTimesVectorBlock
Definition KohnShamDFTBaseOperator.h:277
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_dstNonLocalTempSinglePrec
Definition KohnShamDFTBaseOperator.h:354
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_dstNonLocalTemp
Definition KohnShamDFTBaseOperator.h:349
dftfe::utils::MemoryStorage< double, memorySpace > d_VeffJxW
Definition KohnShamDFTBaseOperator.h:287
const dftfe::uInt d_densityQuadratureID
Definition KohnShamDFTBaseOperator.h:324
dftfe::utils::MemoryStorage< double, memorySpace > d_invJacderExcWithSigmaTimesGradRhoJxW
Definition KohnShamDFTBaseOperator.h:292
dftfe::uInt d_kPointIndex
Definition KohnShamDFTBaseOperator.h:327
dftParameters * d_dftParamsPtr
Definition KohnShamDFTBaseOperator.h:256
dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > d_cellWaveFunctionMatrixDst
Definition KohnShamDFTBaseOperator.h:269
KohnShamDFTBaseOperator(std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperPtrHost, 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::pseudopotentialBaseClass< dataTypes::number, memorySpace > > pseudopotentialClassPtr, std::shared_ptr< excManager< memorySpace > > excManagerPtr, dftParameters *dftParamsPtr, const dftfe::uInt densityQuadratureID, const dftfe::uInt lpspQuadratureID, const dftfe::uInt feOrderPlusOneQuadratureID, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain)
std::shared_ptr< excManager< memorySpace > > d_excManagerPtr
Definition KohnShamDFTBaseOperator.h:255
const dftfe::uInt d_feOrderPlusOneQuadratureID
Definition KohnShamDFTBaseOperator.h:326
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, memorySpace > > d_basisOperationsPtr
Definition KohnShamDFTBaseOperator.h:246
const dftfe::utils::MemoryStorage< double, memorySpace > & getSqrtMassVector()
std::shared_ptr< dftUtils::constraintMatrixInfo< memorySpace > > inverseSqrtMassVectorScaledConstraintsNoneDataInfoPtr
Definition KohnShamDFTBaseOperator.h:308
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_tempBlockVectorOverlapInvXSinglePrec
Definition KohnShamDFTBaseOperator.h:284
dftfe::uInt d_spinIndex
Definition KohnShamDFTBaseOperator.h:328
const dftfe::uInt this_mpi_process
Definition KohnShamDFTBaseOperator.h:334
dftfe::utils::MemoryStorage< dataTypes::numberFP32, memorySpace > d_cellWaveFunctionMatrixDstSinglePrec
Definition KohnShamDFTBaseOperator.h:274
dftfe::uInt d_cellsBlockSizeHX
Definition KohnShamDFTBaseOperator.h:336
std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > d_pseudopotentialNonLocalOperator
Definition KohnShamDFTBaseOperator.h:233
std::vector< double > d_kPointWeights
Definition KohnShamDFTBaseOperator.h:312
std::vector< dftfe::utils::MemoryStorage< dataTypes::numberFP32, memorySpace > > d_cellHamiltonianMatrixSinglePrec
Definition KohnShamDFTBaseOperator.h:261
std::vector< double > d_kPointCoordinates
Definition KohnShamDFTBaseOperator.h:310
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > d_mapNodeIdToProcId
Definition KohnShamDFTBaseOperator.h:356
dealii::TimerOutput computing_timer
Definition KohnShamDFTBaseOperator.h:342
const dftfe::utils::MemoryStorage< double, memorySpace > & getInverseSqrtMassVector()
std::shared_ptr< dftfe::pseudopotentialBaseClass< dataTypes::number, memorySpace > > d_pseudopotentialClassPtr
Definition KohnShamDFTBaseOperator.h:254
dealii::ConditionalOStream pcout
Definition KohnShamDFTBaseOperator.h:339
const dftfe::uInt d_lpspQuadratureID
Definition KohnShamDFTBaseOperator.h:325
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)
Computing Y = scalarHX*HM^{-1}X + scalarX*X + scalarY*Y for a given X and Y in reduced precision.
void overlapInverseMatrixTimesX(dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > &src, const double scalarOinvX, const double scalarY, const double scalarX, dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > &dst)
Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in Reduced precision.
void overlapSqrtInverseMatrixTimesX(dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &src, const double scalarOinvX, const double scalarY, const double scalarX, dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &dst)
Computing Y = scalarOinvX*O^{-1/2}X + scalarX*X + scalarY*Y for a given X and Y in full precision.
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)
Computing Y = scalarHX*M^{-1}HX + scalarX*X + scalarY*Y for a given X and Y in full precision.
void overlapInverseMatrixTimesX(dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &src, const double scalarOinvX, const double scalarY, const double scalarX, dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &dst)
Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in full precision.
void overlapMatrixTimesX(dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &src, const double scalarOX, const double scalarY, const double scalarX, dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > &dst, const bool useApproximateMatrixEntries=true)
Computing Y = scalarOX*OX + scalarX*X + scalarY*Y for a given X and Y in full precision.
KohnShamDFTStandardEigenOperator(std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > BLASWrapperPtr, std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > BLASWrapperPtrHost, 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::pseudopotentialBaseClass< dataTypes::number, memorySpace > > pseudopotentialClassPtr, std::shared_ptr< excManager< memorySpace > > excManagerPtr, dftParameters *dftParamsPtr, const dftfe::uInt densityQuadratureID, const dftfe::uInt lpspQuadratureID, const dftfe::uInt feOrderPlusOneQuadratureID, const MPI_Comm &mpi_comm_parent, const MPI_Comm &mpi_comm_domain)
Definition FEBasisOperations.h:87
Namespace which declares the input parameters and the functions to parse them from the input paramete...
Definition dftParameters.h:36
Definition excManager.h:28
Definition BLASWrapper.h:35
An class template to encapsulate a MultiVector. A MultiVector is a collection of vectors belonging t...
Definition MultiVector.h:127
Definition pseudopotentialBaseClass.h:60
double number
Definition dftfeDataTypes.h:41
float numberFP32
Definition dftfeDataTypes.h:42
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10