DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftfe::KohnShamHamiltonianOperator< memorySpace > Class Template Reference

#include <KohnShamHamiltonianOperator.h>

Inheritance diagram for dftfe::KohnShamHamiltonianOperator< memorySpace >:
dftfe::operatorDFTClass< memorySpace >

Public Member Functions

 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< memorySpace > > 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 ()
 
void resetKohnShamOp ()
 
const MPI_Comm & getMPICommunicatorDomain ()
 
dftUtils::constraintMatrixInfo< dftfe::utils::MemorySpace::HOST > * getOverloadedConstraintMatrixHost () const
 
dftUtils::constraintMatrixInfo< memorySpace > * getOverloadedConstraintMatrix () const
 
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)
 
void computeVEff (std::shared_ptr< AuxDensityMatrix< memorySpace > > auxDensityXCRepresentationPtr, const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > &phiValues, const unsigned int spinIndex=0)
 Computes effective potential involving exchange-correlation functionals.
 
void setVEff (const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > &vKS_quadValues, const unsigned int spinIndex)
 Sets the V-eff potential.
 
void computeVEffExternalPotCorr (const std::map< dealii::CellId, std::vector< double > > &externalPotCorrValues)
 
void computeVEffPrime (std::shared_ptr< AuxDensityMatrix< memorySpace > > 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)
 
void reinitkPointSpinIndex (const unsigned int kPointIndex, const unsigned int spinIndex)
 sets the data member to appropriate kPoint and spin Index
 
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)
 Computing Y = scalarHX*HX + scalarX*X + scalarY*Y for a given X and Y in full precision.
 
void HXWithLowdinOrthonormalisedInput (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)
 Computing Y = scalarHX*M^{1/2}HM^{1/2}X + scalarX*X + scalarY*Y for a given X and Y in full precision. Used for TD-DFT and Inverse DFT calc.
 
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.
 
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 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::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 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 setVEffExternalPotCorrToZero ()
 
- Public Member Functions inherited from dftfe::operatorDFTClass< memorySpace >
virtual ~operatorDFTClass ()
 Destructor.
 

Private Attributes

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< memorySpace > > d_excManagerPtr
 
dftParametersd_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::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_tempBlockVectorOverlapInvX
 
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_tempBlockVectorOverlapInvXSinglePrec
 
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
 
std::shared_ptr< dftUtils::constraintMatrixInfo< memorySpace > > inverseMassVectorScaledConstraintsNoneDataInfoPtr
 
std::shared_ptr< dftUtils::constraintMatrixInfo< memorySpace > > inverseSqrtMassVectorScaledConstraintsNoneDataInfoPtr
 
std::vector< double > d_kPointCoordinates
 
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
 
unsigned int d_nOMPThreads
 
dealii::ConditionalOStream pcout
 
dealii::TimerOutput computing_timer
 
std::shared_ptr< hubbard< dataTypes::number, memorySpace > > d_hubbardClassPtr
 
bool d_useHubbard
 
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_srcNonLocalTemp
 
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > d_dstNonLocalTemp
 
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_srcNonLocalTempSinglePrec
 
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > d_dstNonLocalTempSinglePrec
 
dftfe::utils::MemoryStorage< dftfe::global_size_type, memorySpace > d_mapNodeIdToProcId
 

Constructor & Destructor Documentation

◆ KohnShamHamiltonianOperator()

template<dftfe::utils::MemorySpace memorySpace>
dftfe::KohnShamHamiltonianOperator< memorySpace >::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< memorySpace > > 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 )

Member Function Documentation

◆ computeCellHamiltonianMatrix()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::computeCellHamiltonianMatrix ( const bool onlyHPrimePartForFirstOrderDensityMatResponse = false)

◆ computeCellHamiltonianMatrixExtPotContribution()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::computeCellHamiltonianMatrixExtPotContribution ( )

◆ computeVEff()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::computeVEff ( std::shared_ptr< AuxDensityMatrix< memorySpace > > auxDensityXCRepresentationPtr,
const dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > & phiValues,
const unsigned int spinIndex = 0 )

Computes effective potential involving exchange-correlation functionals.

Parameters
auxDensityMatrixRepresentationcore plus valence electron-density
phiValueselectrostatic potential arising both from electron-density and nuclear charge

◆ computeVEffExternalPotCorr()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::computeVEffExternalPotCorr ( const std::map< dealii::CellId, std::vector< double > > & externalPotCorrValues)

◆ computeVEffPrime()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::computeVEffPrime ( std::shared_ptr< AuxDensityMatrix< memorySpace > > 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 )

◆ getInverseSqrtMassVector()

template<dftfe::utils::MemorySpace memorySpace>
const dftfe::utils::MemoryStorage< double, memorySpace > & dftfe::KohnShamHamiltonianOperator< memorySpace >::getInverseSqrtMassVector ( )
virtual

◆ getMPICommunicatorDomain()

template<dftfe::utils::MemorySpace memorySpace>
const MPI_Comm & dftfe::KohnShamHamiltonianOperator< memorySpace >::getMPICommunicatorDomain ( )
virtual

◆ getOverloadedConstraintMatrix()

template<dftfe::utils::MemorySpace memorySpace>
dftUtils::constraintMatrixInfo< memorySpace > * dftfe::KohnShamHamiltonianOperator< memorySpace >::getOverloadedConstraintMatrix ( ) const
inlinevirtual

◆ getOverloadedConstraintMatrixHost()

template<dftfe::utils::MemorySpace memorySpace>
dftUtils::constraintMatrixInfo< dftfe::utils::MemorySpace::HOST > * dftfe::KohnShamHamiltonianOperator< memorySpace >::getOverloadedConstraintMatrixHost ( ) const
virtual

◆ getScratchFEMultivector()

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > & dftfe::KohnShamHamiltonianOperator< memorySpace >::getScratchFEMultivector ( const unsigned int numVectors,
const unsigned int index )
virtual

◆ getScratchFEMultivectorSinglePrec()

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > & dftfe::KohnShamHamiltonianOperator< memorySpace >::getScratchFEMultivectorSinglePrec ( const unsigned int numVectors,
const unsigned int index )
virtual

◆ getSqrtMassVector()

template<dftfe::utils::MemorySpace memorySpace>
const dftfe::utils::MemoryStorage< double, memorySpace > & dftfe::KohnShamHamiltonianOperator< memorySpace >::getSqrtMassVector ( )
virtual

◆ HX()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::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 )
virtual

Computing Y = scalarHX*HX + scalarX*X + scalarY*Y for a given X and Y in full precision.

Parameters
srcX vector
scalarHXscalar for HX
scalarYscalar for Y
scalarXscalar for X
dstY vector
onlyHPrimePartForFirstOrderDensityMatResponseflag to compute only HPrime part for first order density matrix response

Implements dftfe::operatorDFTClass< memorySpace >.

◆ HXCheby() [1/2]

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::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 )
virtual

Computing Y = scalarHX*M^{-1}HX + scalarX*X + scalarY*Y for a given X and Y in full precision.

Parameters
srcX vector
scalarHXscalar for HX
scalarYscalar for Y
scalarXscalar for X
dstY vector
onlyHPrimePartForFirstOrderDensityMatResponseflag to compute only HPrime part for first order density matrix response
skip1flag to skip extraction
skip2flag to skip nonLoal All Reduce
skip3flag to skip local HX and Assembly

Implements dftfe::operatorDFTClass< memorySpace >.

◆ HXCheby() [2/2]

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::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 )
virtual

Computing Y = scalarHX*HM^{-1}X + scalarX*X + scalarY*Y for a given X and Y in reduced precision.

Parameters
srcX vector
scalarHXscalar for HX
scalarYscalar for Y
scalarXscalar for X
dstY vector
onlyHPrimePartForFirstOrderDensityMatResponseflag to compute only HPrime part for first order density matrix response
skip1flag to skip extraction
skip2flag to skip nonLoal All Reduce
skip3flag to skip local HX and Assembly

Implements dftfe::operatorDFTClass< memorySpace >.

◆ HXWithLowdinOrthonormalisedInput()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::HXWithLowdinOrthonormalisedInput ( 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 )
virtual

Computing Y = scalarHX*M^{1/2}HM^{1/2}X + scalarX*X + scalarY*Y for a given X and Y in full precision. Used for TD-DFT and Inverse DFT calc.

Parameters
srcX vector
scalarHXscalar for HX
scalarYscalar for Y
scalarXscalar for X
dstY vector
onlyHPrimePartForFirstOrderDensityMatResponseflag to compute only HPrime part for first order density matrix response

Implements dftfe::operatorDFTClass< memorySpace >.

◆ init()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::init ( const std::vector< double > & kPointCoordinates,
const std::vector< double > & kPointWeights )
virtual

◆ overlapInverseMatrixTimesX() [1/2]

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::overlapInverseMatrixTimesX ( dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > & src,
const double scalarOinvX,
const double scalarY,
const double scalarX,
dftfe::linearAlgebra::MultiVector< dataTypes::number, memorySpace > & dst )
virtual

Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in full precision.

Parameters
srcX vector
scalarOinvXscalar for O^{-1}X
scalarYscalar for Y
scalarXscalar for X
dstY vector

Implements dftfe::operatorDFTClass< memorySpace >.

◆ overlapInverseMatrixTimesX() [2/2]

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::overlapInverseMatrixTimesX ( dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > & src,
const double scalarOinvX,
const double scalarY,
const double scalarX,
dftfe::linearAlgebra::MultiVector< dataTypes::numberFP32, memorySpace > & dst )
virtual

Computing Y = scalarOinvX*O^{-1}X + scalarX*X + scalarY*Y for a given X and Y in Reduced precision.

Parameters
srcX vector
scalarOinvXscalar for O^{-1}X
scalarYscalar for Y
scalarXscalar for X
dstY vector

Implements dftfe::operatorDFTClass< memorySpace >.

◆ overlapMatrixTimesX()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::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 )
virtual

Computing Y = scalarOX*OX + scalarX*X + scalarY*Y for a given X and Y in full precision.

Parameters
srcX vector
scalarHXscalar for OX
scalarYscalar for Y
scalarXscalar for X
dstY vector
useApproximateMatrixEntriesflag to use approximate overlap matrix

Implements dftfe::operatorDFTClass< memorySpace >.

◆ reinitkPointSpinIndex()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::reinitkPointSpinIndex ( const unsigned int kPointIndex,
const unsigned int spinIndex )
virtual

sets the data member to appropriate kPoint and spin Index

Parameters
kPointIndexk-point Index to set

Implements dftfe::operatorDFTClass< memorySpace >.

◆ reinitNumberWavefunctions()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::reinitNumberWavefunctions ( const unsigned int numWfc)
virtual

◆ resetExtPotHamFlag()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::resetExtPotHamFlag ( )

◆ resetKohnShamOp()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::resetKohnShamOp ( )

◆ setVEff()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::setVEff ( const std::vector< dftfe::utils::MemoryStorage< double, dftfe::utils::MemorySpace::HOST > > & vKS_quadValues,
const unsigned int spinIndex )

Sets the V-eff potential.

Parameters
vKS_quadValuesthe input V-KS values stored at the quadrature points
spinIndexspin index

◆ setVEffExternalPotCorrToZero()

template<dftfe::utils::MemorySpace memorySpace>
void dftfe::KohnShamHamiltonianOperator< memorySpace >::setVEffExternalPotCorrToZero ( )

Member Data Documentation

◆ computing_timer

template<dftfe::utils::MemorySpace memorySpace>
dealii::TimerOutput dftfe::KohnShamHamiltonianOperator< memorySpace >::computing_timer
private

◆ d_basisOperationsPtr

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr< dftfe::basis::FEBasisOperations<dataTypes::number, double, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_basisOperationsPtr
private

◆ d_basisOperationsPtrHost

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr< dftfe::basis::FEBasisOperations<dataTypes::number, double, dftfe::utils::MemorySpace::HOST> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_basisOperationsPtrHost
private

◆ d_BLASWrapperPtr

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<dftfe::linearAlgebra::BLASWrapper<memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_BLASWrapperPtr
private

◆ d_cellHamiltonianMatrix

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dftfe::utils::MemoryStorage<dataTypes::number, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellHamiltonianMatrix
private

◆ d_cellHamiltonianMatrixExtPot

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<double, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellHamiltonianMatrixExtPot
private

◆ d_cellHamiltonianMatrixSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dftfe::utils::MemoryStorage<dataTypes::numberFP32, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellHamiltonianMatrixSinglePrec
private

◆ d_cellsBlockSizeHamiltonianConstruction

template<dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellsBlockSizeHamiltonianConstruction
private

◆ d_cellsBlockSizeHX

template<dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellsBlockSizeHX
private

◆ d_cellWaveFunctionMatrixDst

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<dataTypes::number, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellWaveFunctionMatrixDst
private

◆ d_cellWaveFunctionMatrixDstSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<dataTypes::numberFP32, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellWaveFunctionMatrixDstSinglePrec
private

◆ d_cellWaveFunctionMatrixSrc

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<dataTypes::number, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellWaveFunctionMatrixSrc
private

◆ d_cellWaveFunctionMatrixSrcSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<dataTypes::numberFP32, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_cellWaveFunctionMatrixSrcSinglePrec
private

◆ d_densityQuadratureID

template<dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_densityQuadratureID
private

◆ d_dftParamsPtr

template<dftfe::utils::MemorySpace memorySpace>
dftParameters* dftfe::KohnShamHamiltonianOperator< memorySpace >::d_dftParamsPtr
private

◆ d_dstNonLocalTemp

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_dstNonLocalTemp
private

◆ d_dstNonLocalTempSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_dstNonLocalTempSinglePrec
private

◆ d_excManagerPtr

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<excManager<memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_excManagerPtr
private

◆ d_feOrderPlusOneQuadratureID

template<dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_feOrderPlusOneQuadratureID
private

◆ d_HamiltonianIndex

template<dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_HamiltonianIndex
private

◆ d_hubbardClassPtr

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<hubbard<dataTypes::number, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_hubbardClassPtr
private

◆ d_invJacderExcWithSigmaTimesGradRhoJxW

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<double, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_invJacderExcWithSigmaTimesGradRhoJxW
private

◆ d_invJacKPointTimesJxW

template<dftfe::utils::MemorySpace memorySpace>
std::vector<dftfe::utils::MemoryStorage<double, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_invJacKPointTimesJxW
private

◆ d_isExternalPotCorrHamiltonianComputed

template<dftfe::utils::MemorySpace memorySpace>
bool dftfe::KohnShamHamiltonianOperator< memorySpace >::d_isExternalPotCorrHamiltonianComputed
private

◆ d_kPointCoordinates

template<dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_kPointCoordinates
private

◆ d_kPointIndex

template<dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_kPointIndex
private

◆ d_kPointWeights

template<dftfe::utils::MemorySpace memorySpace>
std::vector<double> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_kPointWeights
private

◆ d_lpspQuadratureID

template<dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_lpspQuadratureID
private

◆ d_mapNodeIdToProcId

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<dftfe::global_size_type, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_mapNodeIdToProcId
private

◆ d_mpiCommDomain

template<dftfe::utils::MemorySpace memorySpace>
const MPI_Comm dftfe::KohnShamHamiltonianOperator< memorySpace >::d_mpiCommDomain
private

◆ d_mpiCommParent

template<dftfe::utils::MemorySpace memorySpace>
const MPI_Comm dftfe::KohnShamHamiltonianOperator< memorySpace >::d_mpiCommParent
private

◆ d_nOMPThreads

template<dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_nOMPThreads
private

◆ d_numVectorsInternal

template<dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_numVectorsInternal
private

◆ d_oncvClassPtr

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<dftfe::oncvClass<dataTypes::number, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_oncvClassPtr
private

◆ d_ONCVnonLocalOperator

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr< AtomicCenteredNonLocalOperator<dataTypes::number, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_ONCVnonLocalOperator
private

◆ d_ONCVnonLocalOperatorSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr< AtomicCenteredNonLocalOperator<dataTypes::numberFP32, memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::d_ONCVnonLocalOperatorSinglePrec
private

◆ d_ONCVNonLocalProjectorTimesVectorBlock

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_ONCVNonLocalProjectorTimesVectorBlock
private

◆ d_ONCVNonLocalProjectorTimesVectorBlockSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_ONCVNonLocalProjectorTimesVectorBlockSinglePrec
private

◆ d_spinIndex

template<dftfe::utils::MemorySpace memorySpace>
unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::d_spinIndex
private

◆ d_srcNonLocalTemp

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_srcNonLocalTemp
private

◆ d_srcNonLocalTempSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_srcNonLocalTempSinglePrec
private

◆ d_tempBlockVectorOverlapInvX

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::number, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_tempBlockVectorOverlapInvX
private

◆ d_tempBlockVectorOverlapInvXSinglePrec

template<dftfe::utils::MemorySpace memorySpace>
dftfe::linearAlgebra::MultiVector<dataTypes::numberFP32, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_tempBlockVectorOverlapInvXSinglePrec
private

◆ d_useHubbard

template<dftfe::utils::MemorySpace memorySpace>
bool dftfe::KohnShamHamiltonianOperator< memorySpace >::d_useHubbard
private

◆ d_VeffExtPotJxW

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<double, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_VeffExtPotJxW
private

◆ d_VeffJxW

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<double, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::d_VeffJxW
private

◆ inverseMassVectorScaledConstraintsNoneDataInfoPtr

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<dftUtils::constraintMatrixInfo<memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::inverseMassVectorScaledConstraintsNoneDataInfoPtr
private

◆ inverseSqrtMassVectorScaledConstraintsNoneDataInfoPtr

template<dftfe::utils::MemorySpace memorySpace>
std::shared_ptr<dftUtils::constraintMatrixInfo<memorySpace> > dftfe::KohnShamHamiltonianOperator< memorySpace >::inverseSqrtMassVectorScaledConstraintsNoneDataInfoPtr
private

◆ n_mpi_processes

template<dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::n_mpi_processes
private

◆ pcout

template<dftfe::utils::MemorySpace memorySpace>
dealii::ConditionalOStream dftfe::KohnShamHamiltonianOperator< memorySpace >::pcout
private

◆ tempHamMatrixImagBlock

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<double, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::tempHamMatrixImagBlock
private

◆ tempHamMatrixRealBlock

template<dftfe::utils::MemorySpace memorySpace>
dftfe::utils::MemoryStorage<double, memorySpace> dftfe::KohnShamHamiltonianOperator< memorySpace >::tempHamMatrixRealBlock
private

◆ this_mpi_process

template<dftfe::utils::MemorySpace memorySpace>
const unsigned int dftfe::KohnShamHamiltonianOperator< memorySpace >::this_mpi_process
private

The documentation for this class was generated from the following file: