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 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 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
 
 PPEXSIData (const PEXSI::GridType *g, Int nprow, Int npcol)
 
void Solve (Int numPole, Real temperature, Real numElectronExact, Real gap, Real deltaE, Real &mu, Real &muMin, Real &muMax, const DistSparseMatrix< Real > &HMat, const DistSparseMatrix< Real > &SMat, Int muMaxIter, Real numElectronTolerance, std::string ColPerm, Int numProcSymbFact, bool isFreeEnergyDensityMatrix, bool isEnergyDensityMatrix, bool isDerivativeTMatrix, std::vector< Real > &muList, std::vector< Real > &numElectronList, std::vector< Real > &numElectronDrvMuList, bool &isConverged)
 Solve is the main subroutine for PPEXSI. More...
 
DistSparseMatrix< Real > & DensityMatrix ()
 DensityMatrix returns the single particle density matrix.
 
DistSparseMatrix< Real > & FreeEnergyDensityMatrix ()
 DensityMatrix returns the free energy density matrix.
 
DistSparseMatrix< Real > & EnergyDensityMatrix ()
 DensityMatrix returns the energy density matrix for computing the force.
 
Real CalculateNumElectron (const DistSparseMatrix< Real > &SMat)
 CalculateNumElectron computes the number of electrons given the current density matrix. More...
 
Real CalculateNumElectronDrvMu (const DistSparseMatrix< Real > &SMat)
 CalculateNumElectronDrvMu computes the derivative of the number of electrons with respect to the chemical potential. More...
 
Real CalculateNumElectronDrvT (const DistSparseMatrix< Real > &SMat)
 CalculateNumElectronDrvT computes the derivative of the number of electrons with respect to the temperature (1/beta, in atomic unit). More...
 
Real CalculateTotalEnergy (const DistSparseMatrix< Real > &HMat)
 CalculateTotalEnergy computes the total energy (band energy part only). More...
 
Real CalculateFreeEnergy (const DistSparseMatrix< Real > &SMat)
 CalculateFreeEnergy computes the total Helmholtz free energy (band energy part only). More...
 
Real CalculateForce (const DistSparseMatrix< Real > &HDerivativeMat, const DistSparseMatrix< Real > &SDerivativeMat)
 CalculateFreeEnergy computes the force, including the Hellman-Feynman force and the Pulay force. More...
 
Real EstimateZeroTemperatureChemicalPotential (Real temperature, Real mu, const DistSparseMatrix< Real > &SMat)
 Estimate the chemical potential at zero temperature using quadratic approximation. More...
 
void CalculateNegativeInertia (const std::vector< Real > &shiftVec, std::vector< Real > &inertiaVec, const DistSparseMatrix< Real > &HMat, const DistSparseMatrix< Real > &SMat, std::string ColPerm, Int numProcSymbFact)
 Compute the negative inertia (the number of eigenvalues below a shift). More...
 
void CalculateNegativeInertiaReal (const std::vector< Real > &shiftVec, std::vector< Real > &inertiaVec, const DistSparseMatrix< Real > &HMat, const DistSparseMatrix< Real > &SMat, std::string ColPerm, Int numProcSymbFact)
 Compute the negative inertia (the number of eigenvalues below a shift) with real arithmetic which directly using SuperLU.
 
void EstimateSpectralRadius (Int method, const DistSparseMatrix< Real > &HMat, const DistSparseMatrix< Real > &SMat, const NumVec< Real > &v0, Real tol, Int maxNumIter, Int &numIter, Real &sigma)
 Estimate the spectral radius of the matrix stencil (H,S). More...
 

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.
Real PEXSI::PPEXSIData::CalculateForce ( const DistSparseMatrix< Real > &  HDerivativeMat,
const DistSparseMatrix< Real > &  SDerivativeMat 
)

CalculateFreeEnergy computes the force, including the Hellman-Feynman force and the Pulay force.

Parameters
[in]HDerivativeMatDerivative of the Hamilotian matrix with respect to the atomic position R_{i,j}, i = 1, ..., natoms, j=1,2,3.
[in]HDerivativeMatDerivative of the overlap matrix with respect to the atomic position R_{i,j}, i = 1, ..., natoms, j=1,2,3.
Returns
The force f_{i,j} with respect to the atomic position R_{i,j}, i = 1, ..., natoms, j=1,2,3.
Real PEXSI::PPEXSIData::CalculateFreeEnergy ( const DistSparseMatrix< Real > &  SMat)

CalculateFreeEnergy computes the total Helmholtz free energy (band energy part only).

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.

Parameters
[in]HMatHamilotian matrix.
Returns
The Helmholtz free energy Tr[rho_f H]
void PEXSI::PPEXSIData::CalculateNegativeInertia ( const std::vector< Real > &  shiftVec,
std::vector< Real > &  inertiaVec,
const DistSparseMatrix< Real > &  HMat,
const DistSparseMatrix< Real > &  SMat,
std::string  ColPerm,
Int  numProcSymbFact 
)

Compute the negative inertia (the number of eigenvalues below a shift).

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]ColPermPermutation method used for SuperLU_DIST
[in]numProcSymbFactNumber of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH.
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.
Real PEXSI::PPEXSIData::CalculateNumElectron ( const DistSparseMatrix< Real > &  SMat)

CalculateNumElectron computes the number of electrons given the current density matrix.

Parameters
[in]SMatoverlap matrix.
Returns
The number of electrons Tr[ S]
Real PEXSI::PPEXSIData::CalculateNumElectronDrvMu ( const DistSparseMatrix< Real > &  SMat)

CalculateNumElectronDrvMu computes the derivative of the number of electrons with respect to the chemical potential.

Parameters
[in]SMatoverlap matrix.
Returns
The derivative of the number of electrons Tr[f'(H-) S]
Real PEXSI::PPEXSIData::CalculateNumElectronDrvT ( const DistSparseMatrix< Real > &  SMat)

CalculateNumElectronDrvT computes the derivative of the number of electrons with respect to the temperature (1/beta, in atomic unit).

Parameters
[in]SMatoverlap matrix.
Returns
The derivative of the number of electrons Tr[f'(H-) S]
Real PEXSI::PPEXSIData::CalculateTotalEnergy ( const DistSparseMatrix< Real > &  HMat)

CalculateTotalEnergy computes the total energy (band energy part only).

Parameters
[in]HMatHamilotian matrix.
Returns
The total energy Tr[ H ].
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.

void PEXSI::PPEXSIData::EstimateSpectralRadius ( Int  method,
const DistSparseMatrix< Real > &  HMat,
const DistSparseMatrix< Real > &  SMat,
const NumVec< Real > &  v0,
Real  tol,
Int  maxNumIter,
Int &  numIter,
Real &  sigma 
)

Estimate the spectral radius of the matrix stencil (H,S).

The current method uses the locally optimal conjugate gradient (LOCG) method.

Note
Rigorously speaking we are not computing the spectral radius, but the largest generalized eigenvalue. For standard Kohn-Sham problem with pseudopotential this should be OK. But this should be changed later, at least as a safe guard, but also may have practical relavance in the case of e.g. all electron calculation.
Parameters
[in]methodMethod for estimating the spectral radius
  • = 0 : Locally optimal conjugate gradient (LOCG) method (default). This works for both S being an identity matrix or not. This is the only implemented one currently.
[in]HMatHamiltonian matrix saved in distributed compressed sparse column format. See DistSparseMatrix.
[in]SMatOverlap matrix saved in distributed compressed sparse column format. See DistSparseMatrix. S can be an identity matrix.
[in]v0Initial starting vector.

If v0.m() == 0, then a random vector is used as the initial starting vector.

Parameters
[in]tolRelative tolerance for estimating sigma.
[in]maxNumIterMaximum number of iterations.
[out]numIterThe number of iterations.
[out]sigmaThe estimated spectral radius.
Real PEXSI::PPEXSIData::EstimateZeroTemperatureChemicalPotential ( Real  temperature,
Real  mu,
const DistSparseMatrix< Real > &  SMat 
)

Estimate the chemical potential at zero temperature using quadratic approximation.

Parameters
[in]temperatureTemperature (in the unit of K), usually high (~3000K)
[in]muComputed chemical potential at high temperature
[in]SMatOverlap matrix
Returns
Estimated chemical potential at zero temperature.
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::Solve ( Int  numPole,
Real  temperature,
Real  numElectronExact,
Real  gap,
Real  deltaE,
Real &  mu,
Real &  muMin,
Real &  muMax,
const DistSparseMatrix< Real > &  HMat,
const DistSparseMatrix< Real > &  SMat,
Int  muMaxIter,
Real  numElectronTolerance,
std::string  ColPerm,
Int  numProcSymbFact,
bool  isFreeEnergyDensityMatrix,
bool  isEnergyDensityMatrix,
bool  isDerivativeTMatrix,
std::vector< Real > &  muList,
std::vector< Real > &  numElectronList,
std::vector< Real > &  numElectronDrvMuList,
bool &  isConverged 
)

Solve is the main subroutine for PPEXSI.

Compute 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 the unit of K)
[in]numElectronExactExact number of electrons
[in]gapBand gap (in the unit of au)
[in]deltaEUpperbound of the spectrum width
[in,out]muInput: Initial guess of chemical potential (in the unit of au). Output: The final update of the chemical potential, for which the number of electrons has NOT been computed. Note: The output mu is NOT the same as last element in muList. muList[end] corresponds to numElectronList[end].
[in,out]muMinInput: Initial guess for the lower bound of the chemical potential. Output: Lower bound of the chemical potential.
[in,out]muMaxInput: Initial guess for the upper bound of the chemical potential. Output: Upper bound of the chemical potential.
[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.
[in]muMaxIterMaximum iteration number for chemical potential
[in]numElectronToleranceStopping criterion for the mu iteration
[in]ColPermPermutation method used for SuperLU_DIST
[in]numProcSymbFactNumber of processors used for parallel symbolic factorization and PARMETIS/PT-SCOTCH.
[in]isFreeEnergyDensityMatrixWhether to compute the Helmholtz free energy matrix
[in]isEnergyDensityMatrixWhether to compute the energy density matrix for force
[in]isDerivativeTMatrixWhether to compute the derivative of the single particle density matrix with respect to the temperature.
[out]muListConvergence history of the chemical potential
[out]numElectronListConvergence history of the number of electrons
[out]numElectronDrvMuListConvergence history of the derivative of electrons with respect to the chemical potential.
[out]isConvergedWhether PEXSI has converged.
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::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.

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