ParGeMSLR
|
Class of parallel CSR matrices. More...
#include <parallel_csr_matrix.hpp>
Public Member Functions | |
ParallelCsrMatrixClass () | |
The constructor of ParallelCsrMatrixClass. More... | |
ParallelCsrMatrixClass (const ParallelCsrMatrixClass< T > &mat) | |
The copy constructor of ParallelCsrMatrixClass. More... | |
ParallelCsrMatrixClass (ParallelCsrMatrixClass< T > &&mat) | |
The move constructor of ParallelCsrMatrixClass. More... | |
ParallelCsrMatrixClass< T > & | operator= (const ParallelCsrMatrixClass< T > &mat) |
The = operator of ParallelCsrMatrixClass. More... | |
ParallelCsrMatrixClass< T > & | operator= (ParallelCsrMatrixClass< T > &&mat) |
The = operator of ParallelCsrMatrixClass. More... | |
virtual | ~ParallelCsrMatrixClass () |
The destructor of ParallelCsrMatrixClass. More... | |
int | Setup (int nrow_local, int ncol_local, parallel_log &parlog) |
Free the current parallel matrix, and set up new paralle matrix. More... | |
int | Setup (int nrow_local, long int nrow_start, long int nrow_global, int ncol_local, long int ncol_start, long int ncol_global, parallel_log &parlog) |
Free the current parallel matrix, and set up new paralle matrix. The default nnz is 0. More... | |
virtual int | Clear () |
Free the current matrix. More... | |
virtual int | GetDataLocation () const |
Get the data location of the matrix. More... | |
virtual int | GetNumRowsLocal () const |
Get the local number of rows of the matrix. More... | |
virtual int | GetNumColsLocal () const |
Get the local number of columns of the matrix. More... | |
long int | GetNumRowsGlobal () const |
Get the global number of rows of the matrix. More... | |
long int | GetNumColsGlobal () const |
Get the global number of columns of the matrix. More... | |
long int | GetRowStartGlobal () const |
Get the global number of the first row of the matrix. More... | |
long int | GetColStartGlobal () const |
Get the global number of the first column of the matrix. More... | |
virtual long int | GetNumNonzeros () const |
Get the number of nonzeros in this matrix. More... | |
int | GetGraphArrays (vector_long &vtxdist, vector_long &xadj, vector_long &adjncy) |
Get the vector arrays for the graph, used as inputs for ParMETIS. Diagonal entry is removed. More... | |
CsrMatrixClass< T > & | GetDiagMat () |
Get the diagonal matrix. More... | |
CsrMatrixClass< T > & | GetOffdMat () |
Get the off-diagonal matrix. More... | |
vector_long & | GetOffdMap () |
Get the off-diagonal column map vector. More... | |
bool & | IsOffdMapSorted () |
Tell if the off-diagonal map vector is sorted. More... | |
int | SetOffdMatSorted (bool is_offd_map_sorted) |
Set the off-diagonal map vector to sorted/unsorted. More... | |
int | SortOffdMap () |
Sort the off-diagonal map vector. More... | |
int | SetupVectorPtrStr (ParallelVectorClass< T > &vec) |
Update the structure of a parallel vector to have same distribution as the row distribution of this matrix. More... | |
virtual int | Eye () |
Create an indentity matrix. More... | |
virtual int | Fill (const T &v) |
Fill the matrix with constant value. More... | |
virtual int | MoveData (const int &location) |
Move the data to another memory location. More... | |
virtual int | Scale (const T &alpha) |
Scale the current csr matrix. More... | |
int | GetDiagScale (T &scale) |
Get the value leads to max diagonal no. More... | |
int | Transpose () |
Transpost the current matrix. More... | |
int | Transpose (ParallelCsrMatrixClass< T > &AT) |
Transpost the current matrix. More... | |
int | SetupMatvec () |
Setup the matvec_helper for parallel matvec. More... | |
int | SetupMatvecStart () |
Setup the matvec_helper for parallel matvec. More... | |
int | SetupMatvecOver () |
Finishing up the setup phase of matvec_helper for parallel matvec. More... | |
virtual int | MatVec (char trans, const T &alpha, const VectorClass< T > &x, const T &beta, VectorClass< T > &y) |
In place parallel csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More... | |
virtual int | MatVec (char trans, const T &alpha, const VectorClass< T > &x, const T &beta, const VectorClass< T > &y, VectorClass< T > &z) |
In place Matrix-Vector product ==> z := alpha*A*x + beta*y, or z := alpha*A'*x + beta*y. Note that z should not equal to y. More... | |
int | MatVecOffd (char trans, const T &alpha, const VectorClass< T > &x, const T &beta, VectorClass< T > &y) |
In place parallel csr Matrix-Vector product with offdiagonal block only ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More... | |
int | MatMat (const T &alpha, const ParallelCsrMatrixClass< T > &A, char transa, const ParallelCsrMatrixClass< T > &B, char transb, const T &beta) |
Dense Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C where C is this matrix. More... | |
int | SubMatrix (vector_long &rows, vector_long &cols, int location, ParallelCsrMatrixClass< T > &parcsrmat_out) |
Extract parallel submatrix from parallel CSR matrix. More... | |
int | Laplacian (int nx, int ny, int nz, int dx, int dy, int dz, T alphax, T alphay, T alphaz, T shift, parallel_log &parlog, bool rand_perturb=false) |
Generate Laplacian matrix on the host memory, 5-pt for 2D problem and 7-pt for 3D problem. More... | |
int | Helmholtz (int n, int dx, int dy, int dz, T w, parallel_log &parlog) |
Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = f. More... | |
int | LaplacianWithPartition (int nx, int ny, int nz, int dx, int dy, int dz, int d2x, int d2y, int d2z, T alphax, T alphay, T alphaz, T shift, parallel_log &parlog, bool rand_perturb=false) |
Generate Laplacian matrix on the host memory, 5-pt for 2D problem and 7-pt for 3D problem. More... | |
int | HelmholtzWithPartition (int n, int dx, int dy, int dz, int d2x, int d2y, int d2z, T w, parallel_log &parlog) |
Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = f. More... | |
int | ReadFromSingleMMFile (const char *matfile, int idxin, parallel_log &parlog) |
Read matrix on the host memory, matrix market format. More... | |
int | PlotPatternGnuPlot (const char *datafilename, int pttype=0) |
Plot the pattern of the parallel csr matrix using gnuplot. Similar to spy in the MATLAB. Function for testing purpose. More... | |
template<typename T1 > | |
int | SetComplexShift (T1 diagonal_shift) |
Set the diagonal complex shift for some preconditioner options. More... | |
template<typename T1 > | |
int | GetComplexShift (T1 &diagonal_shift) |
Get the diagonal complex shift for some preconditioner options. More... | |
int & | GetSeparatorNumSubdomains () |
Get number of subdomains on the top level if setted in advange. More... | |
vector_int & | GetSeparatorDomi () |
Number of subdomain for each node, -1 is the separator. More... | |
Public Attributes | |
std::conditional< PargemslrIsDoublePrecision< T >::value, double, float >::type | _diagonal_shift |
The diagonal complex shift in some preconditioners. The preconditioner would be built on A+_diagonal_shift*I. More... | |
CommunicationHelperClass | _comm_helper |
The data struxture for the parallel matvec. Stores communication inforamtion, and working buffers. More... | |
CommunicationHelperClass | _trans_comm_helper |
The data struxture for the parallel matvec of A transpose. Stores communication inforamtion, and working buffers. More... | |
Class of parallel CSR matrices.
template pargemslr::ParallelCsrMatrixClass< T >::ParallelCsrMatrixClass | ( | ) |
The constructor of ParallelCsrMatrixClass. The default memory location is the host memory.
pargemslr::ParallelCsrMatrixClass< T >::ParallelCsrMatrixClass | ( | const ParallelCsrMatrixClass< T > & | mat | ) |
The copy constructor of ParallelCsrMatrixClass.
[in] | mat | The target matrix. |
pargemslr::ParallelCsrMatrixClass< T >::ParallelCsrMatrixClass | ( | ParallelCsrMatrixClass< T > && | mat | ) |
The move constructor of ParallelCsrMatrixClass.
[in] | mat | The target matrix. |
|
virtual |
The destructor of ParallelCsrMatrixClass. Simply a call to the free function.
|
virtual |
Free the current matrix.
|
virtual |
Create an indentity matrix. Should setup the diag_mat and offd_mat first.
|
virtual |
Fill the matrix with constant value.
[in] | v | The value to be filled. |
template long int pargemslr::ParallelCsrMatrixClass< T >::GetColStartGlobal | ( | ) | const |
Get the global number of the first column of the matrix.
|
inline |
Get the diagonal complex shift for some preconditioner options.
[in] | diagonal_shift | The complex shift will be 0+shifti. |
|
virtual |
Get the data location of the matrix.
template CsrMatrixClass< complexd > & pargemslr::ParallelCsrMatrixClass< T >::GetDiagMat | ( | ) |
Get the diagonal matrix.
template int pargemslr::ParallelCsrMatrixClass< T >::GetDiagScale | ( | T & | scale | ) |
Scale the current csr matrix.
[in] | scale | The scale. |
template int pargemslr::ParallelCsrMatrixClass< T >::GetGraphArrays | ( | vector_long & | vtxdist, |
vector_long & | xadj, | ||
vector_long & | adjncy | ||
) |
Get the vector arrays for the graph, used as inputs for ParMETIS. Diagonal entry is removed.
[in] | vtxdist | The vertex distribution of size np+1. |
[in] | xadj | The local I in CSR format. |
[in] | adjncy | The local J in CSR format. |
template long int pargemslr::ParallelCsrMatrixClass< T >::GetNumColsGlobal | ( | ) | const |
Get the global number of columns of the matrix.
|
virtual |
Get the local number of columns of the matrix.
|
virtual |
Get the number of nonzeros in this matrix.
template long int pargemslr::ParallelCsrMatrixClass< T >::GetNumRowsGlobal | ( | ) | const |
Get the global number of rows of the matrix.
|
virtual |
Get the local number of rows of the matrix.
template vector_long & pargemslr::ParallelCsrMatrixClass< T >::GetOffdMap | ( | ) |
Get the off-diagonal column map vector.
template CsrMatrixClass< complexd > & pargemslr::ParallelCsrMatrixClass< T >::GetOffdMat | ( | ) |
Get the off-diagonal matrix.
template long int pargemslr::ParallelCsrMatrixClass< T >::GetRowStartGlobal | ( | ) | const |
Get the global number of the first row of the matrix.
|
inline |
Number of subdomain for each node, -1 is the separator.
|
inline |
Get number of subdomains on the top level if setted in advange.
template int pargemslr::ParallelCsrMatrixClass< T >::Helmholtz | ( | int | n, |
int | dx, | ||
int | dy, | ||
int | dz, | ||
T | w, | ||
parallel_log & | parlog | ||
) |
Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = f. Robin boundary condition u'-iwu = 0 is used.
[in] | n | Problem size on each dimension |
[in] | dx | |
[in] | dy | |
[in] | dz | the decomposition of the problem on the z dimension. Should have np == dx*dy*dz. |
[in] | w | The w value. |
[in] | parlog | The Parallel log data structure, if parallel_log._comm == NULL, will use the global one. |
template int pargemslr::ParallelCsrMatrixClass< T >::HelmholtzWithPartition | ( | int | n, |
int | dx, | ||
int | dy, | ||
int | dz, | ||
int | d2x, | ||
int | d2y, | ||
int | d2z, | ||
T | w, | ||
parallel_log & | parlog | ||
) |
Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = f. Robin boundary condition u'-iwu = 0 is used.
[in] | n | Problem size on each dimension |
[in] | dx | |
[in] | dy | |
[in] | dz | the decomposition of the problem on each dimension. Should have np == dx*dy*dz. |
[in] | d2x | |
[in] | d2y | |
[in] | d2z | the decomposition of the problem on each dimension. Should have d2x*d2y*d2z % np == 0. |
[in] | w | The w value. |
[in] | parlog | The Parallel log data structure, if parallel_log._comm == NULL, will use the global one. |
template bool & pargemslr::ParallelCsrMatrixClass< T >::IsOffdMapSorted | ( | ) |
Tell if the off-diagonal map vector is sorted.
template int pargemslr::ParallelCsrMatrixClass< T >::Laplacian | ( | int | nx, |
int | ny, | ||
int | nz, | ||
int | dx, | ||
int | dy, | ||
int | dz, | ||
T | alphax, | ||
T | alphay, | ||
T | alphaz, | ||
T | shift, | ||
parallel_log & | parlog, | ||
bool | rand_perturb = false |
||
) |
Generate Laplacian matrix on the host memory, 5-pt for 2D problem and 7-pt for 3D problem.
[in] | nx | |
[in] | ny | |
[in] | nz | Problem size on each dimension |
[in] | dx | |
[in] | dy | |
[in] | dz | the decomposition of the problem on each dimension. Should have np == dx*dy*dz. |
[in] | alphax | |
[in] | alphay | |
[in] | alphaz | Alpha value on each direction |
[in] | shift | The shift value on diagonal |
[in] | parlog | The Parallel log data structure, if parallel_log._comm == NULL, will use the global one |
[in] | rand_perturb | Set to true to add extra random real pertubation to diagonal entries as abs(shift)*Rand(). We only add shift that might makes the problem harder. |
template int pargemslr::ParallelCsrMatrixClass< T >::LaplacianWithPartition | ( | int | nx, |
int | ny, | ||
int | nz, | ||
int | dx, | ||
int | dy, | ||
int | dz, | ||
int | d2x, | ||
int | d2y, | ||
int | d2z, | ||
T | alphax, | ||
T | alphay, | ||
T | alphaz, | ||
T | shift, | ||
parallel_log & | parlog, | ||
bool | rand_perturb = false |
||
) |
Generate Laplacian matrix on the host memory, 5-pt for 2D problem and 7-pt for 3D problem.
Meth partition is created in this case based on infor in d2x/y/z.
[in] | nx | |
[in] | ny | |
[in] | nz | Problem size on each dimension |
[in] | dx | |
[in] | dy | |
[in] | dz | the decomposition of the processors on each dimension. Should have dx*dy*dz == np. |
[in] | d2x | |
[in] | d2y | |
[in] | d2z | the decomposition of the problem on each dimension. Should have d2x*d2y*d2z % np == 0. |
[in] | alphax | |
[in] | alphay | |
[in] | alphaz | Alpha value on each direction |
[in] | shift | The shift value on diagonal |
[in] | parlog | The Parallel log data structure, if parallel_log._comm == NULL, will use the global one. |
[in] | rand_perturb | Set to true to add extra random real pertubation to diagonal entries as abs(shift)*Rand(). We only add shift that might makes the problem harder. |
template int pargemslr::ParallelCsrMatrixClass< T >::MatMat | ( | const T & | alpha, |
const ParallelCsrMatrixClass< T > & | A, | ||
char | transa, | ||
const ParallelCsrMatrixClass< T > & | B, | ||
char | transb, | ||
const T & | beta | ||
) |
Dense Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C where C is this matrix.
If C is not pre-allocated, C will be allocated such that it has the same location as A.
[in] | alpha | The alpha value. |
[in] | A | The first csr matrix. |
[in] | transa | Whether we transpose A or not. |
[in] | B | The second csr matrix. |
[in] | transb | Whether we transpose B or not. |
[in] | beta | The beta value. |
|
virtual |
In place Matrix-Vector product ==> z := alpha*A*x + beta*y, or z := alpha*A'*x + beta*y. Note that z should not equal to y.
[in] | trans | Whether or not transpose matrix A. |
[in] | alpha | The alpha value. |
[in] | x | The first vector. |
[in] | beta | The beta value. |
[in] | y | The second vector. |
[out] | z | The output vector. |
|
virtual |
In place parallel csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y.
[in] | trans | Whether or not transpose matrix A. |
[in] | alpha | The alpha value. |
[in] | x | The left vector. |
[in] | beta | The beta value. |
[in,out] | y | The product vector. |
template int pargemslr::ParallelCsrMatrixClass< T >::MatVecOffd | ( | char | trans, |
const T & | alpha, | ||
const VectorClass< T > & | x, | ||
const T & | beta, | ||
VectorClass< T > & | y | ||
) |
In place parallel csr Matrix-Vector product with offdiagonal block only ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y.
[in] | trans | Whether or not transpose matrix A. |
[in] | alpha | The alpha value. |
[in] | x | The left vector. |
[in] | beta | The beta value. |
[in,out] | y | The product vector. |
|
virtual |
Move the data to another memory location.
[in] | location | The location move to. |
ParallelCsrMatrixClass< T > & pargemslr::ParallelCsrMatrixClass< T >::operator= | ( | const ParallelCsrMatrixClass< T > & | mat | ) |
The = operator of ParallelCsrMatrixClass.
[in] | mat | The target matrix. |
ParallelCsrMatrixClass< T > & pargemslr::ParallelCsrMatrixClass< T >::operator= | ( | ParallelCsrMatrixClass< T > && | mat | ) |
The = operator of ParallelCsrMatrixClass.
[in] | mat | The target matrix. |
template int pargemslr::ParallelCsrMatrixClass< T >::PlotPatternGnuPlot | ( | const char * | datafilename, |
int | pttype = 0 |
||
) |
Plot the pattern of the parallel csr matrix using gnuplot. Similar to spy in the MATLAB. Function for testing purpose.
[in] | datafilename | The filename of the temp file holding the data. |
template int pargemslr::ParallelCsrMatrixClass< T >::ReadFromSingleMMFile | ( | const char * | matfile, |
int | idxin, | ||
parallel_log & | parlog | ||
) |
Read matrix on the host memory, matrix market format.
[in] | matfile | The file. |
[in] | idxin | The index base of the input file. 0-based or 1-based. |
[in] | parlog | The parallel_log data structure. |
|
virtual |
Scale the current csr matrix.
[in] | alpha | The scale. |
|
inline |
Set the diagonal complex shift for some preconditioner options.
[in] | diagonal_shift | The complex shift will be 0+shifti. |
template int pargemslr::ParallelCsrMatrixClass< T >::SetOffdMatSorted | ( | bool | is_offd_map_sorted | ) |
Set the off-diagonal map vector to sorted/unsorted.
template int pargemslr::ParallelCsrMatrixClass< T >::Setup | ( | int | nrow_local, |
int | ncol_local, | ||
parallel_log & | parlog | ||
) |
Free the current parallel matrix, and set up new paralle matrix.
[in] | nrow_local | The local number of rows. |
[in] | ncol_local | The local number of cols. |
[in] | parlog | The Parallel log data structure, if parallel_log._comm == NULL, will use the global one. |
template int pargemslr::ParallelCsrMatrixClass< T >::Setup | ( | int | nrow_local, |
long int | nrow_start, | ||
long int | nrow_global, | ||
int | ncol_local, | ||
long int | ncol_start, | ||
long int | ncol_global, | ||
parallel_log & | parlog | ||
) |
Free the current parallel matrix, and set up new paralle matrix. The default nnz is 0.
[in] | nrow_local | The local number of rows. |
[in] | nrow_start | The start index of the rows. |
[in] | nrow_global | The global length of the rows. |
[in] | ncol_local | The local number of cols. |
[in] | ncol_start | The start index of the columns. |
[in] | ncol_global | The global length of the columns. |
[in] | parlog | The Parallel log data structure, if parallel_log._comm == NULL, will use the global one. |
template int pargemslr::ParallelCsrMatrixClass< T >::SetupMatvec | ( | ) |
Setup the matvec_helper for parallel matvec.
template int pargemslr::ParallelCsrMatrixClass< T >::SetupMatvecOver | ( | ) |
Finishing up the setup phase of matvec_helper for parallel matvec.
template int pargemslr::ParallelCsrMatrixClass< T >::SetupMatvecStart | ( | ) |
Setup the matvec_helper for parallel matvec.
Non-blocking communication, need to call SetupMatvecOver function to finish up.
template int pargemslr::ParallelCsrMatrixClass< T >::SetupVectorPtrStr | ( | ParallelVectorClass< T > & | vec | ) |
Update the structure of a parallel vector to have same distribution as the row distribution of this matrix.
[in,out] | The | target vector. |
template int pargemslr::ParallelCsrMatrixClass< T >::SortOffdMap | ( | ) |
Sort the off-diagonal map vector.
template int pargemslr::ParallelCsrMatrixClass< T >::SubMatrix | ( | vector_long & | rows, |
vector_long & | cols, | ||
int | location, | ||
ParallelCsrMatrixClass< T > & | parcsrmat_out | ||
) |
Extract parallel submatrix from parallel CSR matrix.
Each processor prvide several local rows/cols to keep (in global indices). No exchange of global information.
[in] | rows | The rows to keep. The permutation in rows also applied. |
[in] | cols | The columns to keep. The permutation in cols also applied. |
[in] | location | The location of the output matrix. |
[out] | parcsrmat_out | The output matrix. |
template int pargemslr::ParallelCsrMatrixClass< T >::Transpose | ( | ) |
Transpost the current matrix.
int pargemslr::ParallelCsrMatrixClass< T >::Transpose | ( | ParallelCsrMatrixClass< T > & | AT | ) |
Transpost the current matrix.
[out] | AT | The output matrix. |
CommunicationHelperClass pargemslr::ParallelCsrMatrixClass< T >::_comm_helper |
The data struxture for the parallel matvec. Stores communication inforamtion, and working buffers.
std::conditional<PargemslrIsDoublePrecision<T>::value, double, float>::type pargemslr::ParallelCsrMatrixClass< T >::_diagonal_shift |
The diagonal complex shift in some preconditioners. The preconditioner would be built on A+_diagonal_shift*I.
CommunicationHelperClass pargemslr::ParallelCsrMatrixClass< T >::_trans_comm_helper |
The data struxture for the parallel matvec of A transpose. Stores communication inforamtion, and working buffers.