Interface subroutines of PEXSI that can be called by C. More...
#include <mpi.h>
#include <stdint.h>
Go to the source code of this file.
Classes | |
struct | PPEXSIOptions |
Structure for the input parameters in DFT calculations. More... | |
Typedefs | |
typedef intptr_t | PPEXSIPlan |
A handle for holding the internal PEXSI data structure. More... | |
Functions | |
void | ReadDistSparseMatrixFormattedHeadInterface (char *filename, int *size, int *nnz, int *nnzLocal, int *numColLocal, MPI_Comm comm) |
Read the sizes of a DistSparseMatrix in formatted form (txt) for allocating memory in C. More... | |
void | ReadDistSparseMatrixFormattedInterface (char *filename, int size, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *nzvalLocal, MPI_Comm comm) |
Reading the data of a formatted DistSparseMatrix. More... | |
void | ReadDistSparseMatrixHeadInterface (char *filename, int *size, int *nnz, int *nnzLocal, int *numColLocal, MPI_Comm comm) |
Read the sizes of a DistSparseMatrix in unformatted form (csc) for allocating memory in C. More... | |
void | ParaReadDistSparseMatrixInterface (char *filename, int size, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *nzvalLocal, MPI_Comm comm) |
Actual reading the data of a DistSparseMatrix using MPI-IO, assuming that the arrays have been allocated outside this subroutine. More... | |
void | PPEXSISetDefaultOptions (PPEXSIOptions *options) |
Set the default options for DFT driver. More... | |
PPEXSIPlan | PPEXSIPlanInitialize (MPI_Comm comm, int numProcRow, int numProcCol, int outputFileIndex, int *info) |
Initialize the PEXSI plan. More... | |
void | PPEXSILoadRealHSMatrix (PPEXSIPlan plan, PPEXSIOptions options, int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *HnzvalLocal, int isSIdentity, double *SnzvalLocal, int *info) |
Load the real H and S matrices into the PEXSI internal data structure. H and S can be either real symmetric or real unsymmetric. More... | |
void | PPEXSILoadComplexHSMatrix (PPEXSIPlan plan, PPEXSIOptions options, int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *HnzvalLocal, int isSIdentity, double *SnzvalLocal, int *info) |
Load the complex H and S matrices into the PEXSI internal data structure. H and S can be either complex symmetric or complex unsymmetric. More... | |
void | PPEXSISymbolicFactorizeRealSymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, int *info) |
Separately perform symbolic factorization to prepare factorization and selected inversion for real arithmetic matrices. More... | |
void | PPEXSISymbolicFactorizeComplexSymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, int *info) |
Separately perform symbolic factorization to prepare factorization and selected inversion for complex arithmetic matrices. More... | |
void | PPEXSIInertiaCountRealMatrix (PPEXSIPlan plan, PPEXSIOptions options, int numShift, double *shiftList, double *inertiaList, int *info) |
Directly compute the negative inertia at a set of shifts for real matrices. More... | |
void | PPEXSIInertiaCountComplexMatrix (PPEXSIPlan plan, PPEXSIOptions options, int numShift, double *shiftList, double *inertiaList, int *info) |
Directly compute the negative inertia at a set of shifts for complex matrices. More... | |
void | PPEXSISymbolicFactorizeRealUnsymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, double *AnzvalLocal, int *info) |
Separately perform symbolic factorization to prepare factorization and selected inversion for real arithmetic matrices. More... | |
void | PPEXSISymbolicFactorizeComplexUnsymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, double *AnzvalLocal, int *info) |
Separately perform symbolic factorization to prepare factorization and selected inversion for complex arithmetic matrices. More... | |
void | PPEXSICalculateFermiOperatorReal (PPEXSIPlan plan, PPEXSIOptions options, double mu, double numElectronExact, double *numElectronPEXSI, double *numElectronDrvMuPEXSI, int *info) |
Directly compute the negative inertia at a set of shifts. More... | |
void | PPEXSICalculateFermiOperatorComplex (PPEXSIPlan plan, PPEXSIOptions options, double mu, double numElectronExact, double *numElectronPEXSI, double *numElectronDrvMuPEXSI, int *info) |
Compute the density matrices and number of electrons for a given chemical potential. Here H and S are complex Hermitian matrices (even though S is real in electronic structure calculation) More... | |
void | PPEXSISelInvRealSymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, double *AnzvalLocal, double *AinvnzvalLocal, int *info) |
Simplified driver interface for computing the selected elements of a real symmetric matrix. More... | |
void | PPEXSISelInvComplexSymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, double *AnzvalLocal, double *AinvnzvalLocal, int *info) |
Simplified driver interface for computing the selected elements of a complex symmetric matrix. More... | |
void | PPEXSISelInvRealUnsymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, double *AnzvalLocal, double *AinvnzvalLocal, int *info) |
Simplified driver interface for computing the selected elements of a real unsymmetric matrix. More... | |
void | PPEXSISelInvComplexUnsymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, double *AnzvalLocal, double *AinvnzvalLocal, int *info) |
Simplified driver interface for computing the selected elements of a complex unsymmetric matrix. More... | |
void | PPEXSIDFTDriver (PPEXSIPlan plan, PPEXSIOptions options, double numElectronExact, double *muPEXSI, double *numElectronPEXSI, double *muMinInertia, double *muMaxInertia, int *numTotalInertiaIter, int *numTotalPEXSIIter, int *info) |
Simplified driver for solving Kohn-Sham DFT. More... | |
void | PPEXSIDFTDriver2 (PPEXSIPlan plan, PPEXSIOptions options, double numElectronExact, double *muPEXSI, double *numElectronPEXSI, double *muMinInertia, double *muMaxInertia, int *numTotalInertiaIter, int *info) |
Simplified driver for solving Kohn-Sham DFT. Version 2. This is still in test phase. More... | |
void | PPEXSIRetrieveRealDFTMatrix (PPEXSIPlan plan, double *DMnzvalLocal, double *EDMnzvalLocal, double *FDMnzvalLocal, double *totalEnergyH, double *totalEnergyS, double *totalFreeEnergy, int *info) |
Retrieve the output matrices after running PPEXSIDFTDriver for real input matrices. More... | |
void | PPEXSIRetrieveComplexDFTMatrix (PPEXSIPlan plan, double *DMnzvalLocal, double *EDMnzvalLocal, double *FDMnzvalLocal, double *totalEnergyH, double *totalEnergyS, double *totalFreeEnergy, int *info) |
Retrieve the output matrices after running PPEXSIDFTDriver for complex input matrices. More... | |
void | PPEXSIPlanFinalize (PPEXSIPlan plan, int *info) |
Release the memory used by PEXSI. More... | |
void | PPEXSIGetPoleDM (double *zshift, double *zweight, int Npole, double temp, double gap, double deltaE, double mu) |
Pole expansion for density matrix (DM). More... | |
void | PPEXSIGetPoleEDM (double *zshift, double *zweight, int Npole, double temp, double gap, double deltaE, double mu) |
Pole expansion for energy density matrix (EDM). More... | |
void | PPEXSIGetPoleFDM (double *zshift, double *zweight, int Npole, double temp, double gap, double deltaE, double mu) |
Pole expansion for free energy density matrix (FDM). More... | |
Interface subroutines of PEXSI that can be called by C.
typedef intptr_t PPEXSIPlan |
A handle for holding the internal PEXSI data structure.
INTEGER*8
data structure, or the INTEGER(C_INTPTR_T)
structure if ISO_C_BINDING is used. void ParaReadDistSparseMatrixInterface | ( | char * | filename, |
int | size, | ||
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | nzvalLocal, | ||
MPI_Comm | comm | ||
) |
Actual reading the data of a DistSparseMatrix using MPI-IO, assuming that the arrays have been allocated outside this subroutine.
This routine can be much faster than reading a DistSparseMatrix sequentially, especially compared to the version using formatted input ReadDistSparseMatrixFormattedInterface.
[in] | filename | (global) Filename for the input matrix. |
[in] | size | (global) Number of rows and columns of the matrix. |
[in] | nnz | (global) Total number of nonzeros. |
[in] | nnzLocal | (local) Number of local nonzeros. |
[in] | numColLocal | (local) Number of local columns. |
[out] | colptrLocal | (local) Dimension: numColLocal+1. Local column pointer in CSC format. |
[out] | rowindLocal | (local) Dimension: nnzLocal. Local row index pointer in CSC format. |
[out] | nzvalLocal | (local) Dimension: nnzLocal. Local nonzero values in CSC format. |
[in] | comm | (global) MPI communicator. |
void PPEXSICalculateFermiOperatorComplex | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double | mu, | ||
double | numElectronExact, | ||
double * | numElectronPEXSI, | ||
double * | numElectronDrvMuPEXSI, | ||
int * | info | ||
) |
Compute the density matrices and number of electrons for a given chemical potential. Here H and S are complex Hermitian matrices (even though S is real in electronic structure calculation)
This can be used as an "expert" interface for solving KSDFT with user-implemented heuristics strategies.
Input/Output
The input parameter options are controlled through the structure PPEXSIOptions. The default value can be obtained through PPEXSISetDefaultOptions.
The input H and S matrices should be given by loading functions (currently it is PPEXSILoadRealSymmetricHSMatrix). The output matrices should be obtained from retrieving functions (currently it is PPEXSIRetrieveRealSymmetricDFTMatrix).
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | numElectronExact | (global) Exact number of electrons, i.e. \(N_e(\mu_{\mathrm{exact}})\). |
[in] | options | (global) Other input parameters for the DFT driver. |
[out] | DMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of density matrix in CSC format. |
[out] | EDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of energy density matrix in CSC format. |
[out] | FDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of free energy density matrix in CSC format. |
[out] | muPEXSI | (global) Chemical potential after the last iteration.
|
[out] | numElectronPEXSI | (global) Number of electrons evaluated at the last step. Note In the case that convergence is not reached within maxPEXSIIter steps, and numElectron does not correspond to the number of electrons evaluated at muPEXSI. |
[out] | muMinInertia | (global) Lower bound for mu after the last inertia count procedure. |
[out] | muMaxInertia | (global) Upper bound for mu after the last inertia count procedure. |
[out] | numTotalInertiaIter | (global) Number of total inertia counting procedure. |
[out] | numTotalPEXSIIter | (global) Number of total PEXSI evaluation procedure. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSICalculateFermiOperatorReal | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double | mu, | ||
double | numElectronExact, | ||
double * | numElectronPEXSI, | ||
double * | numElectronDrvMuPEXSI, | ||
int * | info | ||
) |
Directly compute the negative inertia at a set of shifts.
This can be used as an "expert" interface for solving KSDFT with user-implemented heuristics strategies.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | numShift | (global) Number of shifts. |
[in] | shiftList | (global) The list of shifts. Size: numShift |
[out] | inertiaList | (global) The list of inertia counts (in double precision but are of integer values). Size: numShift |
[out] | info | (local) whether the current processor returns the correct information.
|
This can be used as an "expert" interface for solving KSDFT with user-implemented heuristics strategies.
Input/Output
The input parameter options are controlled through the structure PPEXSIOptions. The default value can be obtained through PPEXSISetDefaultOptions.
The input H and S matrices should be given by loading functions (currently it is PPEXSILoadRealSymmetricHSMatrix). The output matrices should be obtained from retrieving functions (currently it is PPEXSIRetrieveRealSymmetricDFTMatrix).
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | numElectronExact | (global) Exact number of electrons, i.e. \(N_e(\mu_{\mathrm{exact}})\). |
[in] | options | (global) Other input parameters for the DFT driver. |
[out] | DMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of density matrix in CSC format. |
[out] | EDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of energy density matrix in CSC format. |
[out] | FDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of free energy density matrix in CSC format. |
[out] | muPEXSI | (global) Chemical potential after the last iteration.
|
[out] | numElectronPEXSI | (global) Number of electrons evaluated at the last step. Note In the case that convergence is not reached within maxPEXSIIter steps, and numElectron does not correspond to the number of electrons evaluated at muPEXSI. |
[out] | muMinInertia | (global) Lower bound for mu after the last inertia count procedure. |
[out] | muMaxInertia | (global) Upper bound for mu after the last inertia count procedure. |
[out] | numTotalInertiaIter | (global) Number of total inertia counting procedure. |
[out] | numTotalPEXSIIter | (global) Number of total PEXSI evaluation procedure. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSIDFTDriver | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double | numElectronExact, | ||
double * | muPEXSI, | ||
double * | numElectronPEXSI, | ||
double * | muMinInertia, | ||
double * | muMaxInertia, | ||
int * | numTotalInertiaIter, | ||
int * | numTotalPEXSIIter, | ||
int * | info | ||
) |
Simplified driver for solving Kohn-Sham DFT.
This function contains both the inertia counting step for estimating the chemical potential, and the Newton's iteration for updating the chemical potential. Heuristics are built into this routine. Expert users and developers can modify this routine to obtain better heuristics. The implementation of this function contains all the heuristics for a DFT solver.
The input parameter options are controlled through the structure PPEXSIOptions. The default value can be obtained through PPEXSISetDefaultOptions.
Basic strategy of the heuristics
Complex Hermitian case
This file should work for both real symmetric and complex Hermitian H and S matrices. However, the complex Hermitian case require asymmetric PSelInv which will be in the future work.
Input/Output
The input H and S matrices should be given by loading functions (currently it is PPEXSILoadRealSymmetricHSMatrix). The output matrices should be obtained from retrieving functions (currently it is PPEXSIRetrieveRealSymmetricDFTMatrix).
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | numElectronExact | (global) Exact number of electrons, i.e. \(N_e(\mu_{\mathrm{exact}})\). |
[out] | DMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of density matrix in CSC format. |
[out] | EDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of energy density matrix in CSC format. |
[out] | FDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of free energy density matrix in CSC format. |
[out] | muPEXSI | (global) Chemical potential after the last iteration.
|
[out] | numElectronPEXSI | (global) Number of electrons evaluated at the last step. Note In the case that convergence is not reached within maxPEXSIIter steps, and numElectron does not correspond to the number of electrons evaluated at muPEXSI. |
[out] | muMinInertia | (global) Lower bound for mu after the last inertia count procedure. |
[out] | muMaxInertia | (global) Upper bound for mu after the last inertia count procedure. |
[out] | numTotalInertiaIter | (global) Number of total inertia counting procedure. |
[out] | numTotalPEXSIIter | (global) Number of total PEXSI evaluation procedure. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSIDFTDriver2 | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double | numElectronExact, | ||
double * | muPEXSI, | ||
double * | numElectronPEXSI, | ||
double * | muMinInertia, | ||
double * | muMaxInertia, | ||
int * | numTotalInertiaIter, | ||
int * | info | ||
) |
Simplified driver for solving Kohn-Sham DFT. Version 2. This is still in test phase.
This function contains both the inertia counting step for estimating the chemical potential, and strategy for updating the chemical potential WITHOUT Newton's iteration. Heuristics are built into this routine. Expert users and developers can modify this routine to obtain better heuristics. The implementation of this function contains all the heuristics for a DFT solver.
The input parameter options are controlled through the structure PPEXSIOptions. The default value can be obtained through PPEXSISetDefaultOptions.
Basic strategy of the heuristics
Input/Output
The input H and S matrices should be given by loading functions (currently it is PPEXSILoadRealSymmetricHSMatrix). The output matrices should be obtained from retrieving functions (currently it is PPEXSIRetrieveRealSymmetricDFTMatrix).
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | numElectronExact | (global) Exact number of electrons, i.e. \(N_e(\mu_{\mathrm{exact}})\). |
[out] | DMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of density matrix in CSC format. |
[out] | EDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of energy density matrix in CSC format. |
[out] | FDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of free energy density matrix in CSC format. |
[out] | muPEXSI | (global) Chemical potential after the last iteration.
|
[out] | numElectronPEXSI | (global) Number of electrons evaluated at the last step. Note In the case that convergence is not reached within maxPEXSIIter steps, and numElectron does not correspond to the number of electrons evaluated at muPEXSI. |
[out] | muMinInertia | (global) Lower bound for mu after the last inertia count procedure. |
[out] | muMaxInertia | (global) Upper bound for mu after the last inertia count procedure. |
[out] | numTotalInertiaIter | (global) Number of total inertia counting procedure. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSIGetPoleDM | ( | double * | zshift, |
double * | zweight, | ||
int | Npole, | ||
double | temp, | ||
double | gap, | ||
double | deltaE, | ||
double | mu | ||
) |
Pole expansion for density matrix (DM).
[out] | zshift | (global) Quadrature point (shift). |
[out] | zweight | (global) Quadrature weight. |
[in] | Npole | (global) Number of poles. |
[in] | temp | (global) Temperature (unit: au). |
[in] | gap | (global) Energy gap (unit: au). |
[in] | deltaE | (global) Spectral radius (unit: au). |
[in] | mu | (global) Chemical potential (unit: au). |
void PPEXSIGetPoleEDM | ( | double * | zshift, |
double * | zweight, | ||
int | Npole, | ||
double | temp, | ||
double | gap, | ||
double | deltaE, | ||
double | mu | ||
) |
Pole expansion for energy density matrix (EDM).
[out] | zshift | (global) Quadrature point (shift). |
[out] | zweight | (global) Quadrature weight. |
[in] | Npole | (global) Number of poles. |
[in] | temp | (global) Temperature (unit: au). |
[in] | gap | (global) Energy gap (unit: au). |
[in] | deltaE | (global) Spectral radius (unit: au). |
[in] | mu | (global) Chemical potential (unit: au). |
void PPEXSIGetPoleFDM | ( | double * | zshift, |
double * | zweight, | ||
int | Npole, | ||
double | temp, | ||
double | gap, | ||
double | deltaE, | ||
double | mu | ||
) |
Pole expansion for free energy density matrix (FDM).
[out] | zshift | (global) Quadrature point (shift). |
[out] | zweight | (global) Quadrature weight. |
[in] | Npole | (global) Number of poles. |
[in] | temp | (global) Temperature (unit: au). |
[in] | gap | (global) Energy gap (unit: au). |
[in] | deltaE | (global) Spectral radius (unit: au). |
[in] | mu | (global) Chemical potential (unit: au). |
void PPEXSIInertiaCountComplexMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
int | numShift, | ||
double * | shiftList, | ||
double * | inertiaList, | ||
int * | info | ||
) |
Directly compute the negative inertia at a set of shifts for complex matrices.
This can be used as an "expert" interface for solving KSDFT with user-implemented heuristics strategies.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | numShift | (global) Number of shifts. |
[in] | shiftList | (global) The list of shifts. Size: numShift |
[out] | inertiaList | (global) The list of inertia counts (in double precision but are of integer values). Size: numShift |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSIInertiaCountRealMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
int | numShift, | ||
double * | shiftList, | ||
double * | inertiaList, | ||
int * | info | ||
) |
Directly compute the negative inertia at a set of shifts for real matrices.
This can be used as an "expert" interface for solving KSDFT with user-implemented heuristics strategies.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | numShift | (global) Number of shifts. |
[in] | shiftList | (global) The list of shifts. Size: numShift |
[out] | inertiaList | (global) The list of inertia counts (in double precision but are of integer values). Size: numShift |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSILoadComplexHSMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
int | nrows, | ||
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | HnzvalLocal, | ||
int | isSIdentity, | ||
double * | SnzvalLocal, | ||
int * | info | ||
) |
Load the complex H and S matrices into the PEXSI internal data structure. H and S can be either complex symmetric or complex unsymmetric.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | nrows | (global) Number of rows and columns of the matrix. |
[in] | nnz | (global) Total number of nonzeros of H. |
[in] | nnzLocal | (local) Number of local nonzeros of H. |
[in] | numColLocal | (local) Number of local columns for H. |
[in] | colptrLocal | (local) Dimension: numColLocal+1. Local column pointer in CSC format. |
[in] | rowindLocal | (local) Dimension: nnzLocal. Local row index pointer in CSC format. |
[in] | HnzvalLocal | (local) Dimension: 2*nnzLocal. Local nonzero values of H in CSC format. |
[in] | isSIdentity | (global) Whether S is an identity matrix. If so, the variable SnzvalLocal is omitted. |
[in] | SnzvalLocal | (local) Dimension: 2*nnzLocal. Local nonzero value of S in CSC format. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSILoadRealHSMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
int | nrows, | ||
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | HnzvalLocal, | ||
int | isSIdentity, | ||
double * | SnzvalLocal, | ||
int * | info | ||
) |
Load the real H and S matrices into the PEXSI internal data structure. H and S can be either real symmetric or real unsymmetric.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | nrows | (global) Number of rows and columns of the matrix. |
[in] | nnz | (global) Total number of nonzeros of H. |
[in] | nnzLocal | (local) Number of local nonzeros of H. |
[in] | numColLocal | (local) Number of local columns for H. |
[in] | colptrLocal | (local) Dimension: numColLocal+1. Local column pointer in CSC format. |
[in] | rowindLocal | (local) Dimension: nnzLocal. Local row index pointer in CSC format. |
[in] | HnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of H in CSC format. |
[in] | isSIdentity | (global) Whether S is an identity matrix. If so, the variable SnzvalLocal is omitted. |
[in] | SnzvalLocal | (local) Dimension: nnzLocal. Local nonzero value of S in CSC format. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSIPlanFinalize | ( | PPEXSIPlan | plan, |
int * | info | ||
) |
Release the memory used by PEXSI.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[out] | info | (local) whether the current processor returns the correct information.
|
PPEXSIPlan PPEXSIPlanInitialize | ( | MPI_Comm | comm, |
int | numProcRow, | ||
int | numProcCol, | ||
int | outputFileIndex, | ||
int * | info | ||
) |
Initialize the PEXSI plan.
In PEXSI, a matrix is generally referred to as a "pole". The factorization and selected inversion procedure for a pole is computed in parallel using numProcRow * numProcCol
processors.
When only selected inversion (PSelInv) is used, it is recommended to set the mpisize of the communicator comm
to be numProcRow * numProcCol
.
When PEXSI is used to evaluate a large number of inverse matrices such as in the electronic structure calculation, mpisize should be numPole*numProcRow * numProcCol
, where numPole
inverse matrices can be processed in parallel.
[in] | comm | (global) Communicator used for the entire PEXSI procedure. The size of this communicator should be a multiple of npPerPole = numProcRow * numProcCol, which is the total number of processors used by each individual pole. |
[in] | numProcRow | (global) Number of processors in the row communication group for each pole. |
[in] | numProcCol | (global) Number of processors in the column communication group for each pole. |
[in] | outputFileIndex | (local) The index for the PEXSI output file. For instance, if this index is 1, then the corresponding processor will output to the file logPEXSI1 . If the index is negative, then no output file is given. Note Each processor must output to a different file. By default, outputFileIndex can be set as mpirank. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSIRetrieveComplexDFTMatrix | ( | PPEXSIPlan | plan, |
double * | DMnzvalLocal, | ||
double * | EDMnzvalLocal, | ||
double * | FDMnzvalLocal, | ||
double * | totalEnergyH, | ||
double * | totalEnergyS, | ||
double * | totalFreeEnergy, | ||
int * | info | ||
) |
Retrieve the output matrices after running PPEXSIDFTDriver for complex input matrices.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[out] | DMnzvalLocal | (local) Dimension: 2*nnzLocal. Nonzero value of density matrix in CSC format. |
[out] | EDMnzvalLocal | (local) Dimension: 2*nnzLocal. Nonzero value of energy density matrix in CSC format. |
[out] | FDMnzvalLocal | (local) Dimension: 2*nnzLocal. Nonzero value of free energy density matrix in CSC format. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSIRetrieveRealDFTMatrix | ( | PPEXSIPlan | plan, |
double * | DMnzvalLocal, | ||
double * | EDMnzvalLocal, | ||
double * | FDMnzvalLocal, | ||
double * | totalEnergyH, | ||
double * | totalEnergyS, | ||
double * | totalFreeEnergy, | ||
int * | info | ||
) |
Retrieve the output matrices after running PPEXSIDFTDriver for real input matrices.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[out] | DMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of density matrix in CSC format. |
[out] | EDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of energy density matrix in CSC format. |
[out] | FDMnzvalLocal | (local) Dimension: nnzLocal. Nonzero value of free energy density matrix in CSC format. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISelInvComplexSymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double * | AnzvalLocal, | ||
double * | AinvnzvalLocal, | ||
int * | info | ||
) |
Simplified driver interface for computing the selected elements of a complex symmetric matrix.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | AnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of A in CSC format.
|
[out] | AinvnzvalLocal | (local) Dimension: 2*nnzLocal. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISelInvComplexUnsymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double * | AnzvalLocal, | ||
double * | AinvnzvalLocal, | ||
int * | info | ||
) |
Simplified driver interface for computing the selected elements of a complex unsymmetric matrix.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | AnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of A in CSC format.
|
[out] | AinvnzvalLocal | (local) Dimension: 2*nnzLocal. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISelInvRealSymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double * | AnzvalLocal, | ||
double * | AinvnzvalLocal, | ||
int * | info | ||
) |
Simplified driver interface for computing the selected elements of a real symmetric matrix.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | AnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of A in CSC format. |
[out] | AinvnzvalLocal | (local) Dimension: nnzLocal. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISelInvRealUnsymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double * | AnzvalLocal, | ||
double * | AinvnzvalLocal, | ||
int * | info | ||
) |
Simplified driver interface for computing the selected elements of a real unsymmetric matrix.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | AnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of A in CSC format. |
[out] | AinvnzvalLocal | (local) Dimension: nnzLocal. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISetDefaultOptions | ( | PPEXSIOptions * | options | ) |
Set the default options for DFT driver.
All default values assume the input unit (for H) is Rydberg.
[in] | options | (global) Pointer to the options containing input parameters for the driver. |
void PPEXSISymbolicFactorizeComplexSymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
int * | info | ||
) |
Separately perform symbolic factorization to prepare factorization and selected inversion for complex arithmetic matrices.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISymbolicFactorizeComplexUnsymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double * | AnzvalLocal, | ||
int * | info | ||
) |
Separately perform symbolic factorization to prepare factorization and selected inversion for complex arithmetic matrices.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | AnzvalLocal | non zero values required to compute row permutation |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISymbolicFactorizeRealSymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
int * | info | ||
) |
Separately perform symbolic factorization to prepare factorization and selected inversion for real arithmetic matrices.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[out] | info | (local) whether the current processor returns the correct information.
|
void PPEXSISymbolicFactorizeRealUnsymmetricMatrix | ( | PPEXSIPlan | plan, |
PPEXSIOptions | options, | ||
double * | AnzvalLocal, | ||
int * | info | ||
) |
Separately perform symbolic factorization to prepare factorization and selected inversion for real arithmetic matrices.
[in] | plan | (local) The plan holding the internal data structure for the PEXSI data structure. |
[in] | options | (global) Other input parameters for the DFT driver. |
[in] | AnzvalLocal | non zero values required to compute row permutation |
[out] | info | (local) whether the current processor returns the correct information.
|
void ReadDistSparseMatrixFormattedHeadInterface | ( | char * | filename, |
int * | size, | ||
int * | nnz, | ||
int * | nnzLocal, | ||
int * | numColLocal, | ||
MPI_Comm | comm | ||
) |
Read the sizes of a DistSparseMatrix in formatted form (txt) for allocating memory in C.
[in] | filename | (global) Filename for the input matrix. |
[out] | size | (global) Number of rows and columns of the matrix. |
[out] | nnz | (global) Total number of nonzeros. |
[out] | nnzLocal | (local) Number of local nonzeros. |
[out] | numColLocal | (local) Number of local columns. |
[in] | comm | (global) MPI communicator. |
void ReadDistSparseMatrixFormattedInterface | ( | char * | filename, |
int | size, | ||
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | nzvalLocal, | ||
MPI_Comm | comm | ||
) |
Reading the data of a formatted DistSparseMatrix.
This routine assumes that the arrays have been allocated outside this subroutine.
[in] | filename | (global) Filename for the input matrix. |
[in] | size | (global) Number of rows and columns of the matrix. |
[in] | nnz | (global) Total number of nonzeros. |
[in] | nnzLocal | (local) Number of local nonzeros. |
[in] | numColLocal | (local) Number of local columns. |
[out] | colptrLocal | (local) Dimension: numColLocal+1. Local column pointer in CSC format. |
[out] | rowindLocal | (local) Dimension: nnzLocal. Local row index pointer in CSC format. |
[out] | nzvalLocal | (local) Dimension: nnzLocal. Local nonzero values in CSC format. |
[in] | comm | (global) MPI communicator. |
void ReadDistSparseMatrixHeadInterface | ( | char * | filename, |
int * | size, | ||
int * | nnz, | ||
int * | nnzLocal, | ||
int * | numColLocal, | ||
MPI_Comm | comm | ||
) |
Read the sizes of a DistSparseMatrix in unformatted form (csc) for allocating memory in C.
[in] | filename | (global) Filename for the input matrix. |
[out] | size | (global) Number of rows and columns of the matrix. |
[out] | nnz | (global) Total number of nonzeros. |
[out] | nnzLocal | (local) Number of local nonzeros. |
[out] | numColLocal | (local) Number of local columns. |
[in] | comm | (global) MPI communicator. |