ParGeMSLR
|
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... | |
Class of CSR/CSC matrices.
template pargemslr::CsrMatrixClass< T >::CsrMatrixClass | ( | ) |
The constructor of CsrMatrixClass. The default memory location is the host memory.
pargemslr::CsrMatrixClass< T >::CsrMatrixClass | ( | const CsrMatrixClass< T > & | mat | ) |
The copy constructor of CsrMatrixClass.
[in] | mat | The target matrix. |
pargemslr::CsrMatrixClass< T >::CsrMatrixClass | ( | CsrMatrixClass< T > && | mat | ) |
The move constructor of CsrMatrixClass.
[in] | mat | The target matrix. |
|
virtual |
The destructor of CsrMatrixClass. Simply a call to the free function.
|
virtual |
Free the current matrix.
template int pargemslr::CsrMatrixClass< T >::Convert | ( | bool | csr | ) |
Convert between csr and csc.
[in] | csr | Set to true to use csr, otherwise csc. |
|
virtual |
Create an indentity matrix.
|
virtual |
Fill the matrix with constant value.
[in] | v | The value to be filled. |
|
inline |
Get the MPI_comm.
int pargemslr::CsrMatrixClass< T >::GetComplexShift | ( | ComplexValueClass< T1 > & | diagonal_shift | ) |
Get the diagonal complex shift for some preconditioner options.
[out] | diagonal_shift | The complex shift. |
int pargemslr::CsrMatrixClass< T >::GetComplexShift | ( | T1 & | diagonal_shift | ) |
Get the diagonal complex shift for some preconditioner options.
[out] | diagonal_shift | The complex shift. |
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.
[out] | comp_indices | 2D vector, comp_indices[i] is the index of the elements in component i. |
[in,out] | ncomps | On 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). |
template complexd * pargemslr::CsrMatrixClass< T >::GetData | ( | ) | const |
Get the data pointer of the matrix.
|
virtual |
Get the data location of the matrix.
template SequentialVectorClass< complexd > & pargemslr::CsrMatrixClass< T >::GetDataVector | ( | ) |
Get the reference to the data vector.
template const SequentialVectorClass< complexd > & pargemslr::CsrMatrixClass< T >::GetDataVector | ( | ) | const |
Get the reference to the data vector.
template int pargemslr::CsrMatrixClass< T >::GetDiagScale | ( | T & | scale | ) |
Scale the current csr matrix.
[in] | scale | The scale. |
template int * pargemslr::CsrMatrixClass< T >::GetI | ( | ) | const |
Get the I pointer of the matrix.
template IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetIVector | ( | ) |
Get the reference to the I vector.
template const IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetIVector | ( | ) | const |
Get the reference to the I vector.
template int * pargemslr::CsrMatrixClass< T >::GetJ | ( | ) | const |
Get the J pointer of the matrix.
template IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetJVector | ( | ) |
Get the reference to the J vector.
template const IntVectorClass< int > & pargemslr::CsrMatrixClass< T >::GetJVector | ( | ) | const |
Get the reference to the J vector.
|
inline |
Get comm, np, and myid. Get the sequential comm.
[in] | np | The number of processors. |
[in] | myid | The local MPI rank number. |
[in] | comm | The MPI_Comm. |
|
virtual |
Get the local number of columns of the matrix.
|
virtual |
Get the number of nonzeros in this matrix.
|
virtual |
Get the local number of rows of the matrix.
template int pargemslr::CsrMatrixClass< T >::Helmholtz | ( | int | n, |
T | 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.
[in] | n | Problem size. |
[in] | w | The w value. |
template bool & pargemslr::CsrMatrixClass< T >::IsCsr | ( | ) |
Csr or csc.
template const bool & pargemslr::CsrMatrixClass< T >::IsCsr | ( | ) | const |
Csr or csc.
template bool pargemslr::CsrMatrixClass< T >::IsHoldingData | ( | ) | const |
Check if the matrix holding its own data.
template bool & pargemslr::CsrMatrixClass< T >::IsRowSorted | ( | ) |
Tell if the column/row indices in each row/column is sorted.
template const bool & pargemslr::CsrMatrixClass< T >::IsRowSorted | ( | ) | const |
Tell if the column/row indices in each row/column is sorted.
template int pargemslr::CsrMatrixClass< T >::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.
[in] | nx | |
[in] | ny | |
[in] | nz | Problem size on each dimension |
[in] | alphax | |
[in] | alphay | |
[in] | alphaz | Alpha value on each direction |
[in] | shift | The shift value on diagonal |
[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::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.
[in] | alpha | The alpha value. |
[in] | A | The first csr matrix. |
[in] | transa | Whether we transpose A or not. |
[in] | beta | The beta value. |
[in] | B | The second csr matrix. |
[in] | transb | Whether we transpose B or not. |
|
virtual |
In place Matrix-Vector product ==> z := alpha*A*x + beta*y, or z := alpha*A'*x + beta*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 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. |
|
virtual |
Move the data to another memory location.
[in] | location | The location move to. |
CsrMatrixClass< T > & pargemslr::CsrMatrixClass< T >::operator= | ( | const CsrMatrixClass< T > & | mat | ) |
CsrMatrixClass< T > & pargemslr::CsrMatrixClass< T >::operator= | ( | CsrMatrixClass< T > && | mat | ) |
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.
[in] | perm | The permutation vector. Set to NULL to plot(A), otherwise plot(A(perm,perm)). |
[in] | conditiona | First condition. |
[in] | conditionb | Secend condition, only plot when conditiona == conditionb. |
[in] | width | The plot width. |
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.
[in] | rperm | The row permutation vector. |
[in] | cperm | The column permutation vector. Set either rperm or cperm to NULL to spy(A), otherwise spy(A(rperm,cperm)). |
[in] | conditiona | First condition. |
[in] | conditionb | Secend condition, only spy when conditiona == conditionb. |
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.
[in] | datafilename | The filename of the temp file holding the data. |
[in] | rperm | The row permutation vector. |
[in] | cperm | The column permutation vector. Set either rperm or cperm to NULL to spy(A), otherwise spy(A(rperm,cperm)). |
[in] | conditiona | First condition. |
[in] | conditionb | Secend condition, only plot when conditiona == conditionb. |
template int pargemslr::CsrMatrixClass< T >::PushBack | ( | int | col, |
T | v | ||
) |
Insert value at the end of the J and A vector, user need to update the I vector.
Currently only for host memory.
[in] | col | The col/row to be inserted into the csr/csc. |
[in] | v | The value to be inserted. |
template int pargemslr::CsrMatrixClass< T >::ReadFromMMFile | ( | const char * | matfile, |
int | idxin | ||
) |
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. |
|
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::CsrMatrixClass< T >::SetNumNonzeros | ( | ) |
Set the number of nonzeros in this matrix.
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.
[in] | nrows | The number of rows. |
[in] | ncols | The number of cols. |
[in] | nnz | The size reserved in J and A array. |
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.
[in] | nrows | The number of rows. |
[in] | ncols | The number of cols. |
[in] | nnz | The size reserved in J and A array. |
[in] | holdvalue | Set to false to disable the A vector. Only save the pattern of the matrix. |
[in] | setzero | If 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. |
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.
[in] | nrows | The number of rows. |
[in] | ncols | The number of cols. |
[in] | nnz | The size reserved in J and A array. |
[in] | setzero | If set to true, I[0] and the length of J and A are set to 0 so you can use the PushBack function. |
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.
[in] | nrows | The number of rows. |
[in] | ncols | The number of cols. |
[in] | nnz | The size reserved in J and A array. |
[in] | location | The location of the data. |
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.
[in] | nrows | The number of rows. |
[in] | ncols | The number of cols. |
[in] | nnz | The size reserved in J and A array. |
[in] | location | The location of the data. |
[in] | holdvalue | Set to false to disable the A vector. Only save the pattern of the matrix. |
[in] | setzero | If 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. |
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.
[in] | nrows | The number of rows. |
[in] | ncols | The number of cols. |
[in] | nnz | The size reserved in J and A array. |
[in] | location | The location of the data. |
[in] | holdvalue | Set to false to disable the A vector. Only save the pattern of the matrix. |
[in] | setzero | If 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] | iscsr | Set to true for csr matrix, false for csc matrix. |
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.
[in] | nrows | The number of rows. |
[in] | ncols | The number of cols. |
[in] | nnz | The size reserved in J and A array. |
[in] | location | The location of the data. |
[in] | setzero | If 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. |
template int pargemslr::CsrMatrixClass< T >::SetupVectorPtrStr | ( | SequentialVectorClass< T > & | vec | ) |
Update the structure of a vector to have same row permutation.
[out] | vec | The target vector. |
template int pargemslr::CsrMatrixClass< T >::SortRow | ( | ) |
Sort the column/row indices in each row/column.
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.
[in] | row_start | The starting row of that submatrix. |
[in] | col_start | The starting column of that submatrix. |
[in] | num_rows | The number of rows. |
[in] | num_cols | The number fo cols. |
[in] | location | The location of the submatrix. |
[out] | csrmat_out | The submatrix. |
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.
[in] | row_starts | The starting rows of the submatrix. |
[in] | col_starts | The starting columns of the submatrix. |
[in] | rows_ends | The ending rows(excluse) of each row interval. |
[in] | cols_ends | The ending columns(excluse) of each column interval. |
[in] | location | The location of the submatrix. |
[out] | csrmat_out | The submatrix. |
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.
[in] | rows | Rows to keep. |
[in] | cols | Cols to keep. |
[in] | location | The location of the submatrix. |
[out] | csrmat_out | The submatrix. |
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.
[in] | rows | Rows to keep. Vector of length _nrows. |
[in] | cols | Cols to keep. Vector of length _nrows. |
[out] | row_perm | The row permutation vector from the original matrix to the new matrix. |
[out] | col_perm | The col permutation vector from the original matrix to the new matrix. |
[in] | complement | Keep zero (true) or nonzero (false) in vector rows and cols. |
[in] | location | The location of the submatrix. |
[out] | csrmat_out | The submatrix. |
template int pargemslr::CsrMatrixClass< T >::Transpose | ( | ) |
Transpost the current matrix.
int pargemslr::CsrMatrixClass< T >::Transpose | ( | CsrMatrixClass< T > & | AT | ) |
Transpost the current matrix.
[out] | AT | The transpose matrix. |
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.