Interface subroutines of PEXSI that can be called by C. More...
#include <mpi.h>
Go to the source code of this file.
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 | PPEXSISelInvInterface (int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *HnzvalLocal, int isSIdentity, double *SnzvalLocal, double *zShift, int ordering, int npSymbFact, MPI_Comm comm, double *AinvnzvalLocal, int *info) |
Driver interface for computing the selected elements of a matrix (in complex arithmietic). More... | |
void | PPEXSIInertiaCountInterface (int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *HnzvalLocal, int isSIdentity, double *SnzvalLocal, double temperature, double numElectronExact, double muMin0, double muMax0, int numPole, int maxIter, double numElectronTolerance, int ordering, int npPerPole, int npSymbFact, MPI_Comm comm, double *muMinInertia, double *muMaxInertia, double *muLowerEdge, double *muUpperEdge, int *numIter, double *muList, double *numElectronList, int *info) |
Driver interface for estimating the chemical potential using inertia counting procedure. More... | |
void | PPEXSIRawInertiaCountInterface (int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *HnzvalLocal, int isSIdentity, double *SnzvalLocal, double muMin0, double muMax0, int numPole, int ordering, int npPerPole, int npSymbFact, MPI_Comm comm, double *muList, int *numElectronList, int *info) |
Directly compute the negative inertia at a set of shifts. More... | |
void | PPEXSISolveInterface (int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *HnzvalLocal, int isSIdentity, double *SnzvalLocal, double temperature, double numElectronExact, double mu0, double muMin0, double muMax0, double gap, double deltaE, int numPole, int maxIter, double numElectronTolerance, int ordering, int npPerPole, int npSymbFact, MPI_Comm comm, double *DMnzvalLocal, double *EDMnzvalLocal, double *FDMnzvalLocal, double *muPEXSI, double *numElectronPEXSI, double *muMinPEXSI, double *muMaxPEXSI, int *numIter, double *muList, double *numElectronList, double *numElectronDrvList, int *info) |
Main driver routine for solving Kohn-Sham density functional theory using PEXSI. More... | |
void | PPEXSILocalDOSInterface (int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *HnzvalLocal, int isSIdentity, double *SnzvalLocal, double Energy, double eta, int ordering, int npSymbFact, MPI_Comm comm, double *localDOSnzvalLocal, int *info) |
Driver interface for computing the local density of states. More... | |
void | PSelInvRealSymmetricInterface (int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *AnzvalLocal, int ordering, int npSymbFact, MPI_Comm comm, int nprow, int npcol, double *AinvnzvalLocal, int *info) |
Simplified driver interface for computing the selected elements of a real symmetric symmetric. More... | |
void | PSelInvComplexSymmetricInterface (int nrows, int nnz, int nnzLocal, int numColLocal, int *colptrLocal, int *rowindLocal, double *AnzvalLocal, int ordering, int npSymbFact, MPI_Comm comm, int nprow, int npcol, double *AinvnzvalLocal, int *info) |
Simplified driver interface for computing the selected elements of a complex symmetric symmetric. More... | |
Interface subroutines of PEXSI that can be called by C.
The FORTRAN interface are wrappers of the C interface. The header file is not needed here, and the interface routines are given directly in interface.cpp.
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 PPEXSIInertiaCountInterface | ( | int | nrows, |
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | HnzvalLocal, | ||
int | isSIdentity, | ||
double * | SnzvalLocal, | ||
double | temperature, | ||
double | numElectronExact, | ||
double | muMin0, | ||
double | muMax0, | ||
int | numPole, | ||
int | maxIter, | ||
double | numElectronTolerance, | ||
int | ordering, | ||
int | npPerPole, | ||
int | npSymbFact, | ||
MPI_Comm | comm, | ||
double * | muMinInertia, | ||
double * | muMaxInertia, | ||
double * | muLowerEdge, | ||
double * | muUpperEdge, | ||
int * | numIter, | ||
double * | muList, | ||
double * | numElectronList, | ||
int * | info | ||
) |
Driver interface for estimating the chemical potential using inertia counting procedure.
The main purpose of this routine is to obtain two set of intervals.
The first interval is [muMinInertia, muMaxInertia]. This should be a tight interval so that
\[ N_e(\mathrm{muMinInertia}) \le \mathrm{numElectronExact} \le N_e(\mathrm{muMaxInertia}). \]
This interval can be used as a bound interval for PPEXSISolveInterface.
The second interval is [muLowerEdge, muUpperEdge]. This interval is defined by
\[ N_e(\mathrm{muLowerEdge}) = \mathrm{numElectronExact} - eps, \]
\[ N_e(\mathrm{muUpppeEdge}) = \mathrm{numElectronExact} + eps. \]
eps is a small number currently fixed to be 0.1. For metals ( \(N'_e(\mu_{\mathrm{exact}}) > 0\)), it should be expected that muLowerEdge is very close to muUpperEdge. For insulators, ( \(N'_e(\mu_{\mathrm{exact}}) \approx 0\)), we have muUpperEdge > muLowerEdge, with the difference characterizing the band gap.
For either metal or insulator,
\[ \frac{1}{2} ( \mathrm{muLowerEdge} + \mathrm{muUpperEdge} ) \]
is a good guess for the chemical potential to be provided to PPEXSISolveInterface.
These intervals are refined for in an iterative procedure, until
\[ N_e(\mathrm{muMaxInertia}) - N_e(\mathrm{muMinInertia}) < numElectronTolerance. \]
The counting of the number of electrons (eigenvalues) is obtained by counting the number of negative eigenvalues for the matrix \(H-\mu_l S\) for a series of \(\{\mu_l\}\). The algorithm makes use of Sylvester's law of inertia.
Finite temperature effect is also taken into account when computing the number of electrons.
Note
[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. |
[in] | temperature | (global) Temperature, in the same unit as H |
[in] | numElectronExact | (global) Exact number of electrons, i.e. \(N_e(\mu_{\mathrm{exact}})\). |
[in] | muMin0 | (global) Initial guess of lower bound for mu |
[in] | muMax0 | (global) Initial guess of upper bound for mu |
[in] | numPole | (global) Number of shifts (chemical potentials) evaluated at each iteration. |
[in] | maxIter | (global) Maximum number of iterations. |
[in] | numElectronTolerance | (global) Stopping criterion of the iterations. |
[in] | ordering | (global) Ordering strategy for factorization and selected inversion.
|
[in] | npSymbFact | (global) Number of processors for PARMETIS/PT-SCOTCH. Only used if the ordering = 0. |
[in] | comm | (global) MPI communicator. |
[out] | muMinInertia | (global) Lower bound for mu after the inertia count. |
[out] | muMaxInertia | (global) Upper bound for mu after the inertia count. |
[out] | muLowerEdge | (global) muLowerEdge satisfies N_e(muLowerEdge) = numElectronExact - eps. |
[out] | muUpperEdge | (global) muUpperEdge satisfies N_e(muUpperEdge) = numElectronExact + eps. |
[out] | numIter | (global) Number of iterations used before reaching the stopping criterion. |
[out] | muList | (global) Dimension: numIter. The list of shifts after the final iteration. |
[out] | numElectronList | (global) Dimension: numIter. The number of electrons at the each shift. Finite temperature effect is taken into account. For obtaining the number of electrons without the finite temperature effect, see PPEXSIRawInertiaCountInterface. |
[out] | info |
|
void PPEXSILocalDOSInterface | ( | int | nrows, |
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | HnzvalLocal, | ||
int | isSIdentity, | ||
double * | SnzvalLocal, | ||
double | Energy, | ||
double | eta, | ||
int | ordering, | ||
int | npSymbFact, | ||
MPI_Comm | comm, | ||
double * | localDOSnzvalLocal, | ||
int * | info | ||
) |
Driver interface for computing the local density of states.
The local density of states (LDOS) can be evaluated as
\[ n(r;E) = \lim_{\eta->0+} \mathrm{Im} \frac{1}{\pi} \langle r\vert (H-(E+i \eta)I)^{-1} \vert r \rangle. \]
For nonorthogonal basis functions, this routine returns of the matrix
\[ \frac{1}{\pi} \mathrm{Im} ( H - (E+i \eta) S )^{-1}, \]
and the LDOS in the real space can be constructed in the same way that the density is constructed.
Note
[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. |
[in] | Energy | (global) Real part of the shift. |
[in] | eta | (global) Imaginary part of the shift, or the "broadening" parameter \(\eta\). |
[in] | ordering | (global) Ordering strategy for factorization and selected inversion.
|
[in] | npSymbFact | (global) Number of processors for PARMETIS/PT-SCOTCH. Only used if the ordering = 0. |
[in] | comm | (global) MPI communicator. |
[out] | localDOSnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of the selected elements of \(\frac{1}{\pi} \mathrm{Im} ( H - (E+i \eta) S )^{-1}\). |
[out] | info |
|
void PPEXSIRawInertiaCountInterface | ( | int | nrows, |
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | HnzvalLocal, | ||
int | isSIdentity, | ||
double * | SnzvalLocal, | ||
double | muMin0, | ||
double | muMax0, | ||
int | numPole, | ||
int | ordering, | ||
int | npPerPole, | ||
int | npSymbFact, | ||
MPI_Comm | comm, | ||
double * | muList, | ||
int * | numElectronList, | ||
int * | info | ||
) |
Directly compute the negative inertia at a set of shifts.
This is a simplified version of PPEXSIInertiaCountInterface, without the including the finite temperature effect, or iterative refinement of the chemical potential. This routine can be used for evaluating density of states after PPEXSISolveInterface has reached convergence.
See PPEXSIInertiaCountInterface for explanation of parameters.
void PPEXSISelInvInterface | ( | int | nrows, |
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | HnzvalLocal, | ||
int | isSIdentity, | ||
double * | SnzvalLocal, | ||
double * | zShift, | ||
int | ordering, | ||
int | npSymbFact, | ||
MPI_Comm | comm, | ||
double * | AinvnzvalLocal, | ||
int * | info | ||
) |
Driver interface for computing the selected elements of a matrix (in complex arithmietic).
Evaluate the selected elements of \((H - z S)^{-1}\) for a shift \(z\in\mathbb{C}\).
Note
[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. |
[in] | zShift | (global) Dimension: 2. Use 2 real numbers to denote a complex shift.
|
[in] | ordering | (global) Ordering strategy for factorization and selected inversion.
|
[in] | npSymbFact | (global) Number of processors for PARMETIS/PT-SCOTCH. Only used if the ordering = 0. |
[in] | comm | (global) MPI communicator. |
[out] | AinvnzvalLocal | (local) Dimension: 2*nnzLocal. Local nonzero values of the selected elements of \((H - z S)^{-1}\).
|
[out] | info |
|
void PPEXSISolveInterface | ( | int | nrows, |
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | HnzvalLocal, | ||
int | isSIdentity, | ||
double * | SnzvalLocal, | ||
double | temperature, | ||
double | numElectronExact, | ||
double | mu0, | ||
double | muMin0, | ||
double | muMax0, | ||
double | gap, | ||
double | deltaE, | ||
int | numPole, | ||
int | maxIter, | ||
double | numElectronTolerance, | ||
int | ordering, | ||
int | npPerPole, | ||
int | npSymbFact, | ||
MPI_Comm | comm, | ||
double * | DMnzvalLocal, | ||
double * | EDMnzvalLocal, | ||
double * | FDMnzvalLocal, | ||
double * | muPEXSI, | ||
double * | numElectronPEXSI, | ||
double * | muMinPEXSI, | ||
double * | muMaxPEXSI, | ||
int * | numIter, | ||
double * | muList, | ||
double * | numElectronList, | ||
double * | numElectronDrvList, | ||
int * | info | ||
) |
Main driver routine for solving Kohn-Sham density functional theory using PEXSI.
Given a good initial guess of the chemical potential and the bound for the chemical potential, this routine uses the Newton method to find the converged chemical potential for a given number of electrons at certain temperature.
After convergence, this routine also returns the density matrix for evaluating the electron density. It also returns the free energy density matrix for evaluating the (Helmholtz or Mermin) free energy, and the energy density matrix for evaluating the atomic force.
Note
[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. |
[in] | temperature | (global) Temperature, in the same unit as H |
[in] | numElectronExact | (global) Exact number of electrons, i.e. \(N_e(\mu_{\mathrm{exact}})\). |
[in] | mu0 | (global) Initial guess for mu. |
[in] | muMin0 | (global) Initial guess of lower bound for mu. |
[in] | muMax0 | (global) Initial guess of upper bound for mu. |
[in] | gap | (global) A lower bound for the band gap. Currently it is recommended to set this to 0, which works for all systems. |
[in] | deltaE | (global) An upper bound for the spectral radius of \(S^{-1} H\). |
[in] | numPole | (global) Number of poles. |
[in] | maxIter | (global) Maximum number of iterations. |
[in] | numElectronTolerance | (global) Stopping criterion of the iterations. |
[in] | ordering | (global) Ordering strategy for factorization and selected inversion.
|
[in] | npPerPole | (global) Number of processors per pole (ppp). |
[in] | npSymbFact | (global) Number of processors for PARMETIS/PT-SCOTCH. Only used if the ordering = 0. |
[in] | comm | (global) MPI communicator. |
[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 muPEXSI after the last iteration. |
[out] | muMinPEXSI | (global) Lower bound for mu after the last iteration. |
[out] | muMaxPEXSI | (global) Upper bound for mu after the last iteration. |
[out] | numIter | (global) Number of iterations used before reaching the stopping criterion. |
[out] | muList | (global) Dimension: numIter. The history of chemical potential. |
[out] | numElectronList | (global) Dimension: numIter. The history of the number of electrons evaluated at each mu in muList. |
[out] | numElectronDrvList | (global) Dimension: numIter. The history of the derivative of the number of electrons with respect to mu, evaluated at each mu in muList. This is used in the Newton iteration. |
[out] | info |
|
void PSelInvComplexSymmetricInterface | ( | int | nrows, |
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | AnzvalLocal, | ||
int | ordering, | ||
int | npSymbFact, | ||
MPI_Comm | comm, | ||
int | nprow, | ||
int | npcol, | ||
double * | AinvnzvalLocal, | ||
int * | info | ||
) |
Simplified driver interface for computing the selected elements of a complex symmetric symmetric.
[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] | AnzvalLocal | (local) Dimension: 2*nnzLocal. Local nonzero values of A in CSC format.
|
[in] | ordering | (global) Ordering strategy for factorization and selected inversion.
|
[in] | npSymbFact | (global) Number of processors for PARMETIS/PT-SCOTCH. Only used if the ordering = 0. |
[in] | comm | (global) MPI communicator. NOTE: mpisize should be equal to nprow * npcol |
[in] | nprow | (global) Number of row processors. |
[in] | npcol | (global) Number of column processors. |
[out] | AinvnzvalLocal | (local) Dimension: 2*nnzLocal. Local nonzero values of the selected elements of \(A^{-1}\).
|
[out] | info |
|
void PSelInvRealSymmetricInterface | ( | int | nrows, |
int | nnz, | ||
int | nnzLocal, | ||
int | numColLocal, | ||
int * | colptrLocal, | ||
int * | rowindLocal, | ||
double * | AnzvalLocal, | ||
int | ordering, | ||
int | npSymbFact, | ||
MPI_Comm | comm, | ||
int | nprow, | ||
int | npcol, | ||
double * | AinvnzvalLocal, | ||
int * | info | ||
) |
Simplified driver interface for computing the selected elements of a real symmetric symmetric.
[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] | AnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of A in CSC format. |
[in] | ordering | (global) Ordering strategy for factorization and selected inversion.
|
[in] | npSymbFact | (global) Number of processors for PARMETIS/PT-SCOTCH. Only used if the ordering = 0. |
[in] | comm | (global) MPI communicator. NOTE: mpisize should be equal to nprow * npcol |
[in] | nprow | (global) Number of row processors. |
[in] | npcol | (global) Number of column processors. |
[out] | AinvnzvalLocal | (local) Dimension: nnzLocal. Local nonzero values of the selected elements of \(A^{-1}\). |
[out] | info |
|
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. |