DFT-FE 1.1.0-pre
Density Functional Theory With Finite-Elements
Loading...
Searching...
No Matches
dftfe::linearAlgebraOperations Namespace Reference

Contains linear algebra functions used in the implementation of an eigen solver. More...

Namespaces

namespace  internal
 Contains internal functions used in linearAlgebraOperations.
 

Functions

void inverse (double *A, int N)
 Compute inverse of serial matrix using LAPACK LU factorization.
 
void inverse (std::complex< double > *A, int N)
 Compute inverse of serial matrix using LAPACK LU factorization.
 
template<typename T, dftfe::utils::MemorySpace memorySpace>
std::pair< double, double > generalisedLanczosLowerUpperBoundEigenSpectrum (const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > &BLASWrapperPtr, operatorDFTClass< memorySpace > &operatorMatrix, dftfe::linearAlgebra::MultiVector< T, memorySpace > &X, dftfe::linearAlgebra::MultiVector< T, memorySpace > &Y, dftfe::linearAlgebra::MultiVector< T, memorySpace > &Z, dftfe::linearAlgebra::MultiVector< T, memorySpace > &tempVec, const dftParameters &dftParams)
 Calculates an estimate of lower and upper bounds of a matrix using k-step Generalised Lanczos method. Algo is present in PAW PRB paper.
 
template<typename T, dftfe::utils::MemorySpace memorySpace>
void chebyshevFilter (operatorDFTClass< memorySpace > &operatorMatrix, dftfe::linearAlgebra::MultiVector< T, memorySpace > &X, dftfe::linearAlgebra::MultiVector< T, memorySpace > &Y, const unsigned int m, const double a, const double b, const double a0)
 Apply Chebyshev filter to a given subspace.
 
template<typename T1, typename T2, dftfe::utils::MemorySpace memorySpace>
void reformulatedChebyshevFilter (const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > &BLASWrapperPtr, operatorDFTClass< memorySpace > &operatorMatrix, dftfe::linearAlgebra::MultiVector< T1, memorySpace > &X, dftfe::linearAlgebra::MultiVector< T1, memorySpace > &Y, dftfe::linearAlgebra::MultiVector< T2, memorySpace > &Residual, dftfe::linearAlgebra::MultiVector< T2, memorySpace > &ResidualNew, std::vector< double > eigenvalues, const unsigned int m, const double a, const double b, const double a0, const bool approxOverlapMatrix)
 Apply Residual based Chebyshev filter to a given subspace.
 
template<typename T>
void gramSchmidtOrthogonalization (T *X, const unsigned int numberComponents, const unsigned int numberDofs, const MPI_Comm &mpiComm)
 Orthogonalize given subspace using GramSchmidt orthogonalization.
 
unsigned int lowdenOrthogonalization (std::vector< dataTypes::number > &X, const unsigned int numberComponents, const MPI_Comm &mpiComm, const dftParameters &dftParams)
 Orthogonalize given subspace using Lowden orthogonalization for double data-type (serial version using LAPACK)
 
template<typename T>
unsigned int pseudoGramSchmidtOrthogonalization (elpaScalaManager &elpaScala, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *X, const unsigned int numberComponents, const unsigned int numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiCommDomain, const bool useMixedPrec, const dftParameters &dftParams)
 Orthogonalize given subspace using Pseudo-Gram-Schmidt orthogonalization (serial version using LAPACK, parallel version using ScaLAPACK)
 
template<typename T>
void rayleighRitzGEP (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, elpaScalaManager &elpaScala, T *X, const unsigned int numberComponents, const unsigned int numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiCommDomain, std::vector< double > &eigenValues, const bool useMixedPrec, const dftParameters &dftParams)
 Compute Rayleigh-Ritz projection (serial version using LAPACK, parallel version using ScaLAPACK)
 
template<typename T>
void rayleighRitz (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, elpaScalaManager &elpaScala, T *X, const unsigned int numberComponents, const unsigned int numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &interBandGroupComm, const MPI_Comm &mpiCommDomain, std::vector< double > &eigenValues, const dftParameters &dftParams, const bool doCommAfterBandParal=true)
 Compute Rayleigh-Ritz projection (serial version using LAPACK, parallel version using ScaLAPACK)
 
template<typename T>
void computeEigenResidualNorm (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *X, const std::vector< double > &eigenValues, const unsigned int numberComponents, const unsigned int numberDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, std::vector< double > &residualNorm, const dftParameters &dftParams)
 Compute residual norm associated with eigenValue problem of the given operator.
 
template<typename T>
void densityMatrixEigenBasisFirstOrderResponse (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, T *X, const unsigned int N, const unsigned int numberLocalDofs, const MPI_Comm &mpiCommParent, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const std::vector< double > &eigenValues, const double fermiEnergy, std::vector< double > &densityMatDerFermiEnergy, elpaScalaManager &elpaScala, const dftParameters &dftParams)
 Compute first order response in density matrix with respect to perturbation in the Hamiltonian. Perturbation is computed in the eigenbasis.
 
void XtHX (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const unsigned int numberComponents, const unsigned int numberLocalDofs, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, std::vector< dataTypes::number > &ProjHam)
 Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.
 
void XtHX (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const unsigned int numberComponents, const unsigned int numberLocalDofs, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, dftfe::ScaLAPACKMatrix< dataTypes::number > &projHamPar, const bool onlyHPrimePartForFirstOrderDensityMatResponse=false)
 Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.
 
void XtOX (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const unsigned int numberComponents, const unsigned int numberLocalDofs, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, dftfe::ScaLAPACKMatrix< dataTypes::number > &projOverlapPar)
 Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.
 
void XtHXMixedPrec (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const unsigned int totalNumberComponents, const unsigned int singlePrecComponents, const unsigned int numberLocalDofs, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, dftfe::ScaLAPACKMatrix< dataTypes::number > &projHamPar, const bool onlyHPrimePartForFirstOrderDensityMatResponse=false)
 Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.
 
void XtOXMixedPrec (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const unsigned int totalNumberComponents, const unsigned int singlePrecComponents, const unsigned int numberLocalDofs, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, dftfe::ScaLAPACKMatrix< dataTypes::number > &projOverlapPar)
 Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.
 
void XtHXXtOX (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const unsigned int numberComponents, const unsigned int numberLocalDofs, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, dftfe::ScaLAPACKMatrix< dataTypes::number > &projHamPar, dftfe::ScaLAPACKMatrix< dataTypes::number > &projOverlapPar, const bool onlyHPrimePartForFirstOrderDensityMatResponse=false)
 Computes the projection of Hamiltonian and Overlap with only a single extraction. Single extraction will be beneficial in full M, PAW cases.
 
void XtHXXtOXMixedPrec (operatorDFTClass< dftfe::utils::MemorySpace::HOST > &operatorMatrix, const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > &BLASWrapperPtr, const dataTypes::number *X, const unsigned int totalNumberComponents, const unsigned int singlePrecComponents, const unsigned int numberLocalDofs, const std::shared_ptr< const dftfe::ProcessGrid > &processGrid, const MPI_Comm &mpiCommDomain, const MPI_Comm &interBandGroupComm, const dftParameters &dftParams, dftfe::ScaLAPACKMatrix< dataTypes::number > &projHamPar, dftfe::ScaLAPACKMatrix< dataTypes::number > &projOverlapPar, const bool onlyHPrimePartForFirstOrderDensityMatResponse=false)
 Computes the projection of Hamiltonian and Overlap with only a single extraction with mixed precision. Single extraction will be beneficial in full M, PAW cases. THe projected Hamiltonain has full precision along blocks of diagonal and for states greater than Ncore THe projected Overlap will be of full precision along the diagonal.
 

Detailed Description

Contains linear algebra functions used in the implementation of an eigen solver.

Author
Phani Motamarri, Sambit Das

Function Documentation

◆ chebyshevFilter()

template<typename T, dftfe::utils::MemorySpace memorySpace>
void dftfe::linearAlgebraOperations::chebyshevFilter ( operatorDFTClass< memorySpace > & operatorMatrix,
dftfe::linearAlgebra::MultiVector< T, memorySpace > & X,
dftfe::linearAlgebra::MultiVector< T, memorySpace > & Y,
const unsigned int m,
const double a,
const double b,
const double a0 )

Apply Chebyshev filter to a given subspace.

Parameters
[in]operatorMatrixAn object which has access to the given matrix
[in,out]XGiven subspace as a dealii array representing multiple fields as a flattened array. In-place update of the given subspace.
[in]numberComponentsNumber of multiple-fields
[in]mChebyshev polynomial degree
[in]alower bound of unwanted spectrum
[in]bupper bound of unwanted spectrum
[in]a0lower bound of wanted spectrum

◆ computeEigenResidualNorm()

template<typename T>
void dftfe::linearAlgebraOperations::computeEigenResidualNorm ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
T * X,
const std::vector< double > & eigenValues,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
std::vector< double > & residualNorm,
const dftParameters & dftParams )

Compute residual norm associated with eigenValue problem of the given operator.

Parameters
[in]operatorMatrixAn object which has access to the given matrix
[in]XGiven subspace as STL vector of dealii vectors
[in]eigenValueseigenValues of the operator
[in]mpiCommParentparent mpi communicator
[in]mpiCommDomaindomain decomposition communicator
[out]residualNormsof the eigen Value problem

◆ densityMatrixEigenBasisFirstOrderResponse()

template<typename T>
void dftfe::linearAlgebraOperations::densityMatrixEigenBasisFirstOrderResponse ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
T * X,
const unsigned int N,
const unsigned int numberLocalDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const std::vector< double > & eigenValues,
const double fermiEnergy,
std::vector< double > & densityMatDerFermiEnergy,
elpaScalaManager & elpaScala,
const dftParameters & dftParams )

Compute first order response in density matrix with respect to perturbation in the Hamiltonian. Perturbation is computed in the eigenbasis.

◆ generalisedLanczosLowerUpperBoundEigenSpectrum()

template<typename T, dftfe::utils::MemorySpace memorySpace>
std::pair< double, double > dftfe::linearAlgebraOperations::generalisedLanczosLowerUpperBoundEigenSpectrum ( const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > & BLASWrapperPtr,
operatorDFTClass< memorySpace > & operatorMatrix,
dftfe::linearAlgebra::MultiVector< T, memorySpace > & X,
dftfe::linearAlgebra::MultiVector< T, memorySpace > & Y,
dftfe::linearAlgebra::MultiVector< T, memorySpace > & Z,
dftfe::linearAlgebra::MultiVector< T, memorySpace > & tempVec,
const dftParameters & dftParams )

Calculates an estimate of lower and upper bounds of a matrix using k-step Generalised Lanczos method. Algo is present in PAW PRB paper.

Parameters
operatorMatrixAn object which has access to the given matrix
vectA dummy vector
Returns
std::pair<double,double> An estimate of the lower and upper bound of the given matrix

◆ gramSchmidtOrthogonalization()

template<typename T>
void dftfe::linearAlgebraOperations::gramSchmidtOrthogonalization ( T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiComm )

Orthogonalize given subspace using GramSchmidt orthogonalization.

Parameters
[in,out]XGiven subspace as flattened array of multi-vectors. In-place update of the given subspace
[in]numberComponentsNumber of multiple-fields
[in]mpiCommglobal communicator

◆ inverse() [1/2]

void dftfe::linearAlgebraOperations::inverse ( double * A,
int N )

Compute inverse of serial matrix using LAPACK LU factorization.

◆ inverse() [2/2]

void dftfe::linearAlgebraOperations::inverse ( std::complex< double > * A,
int N )

Compute inverse of serial matrix using LAPACK LU factorization.

◆ lowdenOrthogonalization()

unsigned int dftfe::linearAlgebraOperations::lowdenOrthogonalization ( std::vector< dataTypes::number > & X,
const unsigned int numberComponents,
const MPI_Comm & mpiComm,
const dftParameters & dftParams )

Orthogonalize given subspace using Lowden orthogonalization for double data-type (serial version using LAPACK)

Parameters
[in,out]XGiven subspace as flattened array of multi-vectors. In-place update of the given subspace
[in]numberComponentsNumber of multiple-fields
[in]mpiCommglobal communicator
Returns
flag indicating success/failure. 1 for failure, 0 for success

◆ pseudoGramSchmidtOrthogonalization()

template<typename T>
unsigned int dftfe::linearAlgebraOperations::pseudoGramSchmidtOrthogonalization ( elpaScalaManager & elpaScala,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
const bool useMixedPrec,
const dftParameters & dftParams )

Orthogonalize given subspace using Pseudo-Gram-Schmidt orthogonalization (serial version using LAPACK, parallel version using ScaLAPACK)

Parameters
[in,out]XGiven subspace as flattened array of multi-vectors. In-place update of the given subspace
[in]numberComponentsNumber of multiple-fields
[in]mpiCommParentparent communicator
[in]interBandGroupComminterpool communicator for parallelization over band groups
[in]mpiCommdomain decomposition communicator
Returns
flag indicating success/failure. 1 for failure, 0 for success

◆ rayleighRitz()

template<typename T>
void dftfe::linearAlgebraOperations::rayleighRitz ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
elpaScalaManager & elpaScala,
T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
std::vector< double > & eigenValues,
const dftParameters & dftParams,
const bool doCommAfterBandParal = true )

Compute Rayleigh-Ritz projection (serial version using LAPACK, parallel version using ScaLAPACK)

Parameters
[in]operatorMatrixAn object which has access to the given matrix
[in,out]XGiven subspace as flattened array of multi-vectors. In-place rotated subspace
[in]numberComponentsNumber of vectors
[in]mpiCommParentparent mpi communicator
[in]interBandGroupComminterpool communicator for parallelization over band groups
[in]mpiCommDomaindomain decomposition communicator
[out]eigenValuesof the Projected Hamiltonian

◆ rayleighRitzGEP()

template<typename T>
void dftfe::linearAlgebraOperations::rayleighRitzGEP ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
elpaScalaManager & elpaScala,
T * X,
const unsigned int numberComponents,
const unsigned int numberDofs,
const MPI_Comm & mpiCommParent,
const MPI_Comm & interBandGroupComm,
const MPI_Comm & mpiCommDomain,
std::vector< double > & eigenValues,
const bool useMixedPrec,
const dftParameters & dftParams )

Compute Rayleigh-Ritz projection (serial version using LAPACK, parallel version using ScaLAPACK)

Parameters
[in]operatorMatrixAn object which has access to the given matrix
[in,out]XGiven subspace as flattened array of multi-vectors. In-place rotated subspace
[in]numberComponentsNumber of vectors
[in]mpiCommDomainparent communicator
[in]interBandGroupComminterpool communicator for parallelization over band groups
[in]mpiCommDomaindomain decomposition communicator
[out]eigenValuesof the Projected Hamiltonian

◆ reformulatedChebyshevFilter()

template<typename T1, typename T2, dftfe::utils::MemorySpace memorySpace>
void dftfe::linearAlgebraOperations::reformulatedChebyshevFilter ( const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< memorySpace > > & BLASWrapperPtr,
operatorDFTClass< memorySpace > & operatorMatrix,
dftfe::linearAlgebra::MultiVector< T1, memorySpace > & X,
dftfe::linearAlgebra::MultiVector< T1, memorySpace > & Y,
dftfe::linearAlgebra::MultiVector< T2, memorySpace > & Residual,
dftfe::linearAlgebra::MultiVector< T2, memorySpace > & ResidualNew,
std::vector< double > eigenvalues,
const unsigned int m,
const double a,
const double b,
const double a0,
const bool approxOverlapMatrix )

Apply Residual based Chebyshev filter to a given subspace.

Parameters
[in]operatorMatrixAn object which has access to the given matrix
[in,out]XGiven subspace as a dealii array representing multiple fields as a flattened array. In-place update of the given subspace.
[in]eigenvaluesestimate of eigenvalues, usually the eigenvalues from previous pass
[in]mChebyshev polynomial degree
[in]alower bound of unwanted spectrum
[in]bupper bound of unwanted spectrum
[in]a0lower bound of wanted spectrum
[in]approxOverlapMatrixto use approximate overlap matrix while computing initial residual

◆ XtHX() [1/2]

void dftfe::linearAlgebraOperations::XtHX ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
const dataTypes::number * X,
const unsigned int numberComponents,
const unsigned int numberLocalDofs,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
std::vector< dataTypes::number > & ProjHam )

Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.

Parameters
XVector of Vectors containing multi-wavefunction fields
numberComponentsnumber of wavefunctions associated with a given node
ProjMatrixprojected small matrix

◆ XtHX() [2/2]

void dftfe::linearAlgebraOperations::XtHX ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
const dataTypes::number * X,
const unsigned int numberComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr< const dftfe::ProcessGrid > & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projHamPar,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false )

Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.

Parameters
XVector of Vectors containing multi-wavefunction fields
numberComponentsnumber of wavefunctions associated with a given node
processGridtwo-dimensional processor grid corresponding to the parallel projHamPar
projHamParparallel ScaLAPACKMatrix which stores the computed projection of the operation into the given subspace

◆ XtHXMixedPrec()

void dftfe::linearAlgebraOperations::XtHXMixedPrec ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
const dataTypes::number * X,
const unsigned int totalNumberComponents,
const unsigned int singlePrecComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr< const dftfe::ProcessGrid > & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projHamPar,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false )

Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.

Parameters
XVector of Vectors containing multi-wavefunction fields
totalNumberComponentsnumber of wavefunctions associated with a given node
singlePrecComponentsnumber of wavecfuntions starting from the first for which the project Hamiltionian block will be computed in single procession. However the cross blocks will still be computed in double precision.
processGridtwo-dimensional processor grid corresponding to the parallel projHamPar
projHamParparallel ScaLAPACKMatrix which stores the computed projection of the operation into the given subspace

◆ XtHXXtOX()

void dftfe::linearAlgebraOperations::XtHXXtOX ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
const dataTypes::number * X,
const unsigned int numberComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr< const dftfe::ProcessGrid > & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projHamPar,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projOverlapPar,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false )

Computes the projection of Hamiltonian and Overlap with only a single extraction. Single extraction will be beneficial in full M, PAW cases.

Parameters
XVector of Vectors containing multi-wavefunction fields
numberComponentsnumber of wavefunctions associated with a given node
numberLocalDofsnumber of DOFs owned in the current procesor
processGridtwo-dimensional processor grid corresponding to the parallel projHamPar
projHamParparallel ScaLAPACKMatrix which stores the computed projection of the Hamiltonian into the given subspace
projOverlapParparallel ScaLAPACKMatrix which stores the computed projection of the Overlap into the given subspace

◆ XtHXXtOXMixedPrec()

void dftfe::linearAlgebraOperations::XtHXXtOXMixedPrec ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
const dataTypes::number * X,
const unsigned int totalNumberComponents,
const unsigned int singlePrecComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr< const dftfe::ProcessGrid > & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projHamPar,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projOverlapPar,
const bool onlyHPrimePartForFirstOrderDensityMatResponse = false )

Computes the projection of Hamiltonian and Overlap with only a single extraction with mixed precision. Single extraction will be beneficial in full M, PAW cases. THe projected Hamiltonain has full precision along blocks of diagonal and for states greater than Ncore THe projected Overlap will be of full precision along the diagonal.

Parameters
XVector of Vectors containing multi-wavefunction fields
numberComponentsnumber of wavefunctions associated with a given node
numberLocalDofsnumber of DOFs owned in the current procesor
processGridtwo-dimensional processor grid corresponding to the parallel projHamPar
projHamParparallel ScaLAPACKMatrix which stores the computed projection of the Hamiltonian into the given subspace
projOverlapParparallel ScaLAPACKMatrix which stores the computed projection of the Overlap into the given subspace

◆ XtOX()

void dftfe::linearAlgebraOperations::XtOX ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
const dataTypes::number * X,
const unsigned int numberComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr< const dftfe::ProcessGrid > & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projOverlapPar )

Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.

Parameters
XVector of Vectors containing multi-wavefunction fields
numberComponentsnumber of wavefunctions associated with a given node
processGridtwo-dimensional processor grid corresponding to the parallel projHamPar
projHamParparallel ScaLAPACKMatrix which stores the computed projection of the operation into the given subspace

◆ XtOXMixedPrec()

void dftfe::linearAlgebraOperations::XtOXMixedPrec ( operatorDFTClass< dftfe::utils::MemorySpace::HOST > & operatorMatrix,
const std::shared_ptr< dftfe::linearAlgebra::BLASWrapper< dftfe::utils::MemorySpace::HOST > > & BLASWrapperPtr,
const dataTypes::number * X,
const unsigned int totalNumberComponents,
const unsigned int singlePrecComponents,
const unsigned int numberLocalDofs,
const std::shared_ptr< const dftfe::ProcessGrid > & processGrid,
const MPI_Comm & mpiCommDomain,
const MPI_Comm & interBandGroupComm,
const dftParameters & dftParams,
dftfe::ScaLAPACKMatrix< dataTypes::number > & projOverlapPar )

Compute projection of the operator into a subspace spanned by a given orthogonal basis HProjConj=X^{T}*HConj*XConj.

Parameters
XVector of Vectors containing multi-wavefunction fields
totalNumberComponentsnumber of wavefunctions associated with a given node
singlePrecComponentsnumber of wavecfuntions starting from the first for which the project Hamiltionian block will be computed in single procession. However the cross blocks will still be computed in double precision.
processGridtwo-dimensional processor grid corresponding to the parallel projHamPar
projHamParparallel ScaLAPACKMatrix which stores the computed projection of the operation into the given subspace