Main class for parallel PEXSI. More...
#include <ppexsi.hpp>
Public Member Functions | |
PPEXSIData (MPI_Comm comm, Int numProcRow, Int numProcCol, Int outputFileIndex) | |
void | LoadRealSymmetricMatrix (Int nrows, Int nnz, Int nnzLocal, Int numColLocal, Int *colptrLocal, Int *rowindLocal, Real *HnzvalLocal, Int isSIdentity, Real *SnzvalLocal, Int verbosity) |
void | SymbolicFactorizeRealSymmetricMatrix (std::string ColPerm, Int numProcSymbFact, Int verbosity) |
Symbolically factorize the loaded matrices for real arithmetic factorization and selected inversion. More... | |
void | SymbolicFactorizeComplexSymmetricMatrix (std::string ColPerm, Int numProcSymbFact, Int verbosity) |
Symbolically factorize the loaded matrices for complex arithmetic factorization and selected inversion. More... | |
void | SelInvRealSymmetricMatrix (double *AnzvalLocal, Int verbosity, double *AinvnzvalLocal) |
void | SelInvComplexSymmetricMatrix (double *AnzvalLocal, Int verbosity, double *AinvnzvalLocal) |
void | CalculateNegativeInertiaReal (const std::vector< Real > &shiftVec, std::vector< Real > &inertiaVec, Int verbosity) |
Compute the negative inertia (the number of eigenvalues below a shift) for real symmetric matrices. The factorization. More... | |
void | CalculateFermiOperatorReal (Int numPole, Real temperature, Real gap, Real deltaE, Real mu, Real numElectronExact, Real numElectronTolerance, Int verbosity, Real &numElectron, Real &numElectronDrvMu) |
Compute the Fermi operator for a given chemical potential for real symmetric matrices. More... | |
void | DFTDriver (Real numElectronExact, Real temperature, Real gap, Real deltaE, Int numPole, Int isInertiaCount, Int maxPEXSIIter, Real muMin0, Real muMax0, Real mu0, Real muInertiaTolerance, Real muInertiaExpansion, Real muPEXSISafeGuard, Real numElectronPEXSITolerance, Int matrixType, Int isSymbolicFactorize, Int ordering, Int numProcSymbFact, Int verbosity, Real &muPEXSI, Real &numElectronPEXSI, Real &muMinInertia, Real &muMaxInertia, Int &numTotalInertiaIter, Int &numTotalPEXSIIter) |
Main driver for solving KSDFT. | |
const GridType * | GridPole () const |
const DistSparseMatrix< Real > & | RhoRealMat () const |
Density matrix. More... | |
const DistSparseMatrix< Real > & | EnergyDensityRealMat () const |
Energy density matrix. More... | |
const DistSparseMatrix< Real > & | FreeEnergyDensityRealMat () const |
Total Helmholtz free energy matrix (band energy part only). More... | |
Real | TotalEnergyH () const |
Real | TotalEnergyS () const |
Real | TotalFreeEnergy () const |
Main class for parallel PEXSI.
void PEXSI::PPEXSIData::CalculateFermiOperatorReal | ( | Int | numPole, |
Real | temperature, | ||
Real | gap, | ||
Real | deltaE, | ||
Real | mu, | ||
Real | numElectronExact, | ||
Real | numElectronTolerance, | ||
Int | verbosity, | ||
Real & | numElectron, | ||
Real & | numElectronDrvMu | ||
) |
Compute the Fermi operator for a given chemical potential for real symmetric matrices.
This routine also computes the single particle density matrix, the Helmholtz free energy density matrix, and the energy density matrix (for computing the Pulay force) simultaneously. These matrices can be called later via member functions DensityMatrix, FreeEnergyDensityMatrix, EnergyDensityMatrix.
[in] | numPole | Number of poles for the pole expansion |
[in] | temperature | Temperature |
[in] | gap | Band gap |
[in] | deltaE | Upperbound of the spectrum width |
[in] | mu | Initial guess of chemical potential. |
[in] | numElectronExact | Exact number of electrons. |
[in] | numElectronTolerance | Tolerance for the number of electrons. This is just used to discard some poles in the pole expansion. |
[in] | verbosity | The level of output information.
|
[out] | numElectron | The number of electron calculated at mu. |
[out] | numElectronDrvMu | The derivative of the number of electron calculated with respect to the chemical potential at mu. |
void PEXSI::PPEXSIData::CalculateNegativeInertiaReal | ( | const std::vector< Real > & | shiftVec, |
std::vector< Real > & | inertiaVec, | ||
Int | verbosity | ||
) |
Compute the negative inertia (the number of eigenvalues below a shift) for real symmetric matrices. The factorization.
uses real arithemetic factorization routine. This subroutine computes the negative inertia of the matrix
I = H - shift * S
where I is the same as the number of eigenvalues lambda for
H x = lambda S x
with lambda < shift according to the Sylvester's law of inertia.
[in] | shiftVec | Shift vectors. |
[out] | inertiaVec | Negative inertia count, the same size as shiftVec. |
[in] | HMat | Hamiltonian matrix saved in distributed compressed sparse column format. See DistSparseMatrix. |
[in] | SMat | Overlap matrix saved in distributed compressed sparse column format. See DistSparseMatrix. |
Note: If SMat.size == 0, SMat is treated as an identity matrix.
[in] | verbosity | The level of output information.
|
|
inline |
Energy density matrix.
Can be used to estimate the total band energy via Tr[EDM*S] or the force, including the Hellman-Feynman force and the Pulay force.
|
inline |
Total Helmholtz free energy matrix (band energy part only).
The Helmholtz free energy is computed by Tr[rho_f*H]. For more information see Alavi, A., Kohanoff, J., Parrinello, M., & Frenkel, D. (1994). Ab initio molecular dynamics with excited electrons. Physical review letters, 73(19), 2599–2602.
|
inline |
Density matrix.
Can be used to estimate the number of electrons by Tr[DM*S] or the band energy via Tr[DM*H]
void PEXSI::PPEXSIData::SymbolicFactorizeComplexSymmetricMatrix | ( | std::string | ColPerm, |
Int | numProcSymbFact, | ||
Int | verbosity | ||
) |
Symbolically factorize the loaded matrices for complex arithmetic factorization and selected inversion.
The symbolic information is saved internally at luComplexMat_ and PMComplexMat_.
[in] | ColPerm | Permutation method used for SuperLU_DIST |
[in] | numProcSymbFact | Number of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH. |
[in] | verbosity | The level of output information.
|
void PEXSI::PPEXSIData::SymbolicFactorizeRealSymmetricMatrix | ( | std::string | ColPerm, |
Int | numProcSymbFact, | ||
Int | verbosity | ||
) |
Symbolically factorize the loaded matrices for real arithmetic factorization and selected inversion.
The symbolic information is saved internally at luRealMat_ and PMRealMat_.
[in] | ColPerm | Permutation method used for SuperLU_DIST |
[in] | numProcSymbFact | Number of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH. |
[in] | verbosity | The level of output information.
|