PEXSI
 All Classes Namespaces Files Functions Variables Typedefs Pages
Public Member Functions | List of all members
PEXSI::PPEXSIData Class Reference

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 LoadRealUnsymmetricMatrix (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 SymbolicFactorizeRealUnsymmetricMatrix (std::string ColPerm, std::string RowPerm, Int numProcSymbFact, Int Transpose, double *AnzvalLocal, 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 SymbolicFactorizeComplexUnsymmetricMatrix (std::string ColPerm, std::string RowPerm, Int numProcSymbFact, Int Transpose, double *AnzvalLocal, Int verbosity)
 Symbolically factorize the loaded matrices for complex arithmetic factorization and selected inversion. More...
 
void SelInvRealSymmetricMatrix (double *AnzvalLocal, Int verbosity, double *AinvnzvalLocal)
 
void SelInvRealUnsymmetricMatrix (double *AnzvalLocal, Int verbosity, double *AinvnzvalLocal)
 
void SelInvComplexSymmetricMatrix (double *AnzvalLocal, Int verbosity, double *AinvnzvalLocal)
 
void SelInvComplexUnsymmetricMatrix (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.
 
void CalculateFermiOperatorReal2 (Int numPole, Real temperature, Real gap, Real deltaE, Real numElectronExact, Real numElectronTolerance, Real muMinPEXSI, Real muMaxPEXSI, Int verbosity, Real &mu, Real &numElectron, bool &isPEXSIConverged)
 Compute the Fermi operator and derivied quantities. More...
 
void DFTDriver2 (Real numElectronExact, Real temperature, Real gap, Real deltaE, Int numPole, Int isInertiaCount, Real muMin0, Real muMax0, Real mu0, Real muInertiaTolerance, Real muInertiaExpansion, Real numElectronPEXSITolerance, Int matrixType, Int isSymbolicFactorize, Int ordering, Int numProcSymbFact, Int verbosity, Real &muPEXSI, Real &numElectronPEXSI, Real &muMinInertia, Real &muMaxInertia, Int &numTotalInertiaIter)
 Updated main driver for DFT. This reuses the pole expansion and only performs one PEXSI iteration per SCF step.
 
const GridTypeGridPole () 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
 

Detailed Description

Main class for parallel PEXSI.

Member Function Documentation

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.

Parameters
[in]numPoleNumber of poles for the pole expansion
[in]temperatureTemperature
[in]gapBand gap
[in]deltaEUpperbound of the spectrum width
[in]muInitial guess of chemical potential.
[in]numElectronExactExact number of electrons.
[in]numElectronToleranceTolerance for the number of electrons. This is just used to discard some poles in the pole expansion.
[in]verbosityThe level of output information.
  • = 0 : No output.
  • = 1 : Basic output (default)
  • = 2 : Detailed output.
[out]numElectronThe number of electron calculated at mu.
[out]numElectronDrvMuThe derivative of the number of electron calculated with respect to the chemical potential at mu.
void PEXSI::PPEXSIData::CalculateFermiOperatorReal2 ( Int  numPole,
Real  temperature,
Real  gap,
Real  deltaE,
Real  numElectronExact,
Real  numElectronTolerance,
Real  muMinPEXSI,
Real  muMaxPEXSI,
Int  verbosity,
Real &  mu,
Real &  numElectron,
bool &  isPEXSIConverged 
)

Compute the Fermi operator and derivied quantities.

This routine also updates the chemical potential mu by reusing Green's functions but with updated contour.

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.

Parameters
[in]numPoleNumber of poles for the pole expansion
[in]temperatureTemperature
[in]gapBand gap
[in]deltaEUpperbound of the spectrum width
[in]numElectronExactExact number of electrons.
[in]numElectronToleranceTolerance for the number of electrons. This is just used to discard some poles in the pole expansion.
[in]muMinPEXSIMinimum of the interval for searching mu.
[in]muMaxPEXSIMaximum of the interval for searching mu.
[in]verbosityThe level of output information.
  • = 0 : No output.
  • = 1 : Basic output (default)
  • = 2 : Detailed output.
[in,out]muInitial guess of chemical potential. On return it gives the updated chemical potential within the range of [muMinPEXSI, muMaxPEXSI]
[out]numElectronThe number of electron calculated at mu.
[out]isConvergedWhether the update strategy for finding the chemical potential has converged.
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.

Parameters
[in]shiftVecShift vectors.
[out]inertiaVecNegative inertia count, the same size as shiftVec.
[in]HMatHamiltonian matrix saved in distributed compressed sparse column format. See DistSparseMatrix.
[in]SMatOverlap matrix saved in distributed compressed sparse column format. See DistSparseMatrix.

Note: If SMat.size == 0, SMat is treated as an identity matrix.

Parameters
[in]verbosityThe level of output information.
  • = 0 : No output.
  • = 1 : Basic output (default)
  • = 2 : Detailed output.
const DistSparseMatrix<Real>& PEXSI::PPEXSIData::EnergyDensityRealMat ( ) const
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.

const DistSparseMatrix<Real>& PEXSI::PPEXSIData::FreeEnergyDensityRealMat ( ) const
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.

const DistSparseMatrix<Real>& PEXSI::PPEXSIData::RhoRealMat ( ) const
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_.

Parameters
[in]ColPermPermutation method used for SuperLU_DIST
[in]numProcSymbFactNumber of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH.
[in]verbosityThe level of output information.
  • = 0 : No output.
  • = 1 : Basic output (default)
  • = 2 : Detailed output.
void PEXSI::PPEXSIData::SymbolicFactorizeComplexUnsymmetricMatrix ( std::string  ColPerm,
std::string  RowPerm,
Int  numProcSymbFact,
Int  Transpose,
double *  AnzvalLocal,
Int  verbosity 
)

Symbolically factorize the loaded matrices for complex arithmetic factorization and selected inversion.

The symbolic information is saved internally at luComplexMat_ and PMComplexUnsymMat_.

Parameters
[in]ColPermPermutation method used for SuperLU_DIST
[in]RowPermRow Permutation method used for SuperLU_DIST
[in]numProcSymbFactNumber of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH.
[in]TransposeTODO
[in]AnzvalLocalnon zero values for row permutation
[in]verbosityThe level of output information.
  • = 0 : No output.
  • = 1 : Basic output (default)
  • = 2 : Detailed output.
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_.

Parameters
[in]ColPermPermutation method used for SuperLU_DIST
[in]numProcSymbFactNumber of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH.
[in]verbosityThe level of output information.
  • = 0 : No output.
  • = 1 : Basic output (default)
  • = 2 : Detailed output.
void PEXSI::PPEXSIData::SymbolicFactorizeRealUnsymmetricMatrix ( std::string  ColPerm,
std::string  RowPerm,
Int  numProcSymbFact,
Int  Transpose,
double *  AnzvalLocal,
Int  verbosity 
)

Symbolically factorize the loaded matrices for real arithmetic factorization and selected inversion.

The symbolic information is saved internally at luRealMat_ and PMRealMat_.

Parameters
[in]ColPermPermutation method used for SuperLU_DIST
[in]RowPermRow Permutation method used for SuperLU_DIST
[in]numProcSymbFactNumber of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH.
[in]TransposeTODO
[in]AnzvalLocalnon zero values for row permutation
[in]verbosityThe level of output information.
  • = 0 : No output.
  • = 1 : Basic output (default)
  • = 2 : Detailed output.

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