PEXSI
 All Classes Namespaces Files Functions Variables Typedefs Pages
Macros | Functions
interface.cpp File Reference

Interface subroutines of PPEXSI that can be called by both C and FORTRAN. More...

#include "c_pexsi_interface.h"
#include "ppexsi.hpp"
#include "blas.hpp"

Macros

#define iC(fun)   { int ierr=fun; if(ierr!=0) exit(1); }
 
#define iA(expr)   { if((expr)==0) { std::cerr<<"wrong "<<__LINE__<<" in " <<__FILE__<<std::endl; std::cerr.flush(); exit(1); } }
 

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 PPEXSILoadRealSymmetricHSMatrix (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 symmetric H and S matrices into the PEXSI internal data structure. 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 PPEXSIInertiaCountRealSymmetricMatrix (PPEXSIPlan plan, PPEXSIOptions options, int numShift, double *shiftList, double *inertiaList, int *info)
 Directly compute the negative inertia at a set of shifts. More...
 
void PPEXSICalculateFermiOperatorReal (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. 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 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 PPEXSIRetrieveRealSymmetricDFTMatrix (PPEXSIPlan plan, double *DMnzvalLocal, double *EDMnzvalLocal, double *FDMnzvalLocal, double *totalEnergyH, double *totalEnergyS, double *totalFreeEnergy, int *info)
 Retrieve the output matrices after running PPEXSIDFTDriver. More...
 
void PPEXSIPlanFinalize (PPEXSIPlan plan, int *info)
 Release the memory used by PEXSI. More...
 
MPI_Comm f2c_comm (MPI_Fint *Fcomm)
 Internal subroutine to convert FORTRAN communicator to C.
 
void f_read_distsparsematrix_formatted_head (char *filename, int *size, int *nnz, int *nnzLocal, int *numColLocal, MPI_Fint *Fcomm)
 FORTRAN interface for ReadDistSparseMatrixFormattedHeadInterface.
 
void f_read_distsparsematrix_formatted (char *filename, int size, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *nzvalLocal, MPI_Fint *Fcomm)
 FORTRAN interface for ReadDistSparseMatrixFormattedInterface.
 
void f_read_distsparsematrix_head (char *filename, int *size, int *nnz, int *nnzLocal, int *numColLocal, MPI_Fint *Fcomm)
 FORTRAN interface for ReadDistSparseMatrixHeadInterface.
 
void f_para_read_distsparsematrix (char *filename, int size, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *nzvalLocal, MPI_Fint *Fcomm)
 FORTRAN interface for ParaReadDistSparseMatrixInterface.
 
PPEXSIPlan f_ppexsi_plan_initialize (MPI_Fint *Fcomm, int numProcRow, int numProcCol, int outputFileIndex, int *info)
 FORTRAN interface for PPEXSIPlanInitialize.
 

Detailed Description

Interface subroutines of PPEXSI that can be called by both C and FORTRAN.

This file will eventually merge with interface.cpp.

Date
Original: 2013-02-03
Revision: 2014-01-01
Revision: 2014-03-07 Second generation interface.

Function Documentation

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.

Parameters
[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 PPEXSICalculateFermiOperatorReal ( 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.

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).

Parameters
[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.
  • In the case that convergence is reached within maxPEXSIIter steps, the value of muPEXSI is the last mu used to achieve accuracy within numElectronPEXSITolerance.
  • In the case that convergence is not reached within maxPEXSIIter steps, and the update from Newton's iteration does not exceed muPEXSISafeGuard, the value of muPEXSI is the last mu plus the update from Newton's 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.
  • = 0: successful exit.
  • > 0: unsuccessful.
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

  • If isInertiaCount == 1, then the inertia counting procedure is invoked until the chemical potential interval is refined from size (muMax0-muMin0) to muInertiaTolerance, or the estimated band gap is larger than muInertiaTolerance.
  • If the change of step in Newton's iteration is larger than muPEXSISafeGuard, the the inertia counting procedure is invoked again, starting from (muMin0, muMax0). If Newton's iteration fails again, the subroutine returns error message with info = 1.
  • The number of shifts in the inertia count is always automatically chosen to be (mpisize / npPerPole). This minimizes the wall clock time of the inertia counting procedure.

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).

Parameters
[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.
  • In the case that convergence is reached within maxPEXSIIter steps, the value of muPEXSI is the last mu used to achieve accuracy within numElectronPEXSITolerance.
  • In the case that convergence is not reached within maxPEXSIIter steps, and the update from Newton's iteration does not exceed muPEXSISafeGuard, the value of muPEXSI is the last mu plus the update from Newton's 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.
  • = 0: successful exit.
  • > 0: unsuccessful.
void PPEXSIInertiaCountRealSymmetricMatrix ( PPEXSIPlan  plan,
PPEXSIOptions  options,
int  numShift,
double *  shiftList,
double *  inertiaList,
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.

Note
Only input from the processors associated with the first pole is required. The information will be broadcast to the other processors in the communicator.
Parameters
[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.
  • = 0: successful exit.
  • > 0: unsuccessful.
void PPEXSILoadRealSymmetricHSMatrix ( 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 symmetric H and S matrices into the PEXSI internal data structure.

Note
Only input from the processors associated with the first pole is required. The information will be broadcast to the other processors in the communicator.
Parameters
[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.
  • = 0: successful exit.
  • > 0: unsuccessful.
void PPEXSIPlanFinalize ( PPEXSIPlan  plan,
int *  info 
)

Release the memory used by PEXSI.

Parameters
[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.
  • = 0: successful exit.
  • > 0: unsuccessful.
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.

Parameters
[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. 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.
  • = 0: successful exit.
  • > 0: unsuccessful.
Returns
(local) The plan holding the internal data structure for the PEXSI data structure.
void PPEXSIRetrieveRealSymmetricDFTMatrix ( PPEXSIPlan  plan,
double *  DMnzvalLocal,
double *  EDMnzvalLocal,
double *  FDMnzvalLocal,
double *  totalEnergyH,
double *  totalEnergyS,
double *  totalFreeEnergy,
int *  info 
)

Retrieve the output matrices after running PPEXSIDFTDriver.

The output matrices are of real arithmetic.

Parameters
[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.
  • = 0: successful exit.
  • > 0: unsuccessful.
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.

Note
The computation is only performed using the group of processors corresponding to the first pole.
Parameters
[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.
  • Use 2 double for one complex number. This ensures the compatibility with FORTRAN.
  • Real part: AnzvalLocal[2*k]. Imag part: AnzvalLocal[2*k+1].
[out]AinvnzvalLocal(local) Dimension: 2*nnzLocal.
[out]info(local) whether the current processor returns the correct information.
  • = 0: successful exit.
  • > 0: unsuccessful.
Note
The computation is only performed using the group of processors corresponding to the first pole.
Parameters
[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: 2*nnzLocal. Local nonzero values of A in CSC format.
  • Use 2 double for one complex number. This ensures the compatibility with FORTRAN.
  • Real part: AnzvalLocal[2*k]. Imag part: AnzvalLocal[2*k+1].
[out]AinvnzvalLocal(local) Dimension: 2*nnzLocal. Local nonzero values of the selected elements of \(A^{-1}\). The format is the same as of AnzvalLocal.
[out]info(local) whether the current processor returns the correct information.
  • = 0: successful exit.
  • > 0: unsuccessful.
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.

Note
The computation is only performed using the group of processors corresponding to the first pole.
Parameters
[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.
  • = 0: successful exit.
  • > 0: unsuccessful.
void PPEXSISetDefaultOptions ( PPEXSIOptions options)

Set the default options for DFT driver.

All default values assume the input unit (for H) is Rydberg.

Parameters
[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.

Parameters
[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.
  • = 0: successful exit.
  • > 0: unsuccessful.
void PPEXSISymbolicFactorizeRealSymmetricMatrix ( PPEXSIPlan  plan,
PPEXSIOptions  options,
int *  info 
)

Separately perform symbolic factorization to prepare factorization and selected inversion for real arithmetic matrices.

Parameters
[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.
  • = 0: successful exit.
  • > 0: unsuccessful.
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.

Parameters
[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.

Parameters
[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.

Parameters
[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.