DFT-FE 1.1.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 basisOperationsPtr,
49 std::shared_ptr<
51 double,
53 basisOperationsPtrHost,
54 std::shared_ptr<
56 pseudopotentialClassPtr,
57 std::shared_ptr<excManager<memorySpace>> excManagerPtr,
58 dftParameters *dftParamsPtr,
59 const dftfe::uInt densityQuadratureID,
60 const dftfe::uInt lpspQuadratureID,
61 const dftfe::uInt feOrderPlusOneQuadratureID,
62 const MPI_Comm &mpi_comm_parent,
63 const MPI_Comm &mpi_comm_domain);
64
65
66
67 /**
68 * @brief Computing Y = scalarOX*OX + scalarX*X + scalarY*Y for a given X and Y in full precision
69 *
70 * @param src X vector
71 * @param scalarHX scalar for OX
72 * @param scalarY scalar for Y
73 * @param scalarX scalar for X
74 * @param dst Y vector
75 * @param useApproximateMatrixEntries flag to use approximate overlap matrix
76 */
77 void
80 const double scalarOX,
81 const double scalarY,
82 const double scalarX,
84 const bool useApproximateMatrixEntries = true);
85
86
87
88 /**
89 * @brief Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in full precision
90 *
91 * @param src X vector
92 * @param scalarOinvX scalar for O^{-1}X
93 * @param scalarY scalar for Y
94 * @param scalarX scalar for X
95 * @param dst Y vector
96 */
97 void
100 const double scalarOinvX,
101 const double scalarY,
102 const double scalarX,
104
105 /**
106 * @brief Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in Reduced precision
107 *
108 * @param src X vector
109 * @param scalarOinvX scalar for O^{-1}X
110 * @param scalarY scalar for Y
111 * @param scalarX scalar for X
112 * @param dst Y vector
113 */
114 void
117 &src,
118 const double scalarOinvX,
119 const double scalarY,
120 const double scalarX,
122 &dst);
123
124 /**
125 * @brief Computing Y = scalarHX*HM^{-1}X + scalarX*X + scalarY*Y for a given X and Y in reduced precision
126 *
127 * @param src X vector
128 * @param scalarHX scalar for HX
129 * @param scalarY scalar for Y
130 * @param scalarX scalar for X
131 * @param dst Y vector
132 * @param onlyHPrimePartForFirstOrderDensityMatResponse flag to compute only HPrime part for first order density matrix response
133 * @param skip1 flag to skip extraction
134 * @param skip2 flag to skip nonLoal All Reduce
135 * @param skip3 flag to skip local HX and Assembly
136 */
137
138 void
140 memorySpace> &src,
141 const double scalarHX,
142 const double scalarY,
143 const double scalarX,
145 memorySpace> &dst,
146 const bool onlyHPrimePartForFirstOrderDensityMatResponse,
147 const bool skip1,
148 const bool skip2,
149 const bool skip3);
150
151 /**
152 * @brief Computing Y = scalarHX*M^{-1}HX + scalarX*X + scalarY*Y for a given X and Y in full precision
153 *
154 * @param src X vector
155 * @param scalarHX scalar for HX
156 * @param scalarY scalar for Y
157 * @param scalarX scalar for X
158 * @param dst Y vector
159 * @param onlyHPrimePartForFirstOrderDensityMatResponse flag to compute only HPrime part for first order density matrix response
160 * @param skip1 flag to skip extraction
161 * @param skip2 flag to skip nonLoal All Reduce
162 * @param skip3 flag to skip local HX and Assembly
163 */
164
165 void
168 const double scalarHX,
169 const double scalarY,
170 const double scalarX,
172 const bool onlyHPrimePartForFirstOrderDensityMatResponse = false,
173 const bool skip1 = false,
174 const bool skip2 = false,
175 const bool skip3 = false);
176
179
180
181 protected:
204 using KohnShamDFTBaseOperator<memorySpace>::
205 d_pseudopotentialNonLocalProjectorTimesVectorBlockSinglePrec;
209 using KohnShamDFTBaseOperator<memorySpace>::d_VeffJxW;
225 using KohnShamDFTBaseOperator<memorySpace>::d_spinIndex;
238 using KohnShamDFTBaseOperator<memorySpace>::pcout;
241 using KohnShamDFTBaseOperator<memorySpace>::d_useHubbard;
249 };
250} // namespace dftfe
251#endif
std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::numberFP32, memorySpace > > d_pseudopotentialNonLocalOperatorSinglePrec
Definition KohnShamDFTBaseOperator.h:356
dftfe::uInt d_cellsBlockSizeHamiltonianConstruction
Definition KohnShamDFTBaseOperator.h:440
std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > d_BLASWrapperPtr
Definition KohnShamDFTBaseOperator.h:359
std::shared_ptr< hubbard< dataTypes::number, memorySpace > > d_hubbardClassPtr
Definition KohnShamDFTBaseOperator.h:449
dftfe::utils::MemoryStorage< dataTypes::numberFP32, memorySpace > d_cellWaveFunctionMatrixSrcSinglePrec
Definition KohnShamDFTBaseOperator.h:388
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, dftfe::utils::MemorySpace::HOST > > d_basisOperationsPtrHost
Definition KohnShamDFTBaseOperator.h:367
dftfe::utils::MemoryStorage< double, memorySpace > d_VeffExtPotJxW
Definition KohnShamDFTBaseOperator.h:404
dftfe::uInt d_HamiltonianIndex
Definition KohnShamDFTBaseOperator.h:434
void reinitNumberWavefunctions(const dftfe::uInt numWfc)
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_srcNonLocalTemp
Definition KohnShamDFTBaseOperator.h:452
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_srcNonLocalTempSinglePrec
Definition KohnShamDFTBaseOperator.h:457
bool d_isExternalPotCorrHamiltonianComputed
Definition KohnShamDFTBaseOperator.h:435
dftfe::utils::MemoryStorage< double, memorySpace > tempHamMatrixRealBlock
Definition KohnShamDFTBaseOperator.h:426
std::shared_ptr< dftUtils::constraintMatrixInfo< memorySpace > > inverseMassVectorScaledConstraintsNoneDataInfoPtr
Definition KohnShamDFTBaseOperator.h:418
dftfe::utils::MemoryStorage< double, memorySpace > d_cellHamiltonianMatrixExtPot
Definition KohnShamDFTBaseOperator.h:379
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:398
dftfe::uInt d_nOMPThreads
Definition KohnShamDFTBaseOperator.h:443
dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > d_cellWaveFunctionMatrixSrc
Definition KohnShamDFTBaseOperator.h:383
dftfe::utils::MemoryStorage< double, memorySpace > tempHamMatrixImagBlock
Definition KohnShamDFTBaseOperator.h:427
const dftfe::uInt n_mpi_processes
Definition KohnShamDFTBaseOperator.h:438
bool d_useHubbard
Definition KohnShamDFTBaseOperator.h:450
const MPI_Comm d_mpiCommParent
Definition KohnShamDFTBaseOperator.h:436
const MPI_Comm d_mpiCommDomain
Definition KohnShamDFTBaseOperator.h:437
std::vector< dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > > d_cellHamiltonianMatrix
Definition KohnShamDFTBaseOperator.h:375
dftfe::uInt d_numVectorsInternal
Definition KohnShamDFTBaseOperator.h:442
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_pseudopotentialNonLocalProjectorTimesVectorBlock
Definition KohnShamDFTBaseOperator.h:393
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_dstNonLocalTempSinglePrec
Definition KohnShamDFTBaseOperator.h:459
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_dstNonLocalTemp
Definition KohnShamDFTBaseOperator.h:454
dftfe::utils::MemoryStorage< double, memorySpace > d_VeffJxW
Definition KohnShamDFTBaseOperator.h:403
const dftfe::uInt d_densityQuadratureID
Definition KohnShamDFTBaseOperator.h:429
dftfe::utils::MemoryStorage< double, memorySpace > d_invJacderExcWithSigmaTimesGradRhoJxW
Definition KohnShamDFTBaseOperator.h:407
dftfe::uInt d_kPointIndex
Definition KohnShamDFTBaseOperator.h:432
dftParameters * d_dftParamsPtr
Definition KohnShamDFTBaseOperator.h:372
dftfe::utils::MemoryStorage< dataTypes::number, memorySpace > d_cellWaveFunctionMatrixDst
Definition KohnShamDFTBaseOperator.h:385
std::shared_ptr< excManager< memorySpace > > d_excManagerPtr
Definition KohnShamDFTBaseOperator.h:371
const dftfe::uInt d_feOrderPlusOneQuadratureID
Definition KohnShamDFTBaseOperator.h:431
std::shared_ptr< dftfe::basis::FEBasisOperations< dataTypes::number, double, memorySpace > > d_basisOperationsPtr
Definition KohnShamDFTBaseOperator.h:362
const dftfe::utils::MemoryStorage< double, memorySpace > & getSqrtMassVector()
std::shared_ptr< dftUtils::constraintMatrixInfo< memorySpace > > inverseSqrtMassVectorScaledConstraintsNoneDataInfoPtr
Definition KohnShamDFTBaseOperator.h:420
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_tempBlockVectorOverlapInvXSinglePrec
Definition KohnShamDFTBaseOperator.h:400
dftfe::uInt d_spinIndex
Definition KohnShamDFTBaseOperator.h:433
const dftfe::uInt this_mpi_process
Definition KohnShamDFTBaseOperator.h:439
KohnShamDFTBaseOperator(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::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)
dftfe::utils::MemoryStorage< dataTypes::numberFP32, memorySpace > d_cellWaveFunctionMatrixDstSinglePrec
Definition KohnShamDFTBaseOperator.h:390
dftfe::uInt d_cellsBlockSizeHX
Definition KohnShamDFTBaseOperator.h:441
std::shared_ptr< AtomicCenteredNonLocalOperator< dataTypes::number, memorySpace > > d_pseudopotentialNonLocalOperator
Definition KohnShamDFTBaseOperator.h:341
std::vector< double > d_kPointWeights
Definition KohnShamDFTBaseOperator.h:424
std::vector< dftfe::utils::MemoryStorage< dataTypes::numberFP32, memorySpace > > d_cellHamiltonianMatrixSinglePrec
Definition KohnShamDFTBaseOperator.h:377
std::vector< double > d_kPointCoordinates
Definition KohnShamDFTBaseOperator.h:422
dftfe::utils::MemoryStorage< dftfe::uInt, memorySpace > d_mapNodeIdToProcId
Definition KohnShamDFTBaseOperator.h:461
dealii::TimerOutput computing_timer
Definition KohnShamDFTBaseOperator.h:447
const dftfe::utils::MemoryStorage< double, memorySpace > & getInverseSqrtMassVector()
std::shared_ptr< dftfe::pseudopotentialBaseClass< dataTypes::number, memorySpace > > d_pseudopotentialClassPtr
Definition KohnShamDFTBaseOperator.h:370
dealii::ConditionalOStream pcout
Definition KohnShamDFTBaseOperator.h:444
const dftfe::uInt d_lpspQuadratureID
Definition KohnShamDFTBaseOperator.h:430
KohnShamDFTStandardEigenOperator(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::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)
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 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.
Definition FEBasisOperations.h:85
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:42
float numberFP32
Definition dftfeDataTypes.h:43
@ HOST
Definition MemorySpaceType.h:34
Definition pseudoPotentialToDftfeConverter.cc:34
std::uint32_t uInt
Definition TypeConfig.h:10