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

Class of CSR/CSC matrices. More...

#include <csr_matrix.hpp>

Public Member Functions

 CsrMatrixClass ()
 The constructor of CsrMatrixClass. More...
 
 CsrMatrixClass (const CsrMatrixClass< T > &mat)
 The copy constructor of CsrMatrixClass. More...
 
 CsrMatrixClass (CsrMatrixClass< T > &&mat)
 The move constructor of CsrMatrixClass. More...
 
CsrMatrixClass< T > & operator= (const CsrMatrixClass< T > &mat)
 The = operator of CsrMatrixClass. More...
 
CsrMatrixClass< T > & operator= (CsrMatrixClass< T > &&mat)
 The = operator of CsrMatrixClass. More...
 
virtual ~CsrMatrixClass ()
 The destructor of CsrMatrixClass. More...
 
int Setup (int nrows, int ncols, int nnz)
 Free the current matrix, and allocate memory to create a new matrix. More...
 
int Setup (int nrows, int ncols, int nnz, int location)
 Free the current matrix, and allocate memory to create a new matrix. More...
 
int Setup (int nrows, int ncols, int nnz, bool setzero)
 Free the current matrix, and allocate memory to create a new matrix. More...
 
int Setup (int nrows, int ncols, int nnz, bool holdvalue, bool setzero)
 Free the current matrix, and allocate memory to create a new matrix. More...
 
int Setup (int nrows, int ncols, int nnz, int location, bool setzero)
 Free the current matrix, and allocate memory to create a new matrix. More...
 
int Setup (int nrows, int ncols, int nnz, int location, bool holdvalue, bool setzero)
 Free the current matrix, and allocate memory to create a new matrix. More...
 
int Setup (int nrows, int ncols, int nnz, int location, bool holdvalue, bool setzero, bool iscsr)
 Free the current matrix, and allocate memory to create a new matrix. More...
 
int PushBack (int col, T v)
 Insert value at the end of the J and A vector, user need to update the I vector. More...
 
int SubMatrix (int row_start, int col_start, int num_rows, int num_cols, int location, CsrMatrixClass< T > &csrmat_out)
 Copy data to extract a submatrix of the current matrix. More...
 
int SubMatrix (vector_int &row_starts, vector_int &col_starts, vector_int &row_ends, vector_int &col_ends, int location, CsrMatrixClass< T > &csrmat_out)
 Copy data to extract a submatrix of the current matrix. More...
 
int SubMatrix (vector_int &rows, vector_int &cols, int location, CsrMatrixClass< T > &csrmat_out)
 Copy data to extract a submatrix of the current matrix. More...
 
int SubMatrixNoPerm (vector_int &rows, vector_int &cols, vector_int &row_perm, vector_int &col_perm, bool complement, int location, CsrMatrixClass< T > &csrmat_out)
 Copy data to extract a submatrix of the current matrix. More...
 
int GetConnectedComponents (std::vector< vector_int > &comp_indices, int &ncomps)
 Get the connected components of A + AT, where A is the current matrix. More...
 
int SetupVectorPtrStr (SequentialVectorClass< T > &vec)
 Update the structure of a vector to have same row permutation. 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...
 
virtual long int GetNumNonzeros () const
 Get the number of nonzeros in this matrix. More...
 
int SetNumNonzeros ()
 Set the number of nonzeros in this matrix. More...
 
int * GetI () const
 Get the I pointer of the matrix. More...
 
int * GetJ () const
 Get the J pointer of the matrix. More...
 
T * GetData () const
 Get the data pointer of the matrix. More...
 
IntVectorClass< int > & GetIVector ()
 Get the reference to the I vector. More...
 
const IntVectorClass< int > & GetIVector () const
 Get the reference to the I vector. More...
 
IntVectorClass< int > & GetJVector ()
 Get the reference to the J vector. More...
 
const IntVectorClass< int > & GetJVector () const
 Get the reference to the J vector. More...
 
SequentialVectorClass< T > & GetDataVector ()
 Get the reference to the data vector. More...
 
const SequentialVectorClass< T > & GetDataVector () const
 Get the reference to the data vector. More...
 
bool IsHoldingData () const
 Check if the matrix holding its own data. More...
 
bool & IsCsr ()
 Csr or csc. More...
 
const bool & IsCsr () const
 Csr or csc. More...
 
bool & IsRowSorted ()
 Tell if the column/row indices in each row/column is sorted. More...
 
const bool & IsRowSorted () const
 Tell if the column/row indices in each row/column is sorted. More...
 
int SortRow ()
 Sort the column/row indices in each row/column. 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 Convert (bool csr)
 Convert between csr and csc. More...
 
int Transpose ()
 Transpost the current matrix. More...
 
int Transpose (CsrMatrixClass< T > &AT)
 Transpost the current matrix. More...
 
virtual int MatVec (char trans, const T &alpha, const VectorClass< T > &x, const T &beta, VectorClass< T > &y)
 In place 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. More...
 
int MatMat (const T &alpha, const CsrMatrixClass< T > &A, char transa, const CsrMatrixClass< 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 PlotPattern (int *rperm, int *cperm, int conditiona, int conditionb)
 Plot the pattern of the csr matrix to the terminal output. Similar to spy in the MATLAB. Function for testing purpose. More...
 
int PlotPatternGnuPlot (const char *datafilename, int *rperm, int *cperm, int conditiona, int conditionb)
 Plot the pattern of the csr matrix using gnuplot. Similar to spy in the MATLAB. Function for testing purpose. More...
 
int Plot (int *perm, int conditiona, int conditionb, int width)
 Plot the csr matrix to the terminal output. Function for testing purpose. More...
 
int Laplacian (int nx, int ny, int nz, T alphax, T alphay, T alphaz, T shift, 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, T w)
 Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = 0. More...
 
int ReadFromMMFile (const char *matfile, int idxin)
 Read matrix on the host memory, matrix market format. 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...
 
template<typename T1 >
int GetComplexShift (ComplexValueClass< T1 > &diagonal_shift)
 Get the diagonal complex shift for some preconditioner options. More...
 
int GetMpiInfo (int &np, int &myid, MPI_Comm &comm) const
 Get comm, np, and myid. Get the sequential comm. More...
 
MPI_Comm GetComm () const
 Get the MPI_comm. 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...
 

Detailed Description

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

Class of CSR/CSC matrices.

Constructor & Destructor Documentation

◆ CsrMatrixClass() [1/3]

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

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

◆ CsrMatrixClass() [2/3]

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

The copy constructor of CsrMatrixClass.

Parameters
[in]matThe target matrix.

◆ CsrMatrixClass() [3/3]

template<typename T >
pargemslr::CsrMatrixClass< T >::CsrMatrixClass ( CsrMatrixClass< T > &&  mat)

The move constructor of CsrMatrixClass.

Parameters
[in]matThe target matrix.

◆ ~CsrMatrixClass()

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

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

Member Function Documentation

◆ Clear()

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

Free the current matrix.

Returns
Return error message.

◆ Convert()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Convert ( bool  csr)

Convert between csr and csc.

Parameters
[in]csrSet to true to use csr, otherwise csc.
Returns
Return error message.

◆ Eye()

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

Create an indentity matrix.

Returns
Return error message.

◆ Fill()

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

Fill the matrix with constant value.

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

◆ GetComm()

template<typename T >
MPI_Comm pargemslr::CsrMatrixClass< T >::GetComm ( ) const
inline

Get the MPI_comm.

Returns
Return the MPI_comm.

◆ GetComplexShift() [1/2]

template<typename T >
template<typename T1 >
int pargemslr::CsrMatrixClass< T >::GetComplexShift ( ComplexValueClass< T1 > &  diagonal_shift)

Get the diagonal complex shift for some preconditioner options.

Parameters
[out]diagonal_shiftThe complex shift.
Returns
Return error message.

◆ GetComplexShift() [2/2]

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

Get the diagonal complex shift for some preconditioner options.

Parameters
[out]diagonal_shiftThe complex shift.
Returns
Return error message.

◆ GetConnectedComponents()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::GetConnectedComponents ( std::vector< vector_int > &  comp_indices,
int &  ncomps 
)

Get the connected components of A + AT, where A is the current matrix.

Parameters
[out]comp_indices2D vector, comp_indices[i] is the index of the elements in component i.
[in,out]ncompsOn entry, If max_ncomps<=1, this parameter will be ignored. If max_ncomps>1, when the actual number of connected component is greater than max_ncomps, we will merge several connected components together.
On exit, the number of final components (after merge).
Returns
Return error message.

◆ GetData()

template<typename T >
template complexd * pargemslr::CsrMatrixClass< T >::GetData ( ) const

Get the data pointer of the matrix.

Returns
Return the data pointer.

◆ GetDataLocation()

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

Get the data location of the matrix.

Returns
Return the data location of the matrix.

◆ GetDataVector() [1/2]

template<typename T >
template SequentialVectorClass< complexd > & pargemslr::CsrMatrixClass< T >::GetDataVector ( )

Get the reference to the data vector.

Returns
Return the reference to the data vector.

◆ GetDataVector() [2/2]

template<typename T >
template const SequentialVectorClass< complexd > & pargemslr::CsrMatrixClass< T >::GetDataVector ( ) const

Get the reference to the data vector.

Returns
Return the reference to the data vector.

◆ GetDiagScale()

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

Scale the current csr matrix.

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

◆ GetI()

template<typename T >
template int * pargemslr::CsrMatrixClass< T >::GetI ( ) const

Get the I pointer of the matrix.

Returns
Return the I pointer.

◆ GetIVector() [1/2]

template<typename T >
template IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetIVector ( )

Get the reference to the I vector.

Returns
Return the reference to the I vector.

◆ GetIVector() [2/2]

template<typename T >
template const IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetIVector ( ) const

Get the reference to the I vector.

Returns
Return the reference to the I vector.

◆ GetJ()

template<typename T >
template int * pargemslr::CsrMatrixClass< T >::GetJ ( ) const

Get the J pointer of the matrix.

Returns
Return the J pointer.

◆ GetJVector() [1/2]

template<typename T >
template IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetJVector ( )

Get the reference to the J vector.

Returns
Return the reference to the J vector.

◆ GetJVector() [2/2]

template<typename T >
template const IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetJVector ( ) const

Get the reference to the J vector.

Returns
Return the reference to the J vector.

◆ GetMpiInfo()

template<typename T >
int pargemslr::CsrMatrixClass< T >::GetMpiInfo ( int &  np,
int &  myid,
MPI_Comm &  comm 
) const
inline

Get comm, np, and myid. Get the sequential comm.

Parameters
[in]npThe number of processors.
[in]myidThe local MPI rank number.
[in]commThe MPI_Comm.
Returns
Return error message.

◆ GetNumColsLocal()

template<typename T >
template int pargemslr::CsrMatrixClass< 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::CsrMatrixClass< T >::GetNumNonzeros ( ) const
virtual

Get the number of nonzeros in this matrix.

Returns
Return the number of nonzeros of the matrix.

◆ GetNumRowsLocal()

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

Get the local number of rows of the matrix.

Returns
Return the local number of rows of the matrix.

◆ Helmholtz()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Helmholtz ( int  n,
w 
)

Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = 0. Robin boundary condition u'-iwu = 0 is used.

Parameters
[in]nProblem size.
[in]wThe w value.
Returns
Return error message.

◆ IsCsr() [1/2]

template<typename T >
template bool & pargemslr::CsrMatrixClass< T >::IsCsr ( )

Csr or csc.

Returns
Return true if matrix is csr.

◆ IsCsr() [2/2]

template<typename T >
template const bool & pargemslr::CsrMatrixClass< T >::IsCsr ( ) const

Csr or csc.

Returns
Return true if matrix is csr.

◆ IsHoldingData()

template<typename T >
template bool pargemslr::CsrMatrixClass< T >::IsHoldingData ( ) const

Check if the matrix holding its own data.

Returns
Return true if matrix holding data.

◆ IsRowSorted() [1/2]

template<typename T >
template bool & pargemslr::CsrMatrixClass< T >::IsRowSorted ( )

Tell if the column/row indices in each row/column is sorted.

Returns
Return true if col/row indices are sorted.

◆ IsRowSorted() [2/2]

template<typename T >
template const bool & pargemslr::CsrMatrixClass< T >::IsRowSorted ( ) const

Tell if the column/row indices in each row/column is sorted.

Returns
Return reference to a bollean indicates if the col/row indices are sorted.

◆ Laplacian()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Laplacian ( int  nx,
int  ny,
int  nz,
alphax,
alphay,
alphaz,
shift,
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]alphax
[in]alphay
[in]alphazAlpha value on each direction
[in]shiftThe shift value on diagonal
[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::CsrMatrixClass< T >::MatMat ( const T &  alpha,
const CsrMatrixClass< T > &  A,
char  transa,
const CsrMatrixClass< 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.

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

◆ MatVec() [1/2]

template<typename T >
int pargemslr::CsrMatrixClass< 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.

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::CsrMatrixClass< T >::MatVec ( char  trans,
const T &  alpha,
const VectorClass< T > &  x,
const T &  beta,
VectorClass< T > &  y 
)
virtual

In place csr Matrix-Vector product ==> 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::CsrMatrixClass< 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 >
CsrMatrixClass< T > & pargemslr::CsrMatrixClass< T >::operator= ( const CsrMatrixClass< T > &  mat)

The = operator of CsrMatrixClass.

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

◆ operator=() [2/2]

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

The = operator of CsrMatrixClass.

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

◆ Plot()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Plot ( int *  perm,
int  conditiona,
int  conditionb,
int  width 
)

Plot the csr matrix to the terminal output. Function for testing purpose.

Parameters
[in]permThe permutation vector. Set to NULL to plot(A), otherwise plot(A(perm,perm)).
[in]conditionaFirst condition.
[in]conditionbSecend condition, only plot when conditiona == conditionb.
[in]widthThe plot width.
Returns
Return error message.

◆ PlotPattern()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::PlotPattern ( int *  rperm,
int *  cperm,
int  conditiona,
int  conditionb 
)

Plot the pattern of the csr matrix to the terminal output. Similar to spy in the MATLAB. Function for testing purpose.

Parameters
[in]rpermThe row permutation vector.
[in]cpermThe column permutation vector. Set either rperm or cperm to NULL to spy(A), otherwise spy(A(rperm,cperm)).
[in]conditionaFirst condition.
[in]conditionbSecend condition, only spy when conditiona == conditionb.
Returns
Return error message.

◆ PlotPatternGnuPlot()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::PlotPatternGnuPlot ( const char *  datafilename,
int *  rperm,
int *  cperm,
int  conditiona,
int  conditionb 
)

Plot the pattern of the 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.
[in]rpermThe row permutation vector.
[in]cpermThe column permutation vector. Set either rperm or cperm to NULL to spy(A), otherwise spy(A(rperm,cperm)).
[in]conditionaFirst condition.
[in]conditionbSecend condition, only plot when conditiona == conditionb.
Returns
Return error message.

◆ PushBack()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::PushBack ( int  col,
v 
)

Insert value at the end of the J and A vector, user need to update the I vector.
Currently only for host memory.

Parameters
[in]colThe col/row to be inserted into the csr/csc.
[in]vThe value to be inserted.
Returns
Return error message.

◆ ReadFromMMFile()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::ReadFromMMFile ( const char *  matfile,
int  idxin 
)

Read matrix on the host memory, matrix market format.

Parameters
[in]matfileThe file.
[in]idxinThe index base of the input file. 0-based or 1-based.
Returns
Return error message.

◆ Scale()

template<typename T >
template int pargemslr::CsrMatrixClass< 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::CsrMatrixClass< 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.

◆ SetNumNonzeros()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::SetNumNonzeros ( )

Set the number of nonzeros in this matrix.

Returns
Return the number of nonzeros in this matrix.

◆ Setup() [1/7]

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Setup ( int  nrows,
int  ncols,
int  nnz 
)

Free the current matrix, and allocate memory to create a new matrix. The default location of data is the host memory.
If you want to use the PushBack funciton for J and A, please use other Setup function with setzero option.

Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]nnzThe size reserved in J and A array.
Returns
Return error message.

◆ Setup() [2/7]

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Setup ( int  nrows,
int  ncols,
int  nnz,
bool  holdvalue,
bool  setzero 
)

Free the current matrix, and allocate memory to create a new matrix. The default location of data is the host memory.
If setzero is set to true, _i_vec[0] and the length of J and A will be set to 0, so that user can use the PushBack function for vectors.

Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]nnzThe size reserved in J and A array.
[in]holdvalueSet to false to disable the A vector. Only save the pattern of the matrix.
[in]setzeroIf set to true, _i_vec[0] and the length of J and A is set to 0 so you can use the push_back function.
Returns
Return error message.

◆ Setup() [3/7]

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Setup ( int  nrows,
int  ncols,
int  nnz,
bool  setzero 
)

Free the current matrix, and allocate memory to create a new matrix. The default location of data is the host memory.
If setzero is set to true, _i_vec[0] and the length of J and A will be set to 0, so that user can use the PushBack function for vectors.

Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]nnzThe size reserved in J and A array.
[in]setzeroIf set to true, I[0] and the length of J and A are set to 0 so you can use the PushBack function.
Returns
Return error message.

◆ Setup() [4/7]

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Setup ( int  nrows,
int  ncols,
int  nnz,
int  location 
)

Free the current matrix, and allocate memory to create a new matrix.

Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]nnzThe size reserved in J and A array.
[in]locationThe location of the data.
Returns
Return error message.

◆ Setup() [5/7]

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Setup ( int  nrows,
int  ncols,
int  nnz,
int  location,
bool  holdvalue,
bool  setzero 
)

Free the current matrix, and allocate memory to create a new matrix. The default location of data is the host memory.
If setzero is set to true, _i_vec[0] and the length of J and A will be set to 0, so that user can use the PushBack function for vectors.

Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]nnzThe size reserved in J and A array.
[in]locationThe location of the data.
[in]holdvalueSet to false to disable the A vector. Only save the pattern of the matrix.
[in]setzeroIf set to true, _i_vec[0] and the length of J and A is set to 0 so you can use the push_back function.
Returns
Return error message.

◆ Setup() [6/7]

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Setup ( int  nrows,
int  ncols,
int  nnz,
int  location,
bool  holdvalue,
bool  setzero,
bool  iscsr 
)

Free the current matrix, and allocate memory to create a new matrix. The default location of data is the host memory.
If setzero is set to true, _i_vec[0] and the length of J and A will be set to 0, so that user can use the PushBack function for vectors.

Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]nnzThe size reserved in J and A array.
[in]locationThe location of the data.
[in]holdvalueSet to false to disable the A vector. Only save the pattern of the matrix.
[in]setzeroIf set to true, _i_vec[0] and the length of J and A is set to 0 so you can use the push_back function.
[in]iscsrSet to true for csr matrix, false for csc matrix.
Returns
Return error message.

◆ Setup() [7/7]

template<typename T >
template int pargemslr::CsrMatrixClass< T >::Setup ( int  nrows,
int  ncols,
int  nnz,
int  location,
bool  setzero 
)

Free the current matrix, and allocate memory to create a new matrix. The default location of data is the host memory.
If setzero is set to true, _i_vec[0] and the length of J and A will be set to 0, so that user can use the PushBack function for vectors.

Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]nnzThe size reserved in J and A array.
[in]locationThe location of the data.
[in]setzeroIf set to true, _i_vec[0] and the length of J and A is set to 0 so you can use the push_back function.
Returns
Return error message.

◆ SetupVectorPtrStr()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::SetupVectorPtrStr ( SequentialVectorClass< T > &  vec)

Update the structure of a vector to have same row permutation.

Parameters
[out]vecThe target vector.
Returns
Return error message.

◆ SortRow()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::SortRow ( )

Sort the column/row indices in each row/column.

Returns
Return error message.

◆ SubMatrix() [1/3]

template<typename T >
int pargemslr::CsrMatrixClass< T >::SubMatrix ( int  row_start,
int  col_start,
int  num_rows,
int  num_cols,
int  location,
CsrMatrixClass< T > &  csrmat_out 
)

Copy data to extract a submatrix of the current matrix.

Parameters
[in]row_startThe starting row of that submatrix.
[in]col_startThe starting column of that submatrix.
[in]num_rowsThe number of rows.
[in]num_colsThe number fo cols.
[in]locationThe location of the submatrix.
[out]csrmat_outThe submatrix.
Returns
Return error message.

◆ SubMatrix() [2/3]

template<typename T >
int pargemslr::CsrMatrixClass< T >::SubMatrix ( vector_int row_starts,
vector_int col_starts,
vector_int row_ends,
vector_int col_ends,
int  location,
CsrMatrixClass< T > &  csrmat_out 
)

Copy data to extract a submatrix of the current matrix.
Each interval [ rows_start[i], rows_end[i] ) gives an interval of rows to keep; interval [ cols_start[i], cols_end[i] ) gives an interval of columns to keep.

Parameters
[in]row_startsThe starting rows of the submatrix.
[in]col_startsThe starting columns of the submatrix.
[in]rows_endsThe ending rows(excluse) of each row interval.
[in]cols_endsThe ending columns(excluse) of each column interval.
[in]locationThe location of the submatrix.
[out]csrmat_outThe submatrix.
Returns
Return error message.

◆ SubMatrix() [3/3]

template<typename T >
int pargemslr::CsrMatrixClass< T >::SubMatrix ( vector_int rows,
vector_int cols,
int  location,
CsrMatrixClass< T > &  csrmat_out 
)

Copy data to extract a submatrix of the current matrix. Rows in the array rows and columns in the array cols will be kept in the order they appear in rows and cols. Can be used to apply permutation.

Parameters
[in]rowsRows to keep.
[in]colsCols to keep.
[in]locationThe location of the submatrix.
[out]csrmat_outThe submatrix.
Returns
Return error message.

◆ SubMatrixNoPerm()

template<typename T >
template int pargemslr::CsrMatrixClass< T >::SubMatrixNoPerm ( vector_int rows,
vector_int cols,
vector_int row_perm,
vector_int col_perm,
bool  complement,
int  location,
CsrMatrixClass< T > &  csrmat_out 
)

Copy data to extract a submatrix of the current matrix.
if complement == false, Nonzero rows and columns in the array rows and cols are kept in their original order.
if complement == true, Zero rows and columns in the array rows and cols are kept in their original order.

Parameters
[in]rowsRows to keep. Vector of length _nrows.
[in]colsCols to keep. Vector of length _nrows.
[out]row_permThe row permutation vector from the original matrix to the new matrix.
[out]col_permThe col permutation vector from the original matrix to the new matrix.
[in]complementKeep zero (true) or nonzero (false) in vector rows and cols.
[in]locationThe location of the submatrix.
[out]csrmat_outThe submatrix.
Returns
Return error message.

◆ Transpose() [1/2]

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

Transpost the current matrix.

Returns
Return error message.

◆ Transpose() [2/2]

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

Transpost the current matrix.

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

Member Data Documentation

◆ _diagonal_shift

template<typename T >
std::conditional<PargemslrIsDoublePrecision<T>::value, double, float>::type pargemslr::CsrMatrixClass< 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.

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