ParGeMSLR
Public Member Functions | Public Attributes | List of all members
pargemslr::ParallelCsrMatrixClass< T > Class Template Reference

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_longGetOffdMap ()
 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_intGetSeparatorDomi ()
 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...
 

Detailed Description

template<typename T>
class pargemslr::ParallelCsrMatrixClass< T >

Class of parallel CSR matrices.

Constructor & Destructor Documentation

◆ ParallelCsrMatrixClass() [1/3]

template<typename T >
template pargemslr::ParallelCsrMatrixClass< T >::ParallelCsrMatrixClass ( )

The constructor of ParallelCsrMatrixClass. The default memory location is the host memory.

◆ ParallelCsrMatrixClass() [2/3]

template<typename T >
pargemslr::ParallelCsrMatrixClass< T >::ParallelCsrMatrixClass ( const ParallelCsrMatrixClass< T > &  mat)

The copy constructor of ParallelCsrMatrixClass.

Parameters
[in]matThe target matrix.

◆ ParallelCsrMatrixClass() [3/3]

The move constructor of ParallelCsrMatrixClass.

Parameters
[in]matThe target matrix.

◆ ~ParallelCsrMatrixClass()

template<typename T >
template pargemslr::ParallelCsrMatrixClass< T >::~ParallelCsrMatrixClass ( )
virtual

The destructor of ParallelCsrMatrixClass. Simply a call to the free function.

Member Function Documentation

◆ Clear()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::Clear ( )
virtual

Free the current matrix.

Returns
Return error message.

◆ Eye()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::Eye ( )
virtual

Create an indentity matrix. Should setup the diag_mat and offd_mat first.

Returns
Return error message.

◆ Fill()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::Fill ( const T &  v)
virtual

Fill the matrix with constant value.

Parameters
[in]vThe value to be filled.
Returns
Return error message.

◆ GetColStartGlobal()

template<typename T >
template long int pargemslr::ParallelCsrMatrixClass< T >::GetColStartGlobal ( ) const

Get the global number of the first column of the matrix.

Returns
Return the global number of the first column of the matrix.

◆ GetComplexShift()

template<typename T >
template<typename T1 >
int pargemslr::ParallelCsrMatrixClass< T >::GetComplexShift ( T1 &  diagonal_shift)
inline

Get the diagonal complex shift for some preconditioner options.

Parameters
[in]diagonal_shiftThe complex shift will be 0+shifti.
Returns
Return error message.

◆ GetDataLocation()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::GetDataLocation ( ) const
virtual

Get the data location of the matrix.

Returns
Return the locatino of the matrix.

◆ GetDiagMat()

template<typename T >
template CsrMatrixClass< complexd > & pargemslr::ParallelCsrMatrixClass< T >::GetDiagMat ( )

Get the diagonal matrix.

Returns
Return the diagonal matrix.

◆ GetDiagScale()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::GetDiagScale ( T &  scale)

Scale the current csr matrix.

Parameters
[in]scaleThe scale.
Returns
Return error message.

◆ GetGraphArrays()

template<typename T >
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.

Parameters
[in]vtxdistThe vertex distribution of size np+1.
[in]xadjThe local I in CSR format.
[in]adjncyThe local J in CSR format.
Returns
Return error message.

◆ GetNumColsGlobal()

template<typename T >
template long int pargemslr::ParallelCsrMatrixClass< T >::GetNumColsGlobal ( ) const

Get the global number of columns of the matrix.

Returns
Return the global number of columns of the matrix.

◆ GetNumColsLocal()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::GetNumColsLocal ( ) const
virtual

Get the local number of columns of the matrix.

Returns
Return the local number of columns of the matrix.

◆ GetNumNonzeros()

template<typename T >
template long int pargemslr::ParallelCsrMatrixClass< T >::GetNumNonzeros ( ) const
virtual

Get the number of nonzeros in this matrix.

Returns
Return the number of nonzeros of the matrix.

◆ GetNumRowsGlobal()

template<typename T >
template long int pargemslr::ParallelCsrMatrixClass< T >::GetNumRowsGlobal ( ) const

Get the global number of rows of the matrix.

Returns
Return the global number of rows of the matrix.

◆ GetNumRowsLocal()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::GetNumRowsLocal ( ) const
virtual

Get the local number of rows of the matrix.

Returns
Return the local number of rows of the matrix.

◆ GetOffdMap()

template<typename T >
template vector_long & pargemslr::ParallelCsrMatrixClass< T >::GetOffdMap ( )

Get the off-diagonal column map vector.

Returns
Return the off-diagonal column map vector.

◆ GetOffdMat()

template<typename T >
template CsrMatrixClass< complexd > & pargemslr::ParallelCsrMatrixClass< T >::GetOffdMat ( )

Get the off-diagonal matrix.

Returns
Return the off-diagonal matrix.

◆ GetRowStartGlobal()

template<typename T >
template long int pargemslr::ParallelCsrMatrixClass< T >::GetRowStartGlobal ( ) const

Get the global number of the first row of the matrix.

Returns
Return the global number of the first row of the matrix.

◆ GetSeparatorDomi()

template<typename T >
vector_int& pargemslr::ParallelCsrMatrixClass< T >::GetSeparatorDomi ( )
inline

Number of subdomain for each node, -1 is the separator.

Returns
Return the number of subdomain for each node, -1 is the separator.

◆ GetSeparatorNumSubdomains()

template<typename T >
int& pargemslr::ParallelCsrMatrixClass< T >::GetSeparatorNumSubdomains ( )
inline

Get number of subdomains on the top level if setted in advange.

Returns
Return the number of subdomains.

◆ Helmholtz()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::Helmholtz ( int  n,
int  dx,
int  dy,
int  dz,
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.

Parameters
[in]nProblem size on each dimension
[in]dx
[in]dy
[in]dzthe decomposition of the problem on the z dimension. Should have np == dx*dy*dz.
[in]wThe w value.
[in]parlogThe Parallel log data structure, if parallel_log._comm == NULL, will use the global one.
Returns
Return error message.

◆ HelmholtzWithPartition()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::HelmholtzWithPartition ( int  n,
int  dx,
int  dy,
int  dz,
int  d2x,
int  d2y,
int  d2z,
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.

Parameters
[in]nProblem size on each dimension
[in]dx
[in]dy
[in]dzthe decomposition of the problem on each dimension. Should have np == dx*dy*dz.
[in]d2x
[in]d2y
[in]d2zthe decomposition of the problem on each dimension. Should have d2x*d2y*d2z % np == 0.
[in]wThe w value.
[in]parlogThe Parallel log data structure, if parallel_log._comm == NULL, will use the global one.
Returns
Return error message.

◆ IsOffdMapSorted()

template<typename T >
template bool & pargemslr::ParallelCsrMatrixClass< T >::IsOffdMapSorted ( )

Tell if the off-diagonal map vector is sorted.

Returns
Return true if off-diagonal map is sorted.

◆ Laplacian()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::Laplacian ( int  nx,
int  ny,
int  nz,
int  dx,
int  dy,
int  dz,
alphax,
alphay,
alphaz,
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.

Parameters
[in]nx
[in]ny
[in]nzProblem size on each dimension
[in]dx
[in]dy
[in]dzthe decomposition of the problem on each dimension. Should have np == dx*dy*dz.
[in]alphax
[in]alphay
[in]alphazAlpha value on each direction
[in]shiftThe shift value on diagonal
[in]parlogThe Parallel log data structure, if parallel_log._comm == NULL, will use the global one
[in]rand_perturbSet 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.
Returns
Return error message.

◆ LaplacianWithPartition()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::LaplacianWithPartition ( int  nx,
int  ny,
int  nz,
int  dx,
int  dy,
int  dz,
int  d2x,
int  d2y,
int  d2z,
alphax,
alphay,
alphaz,
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.

Parameters
[in]nx
[in]ny
[in]nzProblem size on each dimension
[in]dx
[in]dy
[in]dzthe decomposition of the processors on each dimension. Should have dx*dy*dz == np.
[in]d2x
[in]d2y
[in]d2zthe decomposition of the problem on each dimension. Should have d2x*d2y*d2z % np == 0.
[in]alphax
[in]alphay
[in]alphazAlpha value on each direction
[in]shiftThe shift value on diagonal
[in]parlogThe Parallel log data structure, if parallel_log._comm == NULL, will use the global one.
[in]rand_perturbSet 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.
Returns
Return error message.

◆ MatMat()

template<typename T >
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.

Note
TODO: currently not supported.
Parameters
[in]alphaThe alpha value.
[in]AThe first csr matrix.
[in]transaWhether we transpose A or not.
[in]BThe second csr matrix.
[in]transbWhether we transpose B or not.
[in]betaThe beta value.
Returns
Return error message.

◆ MatVec() [1/2]

template<typename T >
int pargemslr::ParallelCsrMatrixClass< T >::MatVec ( char  trans,
const T &  alpha,
const VectorClass< T > &  x,
const T &  beta,
const VectorClass< T > &  y,
VectorClass< T > &  z 
)
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.

Note
TODO: Currently A' is not supported.
Parameters
[in]transWhether or not transpose matrix A.
[in]alphaThe alpha value.
[in]xThe first vector.
[in]betaThe beta value.
[in]yThe second vector.
[out]zThe output vector.
Returns
Return error message.

◆ MatVec() [2/2]

template<typename T >
int pargemslr::ParallelCsrMatrixClass< T >::MatVec ( char  trans,
const T &  alpha,
const VectorClass< T > &  x,
const T &  beta,
VectorClass< T > &  y 
)
virtual

In place parallel csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y.

Note
TODO: Currently A' is not supported.
Parameters
[in]transWhether or not transpose matrix A.
[in]alphaThe alpha value.
[in]xThe left vector.
[in]betaThe beta value.
[in,out]yThe product vector.
Returns
Return error message.

◆ MatVecOffd()

template<typename T >
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.

Parameters
[in]transWhether or not transpose matrix A.
[in]alphaThe alpha value.
[in]xThe left vector.
[in]betaThe beta value.
[in,out]yThe product vector.
Returns
Return error message.

◆ MoveData()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::MoveData ( const int &  location)
virtual

Move the data to another memory location.

Parameters
[in]locationThe location move to.
Returns
Return error message.

◆ operator=() [1/2]

template<typename T >
ParallelCsrMatrixClass< T > & pargemslr::ParallelCsrMatrixClass< T >::operator= ( const ParallelCsrMatrixClass< T > &  mat)

The = operator of ParallelCsrMatrixClass.

Parameters
[in]matThe target matrix.
Returns
Return the matrix.

◆ operator=() [2/2]

template<typename T >
ParallelCsrMatrixClass< T > & pargemslr::ParallelCsrMatrixClass< T >::operator= ( ParallelCsrMatrixClass< T > &&  mat)

The = operator of ParallelCsrMatrixClass.

Parameters
[in]matThe target matrix.
Returns
Return the matrix.

◆ PlotPatternGnuPlot()

template<typename T >
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.

Parameters
[in]datafilenameThe filename of the temp file holding the data.
Returns
Return error message.

◆ ReadFromSingleMMFile()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::ReadFromSingleMMFile ( const char *  matfile,
int  idxin,
parallel_log parlog 
)

Read matrix on the host memory, matrix market format.

Note
The current version doesn't support long int as input.
Parameters
[in]matfileThe file.
[in]idxinThe index base of the input file. 0-based or 1-based.
[in]parlogThe parallel_log data structure.
Returns
Return error message.

◆ Scale()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::Scale ( const T &  alpha)
virtual

Scale the current csr matrix.

Parameters
[in]alphaThe scale.
Returns
Return error message.

◆ SetComplexShift()

template<typename T >
template<typename T1 >
int pargemslr::ParallelCsrMatrixClass< T >::SetComplexShift ( T1  diagonal_shift)
inline

Set the diagonal complex shift for some preconditioner options.

Parameters
[in]diagonal_shiftThe complex shift will be 0+shifti.
Returns
Return error message.

◆ SetOffdMatSorted()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::SetOffdMatSorted ( bool  is_offd_map_sorted)

Set the off-diagonal map vector to sorted/unsorted.

Returns
Return error message.

◆ Setup() [1/2]

template<typename T >
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.

Note
The diag_mat and offd_mat need to be initialized separately.
Parameters
[in]nrow_localThe local number of rows.
[in]ncol_localThe local number of cols.
[in]parlogThe Parallel log data structure, if parallel_log._comm == NULL, will use the global one.
Returns
Return error message.

◆ Setup() [2/2]

template<typename T >
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.

Note
The diag_mat and offd_mat need to be initialized separately.
Parameters
[in]nrow_localThe local number of rows.
[in]nrow_startThe start index of the rows.
[in]nrow_globalThe global length of the rows.
[in]ncol_localThe local number of cols.
[in]ncol_startThe start index of the columns.
[in]ncol_globalThe global length of the columns.
[in]parlogThe Parallel log data structure, if parallel_log._comm == NULL, will use the global one.
Returns
Return error message.

◆ SetupMatvec()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::SetupMatvec ( )

Setup the matvec_helper for parallel matvec.

Returns
Return error message.

◆ SetupMatvecOver()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::SetupMatvecOver ( )

Finishing up the setup phase of matvec_helper for parallel matvec.

Returns
Return error message.

◆ SetupMatvecStart()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::SetupMatvecStart ( )

Setup the matvec_helper for parallel matvec.
Non-blocking communication, need to call SetupMatvecOver function to finish up.

Returns
Return error message.

◆ SetupVectorPtrStr()

template<typename T >
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.

Parameters
[in,out]Thetarget vector.
Returns
Return error message.

◆ SortOffdMap()

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::SortOffdMap ( )

Sort the off-diagonal map vector.

Returns
Return error message.

◆ SubMatrix()

template<typename T >
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.

Parameters
[in]rowsThe rows to keep. The permutation in rows also applied.
[in]colsThe columns to keep. The permutation in cols also applied.
[in]locationThe location of the output matrix.
[out]parcsrmat_outThe output matrix.
Returns
Return error message.

◆ Transpose() [1/2]

template<typename T >
template int pargemslr::ParallelCsrMatrixClass< T >::Transpose ( )

Transpost the current matrix.

Returns
Return error message.

◆ Transpose() [2/2]

template<typename T >
int pargemslr::ParallelCsrMatrixClass< T >::Transpose ( ParallelCsrMatrixClass< T > &  AT)

Transpost the current matrix.

Parameters
[out]ATThe output matrix.
Returns
Return error message.

Member Data Documentation

◆ _comm_helper

template<typename T >
CommunicationHelperClass pargemslr::ParallelCsrMatrixClass< T >::_comm_helper

The data struxture for the parallel matvec. Stores communication inforamtion, and working buffers.

◆ _diagonal_shift

template<typename T >
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.

Note
TODO: if useful, change to private.

◆ _trans_comm_helper

template<typename T >
CommunicationHelperClass pargemslr::ParallelCsrMatrixClass< T >::_trans_comm_helper

The data struxture for the parallel matvec of A transpose. Stores communication inforamtion, and working buffers.

Note
Not yet in use.

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