ParGeMSLR
Functions
matrixops.hpp File Reference

Matrix operations. More...

#include <vector>
#include "../vectors/vector.hpp"
#include "../utils/mmio.hpp"
#include "matrix.hpp"
#include "csr_matrix.hpp"
#include "parallel_csr_matrix.hpp"
#include "coo_matrix.hpp"
#include "dense_matrix.hpp"

Go to the source code of this file.

Functions

PrecisionEnum pargemslr::GetMatrixPrecision (const MatrixClass< int > &mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPrecision (const MatrixClass< long int > &mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPrecision (const MatrixClass< float > &mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPrecision (const MatrixClass< double > &mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPrecision (const MatrixClass< complexs > &mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPrecision (const MatrixClass< complexd > &mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPPrecision (const MatrixClass< int > *mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPPrecision (const MatrixClass< long int > *mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPPrecision (const MatrixClass< float > *mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPPrecision (const MatrixClass< double > *mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPPrecision (const MatrixClass< complexs > *mat)
 Get the precision of a matrix. More...
 
PrecisionEnum pargemslr::GetMatrixPPrecision (const MatrixClass< complexd > *mat)
 Get the precision of a matrix. More...
 
template<typename T >
int pargemslr::DenseMatrixPMatVecTemplate (char trans, int nrows, int ncols, const T &alpha, const T *aa, int ldim, const T *x, const T &beta, T *y)
 In place dense complex Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::DenseMatrixMatVec (const DenseMatrixClass< float > &A, char trans, const float &alpha, const VectorClass< float > &x, const float &beta, VectorClass< float > &y)
 In place dense float Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::DenseMatrixMatVec (const DenseMatrixClass< double > &A, char trans, const double &alpha, const VectorClass< double > &x, const double &beta, VectorClass< double > &y)
 In place dense double Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::DenseMatrixMatVec (const DenseMatrixClass< complexs > &A, char trans, const complexs &alpha, const VectorClass< complexs > &x, const complexs &beta, VectorClass< complexs > &y)
 In place dense single complex Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::DenseMatrixMatVec (const DenseMatrixClass< complexd > &A, char trans, const complexd &alpha, const VectorClass< complexd > &x, const complexd &beta, VectorClass< complexd > &y)
 In place dense double complex Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::DenseMatrixInvertHost (DenseMatrixClass< float > &A)
 Compute the inverese of a general float matrix. More...
 
int pargemslr::DenseMatrixInvertHost (DenseMatrixClass< double > &A)
 Compute the inverese of a general double matrix. More...
 
int pargemslr::DenseMatrixInvertHost (DenseMatrixClass< complexs > &A)
 Compute the inverese of a general single complex matrix. More...
 
int pargemslr::DenseMatrixInvertHost (DenseMatrixClass< complexd > &A)
 Compute the inverese of a general double complex matrix. More...
 
int pargemslr::DenseMatrixInvertUpperTriangularHost (DenseMatrixClass< float > &A)
 Compute the inverese of a upper triangular float matrix. More...
 
int pargemslr::DenseMatrixInvertUpperTriangularHost (DenseMatrixClass< double > &A)
 Compute the inverese of a upper triangular double matrix. More...
 
int pargemslr::DenseMatrixInvertUpperTriangularHost (DenseMatrixClass< complexs > &A)
 Compute the inverese of a upper triangular single complex matrix. More...
 
int pargemslr::DenseMatrixInvertUpperTriangularHost (DenseMatrixClass< complexd > &A)
 Compute the inverese of a upper triangular double complex matrix. More...
 
int pargemslr::DenseMatrixQRDecompositionHost (DenseMatrixClass< float > &A, DenseMatrixClass< float > &Q)
 Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn. More...
 
int pargemslr::DenseMatrixQRDecompositionHost (DenseMatrixClass< double > &A, DenseMatrixClass< double > &Q)
 Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn. More...
 
int pargemslr::DenseMatrixQRDecompositionHost (DenseMatrixClass< complexs > &A, DenseMatrixClass< complexs > &Q)
 Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn. More...
 
int pargemslr::DenseMatrixQRDecompositionHost (DenseMatrixClass< complexd > &A, DenseMatrixClass< complexd > &Q)
 Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< float > &A, DenseMatrixClass< float > &Q)
 Transform a dense float matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< float > &A, int start, int end, DenseMatrixClass< float > &Q)
 Transform a dense float matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< double > &A, DenseMatrixClass< double > &Q)
 Transform a dense double matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< double > &A, int start, int end, DenseMatrixClass< double > &Q)
 Transform a dense double matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< complexs > &A, DenseMatrixClass< complexs > &Q)
 Transform a dense single complex matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< complexs > &A, int start, int end, DenseMatrixClass< complexs > &Q)
 Transform a dense single complex matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< complexd > &A, DenseMatrixClass< complexd > &Q)
 Transform a dense double complex matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixHessDecompositionHost (DenseMatrixClass< complexd > &A, int start, int end, DenseMatrixClass< complexd > &Q)
 Transform a dense double complex matrix A into hessenberg matrix Q^HAQ = Hess. More...
 
int pargemslr::DenseMatrixRealHessSchurDecompositionHost (DenseMatrixClass< float > &A, DenseMatrixClass< float > &Q, vector_seq_float &wr, vector_seq_float &wi)
 Compute the Schur decomposition of a float dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixRealHessSchurDecompositionHost (DenseMatrixClass< float > &A, int start, int end, DenseMatrixClass< float > &Q, vector_seq_float &wr, vector_seq_float &wi)
 Compute the Schur decomposition of a float dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixRealHessSchurDecompositionHost (DenseMatrixClass< double > &A, DenseMatrixClass< double > &Q, vector_seq_double &wr, vector_seq_double &wi)
 Compute the Schur decomposition of a double dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixRealHessSchurDecompositionHost (DenseMatrixClass< double > &A, int start, int end, DenseMatrixClass< double > &Q, vector_seq_double &wr, vector_seq_double &wi)
 Compute the Schur decomposition of a double dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixComplexHessSchurDecompositionHost (DenseMatrixClass< complexs > &A, DenseMatrixClass< complexs > &Q, vector_seq_complexs &w)
 Compute the Schur decomposition of a single complex dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixComplexHessSchurDecompositionHost (DenseMatrixClass< complexs > &A, int start, int end, DenseMatrixClass< complexs > &Q, vector_seq_complexs &w)
 Compute the Schur decomposition of a single complex dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixComplexHessSchurDecompositionHost (DenseMatrixClass< complexd > &A, DenseMatrixClass< complexd > &Q, vector_seq_complexd &w)
 Compute the Schur decomposition of a double complex dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixComplexHessSchurDecompositionHost (DenseMatrixClass< complexd > &A, int start, int end, DenseMatrixClass< complexd > &Q, vector_seq_complexd &w)
 Compute the Schur decomposition of a double complex dense hessenberg matrix A = QUQ^H. More...
 
int pargemslr::DenseMatrixRealOrderSchur (DenseMatrixClass< float > &A, DenseMatrixClass< float > &Q, vector_seq_float &wr, vector_seq_float &wi, vector_int &select)
 Reorder the float schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixRealOrderSchur (DenseMatrixClass< double > &A, DenseMatrixClass< double > &Q, vector_seq_double &wr, vector_seq_double &wi, vector_int &select)
 Reorder the double schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixComplexOrderSchur (DenseMatrixClass< complexs > &A, DenseMatrixClass< complexs > &Q, vector_seq_complexs &w, vector_int &select)
 Reorder the single complex schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixComplexOrderSchur (DenseMatrixClass< complexd > &A, DenseMatrixClass< complexd > &Q, vector_seq_complexd &w, vector_int &select)
 Reorder the double complex schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixRealOrderSchurClusters (DenseMatrixClass< float > &A, DenseMatrixClass< float > &Q, vector_seq_float &wr, vector_seq_float &wi, vector_int &clusters)
 Reorder the float schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixRealOrderSchurClusters (DenseMatrixClass< double > &A, DenseMatrixClass< double > &Q, vector_seq_double &wr, vector_seq_double &wi, vector_int &clusters)
 Reorder the double schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixComplexOrderSchurClusters (DenseMatrixClass< complexs > &A, DenseMatrixClass< complexs > &Q, vector_seq_complexs &w, vector_int &clusters)
 Reorder the single complex schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixComplexOrderSchurClusters (DenseMatrixClass< complexd > &A, DenseMatrixClass< complexd > &Q, vector_seq_complexd &w, vector_int &clusters)
 Reorder the double complex schur decomposition, puch selected eigenvalues to the leading part. More...
 
int pargemslr::DenseMatrixRealHessEigenDecompositionHost (DenseMatrixClass< float > &A, DenseMatrixClass< float > &Q, vector_seq_float &wr, vector_seq_float &wi)
 Compute the Eigen decomposition of float dense hessenberg matrix A, AQ = QD. A is not modified. More...
 
int pargemslr::DenseMatrixRealHessEigenDecompositionHost (DenseMatrixClass< double > &A, DenseMatrixClass< double > &Q, vector_seq_double &wr, vector_seq_double &wi)
 Compute the Eigen decomposition of double dense hessenberg matrix A, AQ = QD. A is not modified. More...
 
int pargemslr::DenseMatrixComplexHessEigenDecompositionHost (DenseMatrixClass< complexs > &A, DenseMatrixClass< complexs > &Q, vector_seq_complexs &w)
 Compute the Eigen decomposition of single complex dense hessenberg matrix A, AQ = QD. A is not modified. More...
 
int pargemslr::DenseMatrixComplexHessEigenDecompositionHost (DenseMatrixClass< complexd > &A, DenseMatrixClass< complexd > &Q, vector_seq_complexd &w)
 Compute the Eigen decomposition of double complex dense hessenberg matrix A, AQ = QD. A is not modified. More...
 
template<typename T >
int pargemslr::DenseMatrixMatMatTemplate (const T &alpha, const DenseMatrixClass< T > &A, char transa, const DenseMatrixClass< T > &B, char transb, const T &beta, DenseMatrixClass< T > &C)
 Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C. More...
 
int pargemslr::DenseMatrixMatMat (const float &alpha, const DenseMatrixClass< float > &A, char transa, const DenseMatrixClass< float > &B, char transb, const float &beta, DenseMatrixClass< float > &C)
 Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C. More...
 
int pargemslr::DenseMatrixMatMat (const double &alpha, const DenseMatrixClass< double > &A, char transa, const DenseMatrixClass< double > &B, char transb, const double &beta, DenseMatrixClass< double > &C)
 Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C. More...
 
int pargemslr::DenseMatrixMatMat (const complexs &alpha, const DenseMatrixClass< complexs > &A, char transa, const DenseMatrixClass< complexs > &B, char transb, const complexs &beta, DenseMatrixClass< complexs > &C)
 Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C. More...
 
int pargemslr::DenseMatrixMatMat (const complexd &alpha, const DenseMatrixClass< complexd > &A, char transa, const DenseMatrixClass< complexd > &B, char transb, const complexd &beta, DenseMatrixClass< complexd > &C)
 Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C. More...
 
template<typename T >
int pargemslr::CsrMatrixPMatVecHostTemplate (const int *ia, const int *ja, const T *aa, int nrows, int ncols, char trans, const T &alpha, const T *x, const T &beta, T *y)
 In place float csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixPMatVecHost (const int *ia, const int *ja, const float *aa, int nrows, int ncols, char trans, const float &alpha, const float *x, const float &beta, float *y)
 In place csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixPMatVecHost (const int *ia, const int *ja, const double *aa, int nrows, int ncols, char trans, const double &alpha, const double *x, const double &beta, double *y)
 In place double csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixPMatVecHost (const int *ia, const int *ja, const complexs *aa, int nrows, int ncols, char trans, const complexs &alpha, const complexs *x, const complexs &beta, complexs *y)
 In place single complex csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixPMatVecHost (const int *ia, const int *ja, const complexd *aa, int nrows, int ncols, char trans, const complexd &alpha, const complexd *x, const complexd &beta, complexd *y)
 In place double complex csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixMatVec (const CsrMatrixClass< float > &A, char trans, const float &alpha, const VectorClass< float > &x, const float &beta, VectorClass< float > &y)
 In place float csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixMatVec (const CsrMatrixClass< double > &A, char trans, const double &alpha, const VectorClass< double > &x, const double &beta, VectorClass< double > &y)
 In place double csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixMatVec (const CsrMatrixClass< complexs > &A, char trans, const complexs &alpha, const VectorClass< complexs > &x, const complexs &beta, VectorClass< complexs > &y)
 In place single complex csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
int pargemslr::CsrMatrixMatVec (const CsrMatrixClass< complexd > &A, char trans, const complexd &alpha, const VectorClass< complexd > &x, const complexd &beta, VectorClass< complexd > &y)
 In place double complex csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y. More...
 
template<int INIDX, int OUTIDX, typename T >
int pargemslr::CsrMatrixP2CscMatrixPHost (int nrows, int ncols, bool copy_data, T *ai, int *ji, int *ii, T *ao, int *jo, int *io)
 Converts a csr matrix to csc matrix. More...
 
template<int INIDX, int OUTIDX, typename T >
int pargemslr::CooMatrixP2CsrMatrixPHost (int nrows, int ncols, int nnz, T *ai, int *ji, int *ii, T *ao, int *jo, int *io)
 Converts a coo matrix to csr matrix. More...
 
template<typename T >
int pargemslr::CsrMatrixTransposeHost (CsrMatrixClass< T > &A, CsrMatrixClass< T > &AT)
 Compute the transpose of a csr/csc matrix. More...
 
template<typename T >
int pargemslr::CsrMatrixAddHost (CsrMatrixClass< T > &A, CsrMatrixClass< T > &B, CsrMatrixClass< T > &C)
 Add sparse matrix in CSR format, C= A + B. More...
 
template<typename T >
int pargemslr::ParallelCsrMatrixTransposeHost (ParallelCsrMatrixClass< T > &A, ParallelCsrMatrixClass< T > &AT)
 Compute the transpose of a parallel csr matrix. More...
 
template<typename T >
int pargemslr::ParallelCsrMatrixAddHost (ParallelCsrMatrixClass< T > &A, ParallelCsrMatrixClass< T > &B, ParallelCsrMatrixClass< T > &C)
 Add sparse matrix in Parallel CSR format, C= A + B. More...
 
template<typename T >
int pargemslr::CsrMatrixMetisKwayHost (CsrMatrixClass< T > &A, int &num_dom, IntVectorClass< int > &map, bool vertexsep, IntVectorClass< int > &sep, int &edgecut, IntVectorClass< int > &perm, IntVectorClass< int > &dom_ptr)
 Kway metis partition. For the vertex seperator, we only supports 2-way vertex seperator. More...
 
template<typename T >
int pargemslr::CsrMatrixMaxMatchingHost (CsrMatrixClass< T > &A, int &nmatch, IntVectorClass< int > &match_row, IntVectorClass< int > &match_col)
 Use Hopcroft–Karp algorithm to find a maximal matching. More...
 
bool pargemslr::CsrMatrixMaxMatchingBfsHost (int nrows, int ncols, int *A_i, int *A_j, IntVectorClass< int > &dist, IntVectorClass< int > &match_row, IntVectorClass< int > &match_col)
 The BFS for maximal matching algorithm. More...
 
bool pargemslr::CsrMatrixMaxMatchingDfsHost (int nrows, int ncols, int *A_i, int *A_j, IntVectorClass< int > &dist, IntVectorClass< int > &match_row, IntVectorClass< int > &match_col, int i)
 The BFS for maximal matching algorithm. More...
 
template<typename T >
int pargemslr::CsrMatrixSortRow (CsrMatrixClass< T > &A)
 Sort columns inside each row of a real csr matrix by ascending order. More...
 
template<typename T >
int pargemslr::CsrSubMatrixAmdHost (CsrMatrixClass< T > &A, vector_int &rowscols, vector_int &perm)
 This function computes the AMD ordering of a submatrix A[rowscols,rowscols] + trans(A[rowscols,rowscols]). More...
 
template<typename T >
int pargemslr::CsrMatrixAmdHost (CsrMatrixClass< T > &A, vector_int &perm)
 This function computes the Amd ordering of a CSR matrix A+trans(A). More...
 
template<typename T >
int pargemslr::CsrSubMatrixNdHost (CsrMatrixClass< T > &A, vector_int &rowscols, vector_int &perm)
 This function computes the ND ordering of a submatrix A[rowscols,rowscols]. More...
 
template<typename T >
int pargemslr::CsrMatrixNdHost (CsrMatrixClass< T > &A, vector_int &perm)
 This function computes the ND ordering of a CSR matrix A. More...
 
template<typename T >
int pargemslr::CsrSubMatrixRcmHost (CsrMatrixClass< T > &A, vector_int &rowscols, vector_int &perm)
 This function computes the RCM ordering of a submatrix A[rowscols,rowscols] + trans(A[rowscols,rowscols]). More...
 
template<typename T >
int pargemslr::CsrMatrixRcmHost (CsrMatrixClass< T > &A, vector_int &perm)
 This function computes the RCM ordering of a CSR matrix A+trans(A). More...
 
template<typename T >
int pargemslr::CsrMatrixRcmRootHost (CsrMatrixClass< T > &G, vector_int &marker, int &root)
 This function find the next root for the RCM ordering of graph G. More...
 
template<typename T >
int pargemslr::CsrMatrixRcmNumberingHost (CsrMatrixClass< T > &G, int root, vector_int &marker, vector_int &perm, int &current_num)
 This function computes the RCM numbering for a connect component start from root. More...
 
template<typename T >
int pargemslr::CsrMatrixRcmPerphnHost (CsrMatrixClass< T > &G, int &root, vector_int &marker)
 This function find an end of the pseudo-peripheral starting from the root. More...
 
template<typename T >
int pargemslr::CsrMatrixRcmBfsHost (CsrMatrixClass< T > &G, int root, vector_int &marker, std::vector< std::vector< int > > &level)
 This function builds the level structure starting from node root. More...
 
int pargemslr::CsrMatrixRcmClearLevelHost (std::vector< std::vector< int > > &level)
 This function frees the level structure. More...
 
int pargemslr::CsrMatrixRcmSwapHost (vector_int &perm, int a, int b)
 This function swap two element in the perm array. More...
 
int pargemslr::CsrMatrixRcmReverseHost (vector_int &perm, int start, int end)
 This function reverse an array between [start,end]. More...
 
int pargemslr::ParmetisKwayHost (vector_long &vtxdist, vector_long &xadj, vector_long &adjncy, long int &num_dom, vector_long &map, parallel_log &parlog)
 Call ParMETIS to partition A based on the adjancency graph of A. More...
 
int pargemslr::ParmetisNodeND (vector_long &vtxdist, vector_long &xadj, vector_long &adjncy, long int &num_dom, vector_long &map, parallel_log &parlog)
 Partition with log(p) levels nd, p is the number of processors. More...
 
template<class VectorType , class MatrixType , typename DataType , typename RealDataType >
int pargemslr::PargemslrSubSpaceIteration (MatrixType &A, int k, int its, DenseMatrixClass< DataType > &V, DenseMatrixClass< DataType > &H, RealDataType tol)
 The subspace iteration procedure on a matrix. Note that we do not check for convergence. Apply a fixed number of iterations. More...
 
template<class VectorType , class MatrixType , typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiNoRestart (MatrixType &A, int mstart, int msteps, DenseMatrixClass< DataType > &V, DenseMatrixClass< DataType > &H, RealDataType tol_orth, RealDataType tol_reorth)
 The standard Arnoldi procedure on a matrix. More...
 
template<class VectorType , class MatrixType , typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiThickRestartNoLock (MatrixType &A, int msteps, int maxits, int rank, int rank2, RealDataType truncate, RealDataType tr_fact, RealDataType tol_eig, RealDataType eig_target_mag, RealDataType eig_truncate, RealDataType(*weight)(ComplexValueClass< RealDataType >), DenseMatrixClass< DataType > &V, DenseMatrixClass< DataType > &H, RealDataType tol_orth, RealDataType tol_reorth)
 The thick-restart Arnoldi without locking. This version should be used when convergence tolerance is high. More...
 
template<typename T >
int pargemslr::PargemslrArnoldiThickRestartChooseEigenValuesReal (DenseMatrixClass< T > &H, DenseMatrixClass< T > &Q, T h_last, T(*weight)(ComplexValueClass< T >), T truncate, int &ncov, int &nicov, int &nsatis, T tol_eig, T eig_target_mag, T eig_truncate, bool &cut, SequentialVectorClass< T > &wr, SequentialVectorClass< T > &wi, vector_int &icov, vector_int &iicov, vector_int &isatis, SequentialVectorClass< T > &dcov, SequentialVectorClass< T > &dicov, SequentialVectorClass< T > &dsatis)
 Pick convergenced eigenvalues and unconvergenced eigenvalues. More...
 
template<typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiThickRestartChooseEigenValuesComplex (DenseMatrixClass< DataType > &H, DenseMatrixClass< DataType > &Q, DataType h_last, RealDataType(*weight)(DataType), RealDataType truncate, int &ncov, int &nicov, int &nsatis, RealDataType tol_eig, RealDataType eig_target_mag, RealDataType eig_truncate, bool &cut, SequentialVectorClass< DataType > &w, vector_int &icov, vector_int &iicov, vector_int &isatis, SequentialVectorClass< RealDataType > &dcov, SequentialVectorClass< RealDataType > &dicov, SequentialVectorClass< RealDataType > &dsatis)
 Pick convergenced eigenvalues and unconvergenced eigenvalues. More...
 
template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiThickRestartBuildThickRestartNewVector (DenseMatrixClass< DataType > &V, DenseMatrixClass< DataType > &H, int m, RealDataType tol_orth, RealDataType tol_reorth, VectorType &v)
 The current V is a invarient subspace, pick a new vector to restart. More...
 
template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrCgs2 (VectorType &w, DenseMatrixClass< DataType > &V, DenseMatrixClass< DataType > &H, RealDataType &t, int k, RealDataType tol_orth)
 The classical Gram-Schmidt with re-orthgonal (always one step of re-orth). More...
 
template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrMgs (VectorType &w, DenseMatrixClass< DataType > &V, DenseMatrixClass< DataType > &H, RealDataType &t, int k, RealDataType tol_orth, RealDataType tol_reorth)
 The modified Gram-Schmidt with re-orthgonal. More...
 
template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrOrthogonal (VectorType &w, DenseMatrixClass< DataType > &V, RealDataType &t, int k, RealDataType tol_orth, RealDataType tol_reorth)
 The modified Gram-Schmidt with re-orthgonal without H matrix. More...
 
template<typename T >
int pargemslr::DenseMatrixPlotHost (DenseMatrixClass< T > &A, int conditiona, int conditionb, int width)
 Plot the dense matrix to the terminal output. Function for testing purpose. More...
 
template<typename T >
int pargemslr::CsrMatrixPlotHost (CsrMatrixClass< T > &A, int *perm, int conditiona, int conditionb, int width)
 Plot the csr matrix to the terminal output. Function for testing purpose. More...
 
int pargemslr::CooMatrixReadFromFile (CooMatrixClass< float > &coo, const char *matfile, int idxin, int idxout)
 Read a coo matrix from matrix marker format file. More...
 
int pargemslr::CooMatrixReadFromFile (CooMatrixClass< double > &coo, const char *matfile, int idxin, int idxout)
 Read a coo matrix from matrix marker format file. More...
 
int pargemslr::CooMatrixReadFromFile (CooMatrixClass< complexs > &coo, const char *matfile, int idxin, int idxout)
 Read a coo matrix from matrix marker format file. More...
 
int pargemslr::CooMatrixReadFromFile (CooMatrixClass< complexd > &coo, const char *matfile, int idxin, int idxout)
 Read a coo matrix from matrix marker format file. More...
 
template<typename T >
int pargemslr::SetupPermutationRKwayRecursive (CsrMatrixClass< T > &A, bool vertexsep, int clvl, int &tlvl, int num_dom, int minsep, int kmin, int kfactor, vector_int &map_v, vector_int &mapptr_v)
 Recursive KWay partition call. More...
 
template<typename T >
int pargemslr::SetupPermutationNDRecursive (CsrMatrixClass< T > &A, bool vertexsep, int clvl, int &tlvl, int minsep, std::vector< std::vector< vector_int > > &level_str)
 Recursive ND partition call. More...
 
int pargemslr::SetupPermutationNDCombineLevels (std::vector< vector_int > &level_stri, int &ndom)
 Compress a certain level of level_str from SetupPermutationNDRecursive into a give number of domains. More...
 
template<typename T >
int pargemslr::ParallelCsrMatrixSetupPermutationParallelRKway (ParallelCsrMatrixClass< T > &A, bool vertexsep, int &nlev, long int ncomp, long int minsep, long int kmin, long int kfactor, vector_int &map_v, vector_int &mapptr_v, bool bj_last)
 The main recursive Kway reordering function. More...
 
int pargemslr::SetupPermutationParallelRKwayRecursive (vector_long &vtxdist, vector_long &xadj, vector_long &adjncy, int clvl, int &tlvl, long int ncomp, long int minsep, long int kmin, long int kfactor, vector_int &map_v, vector_int &mapptr_v, bool bj_last, parallel_log &parlog)
 This is the recursive function to build the mapping information. More...
 
int pargemslr::ParallelRKwayGetSeparator (vector_long &vtxdist, vector_long &xadj, vector_long &adjncy, vector_long &vtxdist_s, vector_long &xadj_s, vector_long &adjncy_s, vector_long &map, int num_dom, vector_int &vtxsep, parallel_log &parlog)
 This function finds the edge separator based on the map information. More...
 
template<typename T >
int pargemslr::SetupPermutationParallelRKwayRecursive2 (ParallelCsrMatrixClass< T > &A, int clvl, int &tlvl, long int ncomp, long int minsep, long int kmin, long int kfactor, vector_int &map_v, vector_int &mapptr_v, bool bj_last, parallel_log &parlog)
 This is the recursive function to build the mapping information, using an approximate vetrex separator (not the minimal). More...
 
template<typename T >
int pargemslr::SetupPermutationParallelKwayVertexSep (ParallelCsrMatrixClass< T > &A, long int &ncomp, vector_int &map_v, vector_long &perm_sep, parallel_log &parlog)
 This function build the k-way partition where k is the power of 2, and provide an approximate vertex saperator. More...
 
template<typename T >
int pargemslr::SetupPermutationParallelKwayVertexSepRecursive (ParallelCsrMatrixClass< T > &A, int clvl, int tlvl, bool &succeed, vector_int &map_v, parallel_log &parlog)
 This recursie call function for building the k-way partition where k is the power of 2, and provide an approximate vertex saperator. More...
 
int pargemslr::ParallelNDGetSeparator (vector_long &vtxdist, vector_long &xadj, vector_long &adjncy, vector_long &map, long int &ndom1, long int &ndom2, long int &edge_cut, parallel_log &parlog)
 This function finds the vertex separator based on the map information. A very simple appraoch, find the edge separator and only keep one side. More...
 
template<typename T >
int pargemslr::ParallelCsrMatrixSetupPermutationParallelND (ParallelCsrMatrixClass< T > &A, bool vertexsep, int &nlev, long int minsep, vector_int &map_v, vector_int &mapptr_v)
 The main ND reordering function using ND in the ParMETIS. More...
 
template<typename T >
int pargemslr::ParallelCsrMatrixSetupIOOrder (ParallelCsrMatrixClass< T > &parcsr_in, vector_int &local_perm, int &nI, CsrMatrixClass< T > &B_mat, CsrMatrixClass< T > &E_mat, CsrMatrixClass< T > &F_mat, ParallelCsrMatrixClass< T > &C_mat, int perm_option, bool perm_c)
 Setup the local order splits the interior and exterior nodes. More...
 

Detailed Description

Matrix operations.
DenseMatrixXxx: functions for dense matrices.
CsrMatrixXxx: functions for csr/csc matrices.

Function Documentation

◆ CooMatrixP2CsrMatrixPHost()

template<int INIDX, int OUTIDX, typename T >
int pargemslr::CooMatrixP2CsrMatrixPHost ( int  nrows,
int  ncols,
int  nnz,
T *  ai,
int *  ji,
int *  ii,
T *  ao,
int *  jo,
int *  io 
)

Converts a coo matrix to csr matrix.

Parameters
[in]INIDXspecifies if COO should be 0/1 index.
[in]OUTIDXspecifies if CSR should be 0/1 index.
[in]nrowsNumber of rows.
[in]ncolsNumber of columns.
[in]nnzNumber of nonzero entries.
[in]aiData vector of the input coo matrix.
[in]jiCol indices vector of the input coo matrix.
[in]iiRow indices vector of the input coo matrix.
[out]aoData vector for the output csr matrix.
[out]joJ vector for the output csr matrix.
[out]ioI vector for the output csr matrix.
Returns
Return error message.

◆ CooMatrixReadFromFile() [1/4]

int pargemslr::CooMatrixReadFromFile ( CooMatrixClass< complexd > &  coo,
const char *  matfile,
int  idxin,
int  idxout 
)

Read a coo matrix from matrix marker format file.

Parameters
[out]cooThe return matrix.
[in]matfileThe file name.
[in]idxinThe index base of the input matrix, 0-based or 1-based.
[in]idxinThe index base of the output matrix, 0-based or 1-based.
Returns
Return error message.

◆ CooMatrixReadFromFile() [2/4]

int pargemslr::CooMatrixReadFromFile ( CooMatrixClass< complexs > &  coo,
const char *  matfile,
int  idxin,
int  idxout 
)

Read a coo matrix from matrix marker format file.

Parameters
[out]cooThe return matrix.
[in]matfileThe file name.
[in]idxinThe index base of the input matrix, 0-based or 1-based.
[in]idxinThe index base of the output matrix, 0-based or 1-based.
Returns
Return error message.

◆ CooMatrixReadFromFile() [3/4]

int pargemslr::CooMatrixReadFromFile ( CooMatrixClass< double > &  coo,
const char *  matfile,
int  idxin,
int  idxout 
)

Read a coo matrix from matrix marker format file.

Parameters
[out]cooThe return matrix.
[in]matfileThe file name.
[in]idxinThe index base of the input matrix, 0-based or 1-based.
[in]idxinThe index base of the output matrix, 0-based or 1-based.
Returns
Return error message.

◆ CooMatrixReadFromFile() [4/4]

int pargemslr::CooMatrixReadFromFile ( CooMatrixClass< float > &  coo,
const char *  matfile,
int  idxin,
int  idxout 
)

Read a coo matrix from matrix marker format file.

Parameters
[out]cooThe return matrix.
[in]matfileThe file name.
[in]idxinThe index base of the input matrix, 0-based or 1-based.
[in]idxinThe index base of the output matrix, 0-based or 1-based.
Returns
Return error message.

◆ CsrMatrixAddHost()

template<typename T >
int pargemslr::CsrMatrixAddHost ( CsrMatrixClass< T > &  A,
CsrMatrixClass< T > &  B,
CsrMatrixClass< T > &  C 
)

Add sparse matrix in CSR format, C= A + B.

Parameters
[in]AMatrix A.
[in]BMatrix B.
[out]CMatrix C = A + B.
Returns
Return error message.

◆ CsrMatrixAmdHost()

template<typename T >
int pargemslr::CsrMatrixAmdHost ( CsrMatrixClass< T > &  A,
vector_int perm 
)

This function computes the Amd ordering of A+trans(A).

Parameters
[in]AThe target CSR matrix.
[out]permThe RCM permutation.
Returns
Return error message.

◆ CsrMatrixMatVec() [1/4]

int pargemslr::CsrMatrixMatVec ( const CsrMatrixClass< complexd > &  A,
char  trans,
const complexd alpha,
const VectorClass< complexd > &  x,
const complexd beta,
VectorClass< complexd > &  y 
)

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

Parameters
[in]AThe matrix.
[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.

◆ CsrMatrixMatVec() [2/4]

int pargemslr::CsrMatrixMatVec ( const CsrMatrixClass< complexs > &  A,
char  trans,
const complexs alpha,
const VectorClass< complexs > &  x,
const complexs beta,
VectorClass< complexs > &  y 
)

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

Parameters
[in]AThe matrix.
[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.

◆ CsrMatrixMatVec() [3/4]

int pargemslr::CsrMatrixMatVec ( const CsrMatrixClass< double > &  A,
char  trans,
const double &  alpha,
const VectorClass< double > &  x,
const double &  beta,
VectorClass< double > &  y 
)

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

Parameters
[in]AThe matrix.
[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.

◆ CsrMatrixMatVec() [4/4]

int pargemslr::CsrMatrixMatVec ( const CsrMatrixClass< float > &  A,
char  trans,
const float &  alpha,
const VectorClass< float > &  x,
const float &  beta,
VectorClass< float > &  y 
)

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

Parameters
[in]AThe matrix.
[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.

◆ CsrMatrixMaxMatchingBfsHost()

bool pargemslr::CsrMatrixMaxMatchingBfsHost ( int  nrows,
int  ncols,
int *  A_i,
int *  A_j,
IntVectorClass< int > &  dist,
IntVectorClass< int > &  match_row,
IntVectorClass< int > &  match_col 
)

The BFS for maximal matching algorithm.

Note
See https://en.wikipedia.org/wiki/Hopcroft%E2%80%93Karp_algorithm for more detail.
Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]A_iThe I in CSR format.
[in]A_jThe J in CSR format.
[in]distThe distance array.
[in]match_rowThe match array for rows.
[in]match_colThe match array for cols.
Returns
Return true if Bfs find new path.

◆ CsrMatrixMaxMatchingDfsHost()

bool pargemslr::CsrMatrixMaxMatchingDfsHost ( int  nrows,
int  ncols,
int *  A_i,
int *  A_j,
IntVectorClass< int > &  dist,
IntVectorClass< int > &  match_row,
IntVectorClass< int > &  match_col,
int  i 
)

The BFS for maximal matching algorithm.

Note
See https://en.wikipedia.org/wiki/Hopcroft%E2%80%93Karp_algorithm for more detail.
Parameters
[in]nrowsThe number of rows.
[in]ncolsThe number of cols.
[in]A_iThe I in CSR format.
[in]A_jThe J in CSR format.
[in]distThe distance array.
[in]match_rowThe match array for rows.
[in]match_colThe match array for cols.
[in]iThe start node.
Returns
Return true if Dfs find new path.

◆ CsrMatrixMaxMatchingHost()

template<typename T >
int pargemslr::CsrMatrixMaxMatchingHost ( CsrMatrixClass< T > &  A,
int &  nmatch,
IntVectorClass< int > &  match_row,
IntVectorClass< int > &  match_col 
)

Use Hopcroft–Karp algorithm to find a maximal matching.

Note
See https://en.wikipedia.org/wiki/Hopcroft%E2%80%93Karp_algorithm for more detail.
Parameters
[in]AMatrix A.
[out]nmatchNumber of matched pairs.
[out]match_rowThe column map vector, match_row[i] = j if i-th row is mapping to j-th col, match_row[i] = -1 if not matched.
[out]match_colThe column map vector, match_col[i] = j if i-th col is mapping to j-th row, match_col[i] = -1 if not matched.
Returns
Return error message.

◆ CsrMatrixMetisKwayHost()

template<typename T >
int pargemslr::CsrMatrixMetisKwayHost ( CsrMatrixClass< T > &  A,
int &  num_dom,
IntVectorClass< int > &  map,
bool  vertexsep,
IntVectorClass< int > &  sep,
int &  edgecut,
IntVectorClass< int > &  perm,
IntVectorClass< int > &  dom_ptr 
)

Kway metis partition. For the vertex seperator, we only supports 2-way vertex seperator.

Note
Only for symmetrix matirx, please compute A = A+A' before calling this function
Parameters
[in]AMatrix A.
[in,out]num_domNumber of domains.
[out]mapThe map vector, map[i] is the domain number of the ith node.
[in]vertexsepSet to true to compute vertex seperator. Only work when num_dom is 2.
[out]sepThe seperator vector, sep[i] != 0 means the ith node is in the seperator.
[out]edgecutThe size of the nodes inside vertexsep.
[out]permThe permutation array.
[out]dom_ptrPointer to the start node of each domain in the new order of length num_dom+1, the last entry is the size of A.
Returns
Return METIS error message.

◆ CsrMatrixNdHost()

template<typename T >
int pargemslr::CsrMatrixNdHost ( CsrMatrixClass< T > &  A,
vector_int perm 
)

This function computes the ND ordering of a CSR matrix A.

Parameters
[in]AThe target CSR matrix.
[out]permThe RCM permutation.
Returns
Return error message.

◆ CsrMatrixP2CscMatrixPHost()

template<int INIDX, int OUTIDX, typename T >
int pargemslr::CsrMatrixP2CscMatrixPHost ( int  nrows,
int  ncols,
bool  copy_data,
T *  ai,
int *  ji,
int *  ii,
T *  ao,
int *  jo,
int *  io 
)

Converts a csr matrix to csc matrix.

Note
The OpenMP version requires fixed schedule for correct result.
Parameters
[in]INIDXspecifies if COO should be 0/1 index.
[in]OUTIDXspecifies if CSR should be 0/1 index.
[in]nrowNumber of rows.
[in]ncolNumber of columns.
[in]copy_dataWhether we copy the data vector.
[in]aiData vector of the input csr matrix.
[in]jiJ vector of the input csr matrix.
[in]iiI vector of the input csr matrix.
[out]aoData vector for the output csc matrix.
[out]joJ vector for the output csc matrix.
[out]ioI vector for the output csc matrix.
Returns
Return error message.

◆ CsrMatrixPlotHost()

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

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

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

◆ CsrMatrixPMatVecHost() [1/4]

int pargemslr::CsrMatrixPMatVecHost ( const int *  ia,
const int *  ja,
const complexd aa,
int  nrows,
int  ncols,
char  trans,
const complexd alpha,
const complexd x,
const complexd beta,
complexd y 
)

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

Parameters
[in]iaThe column pointer vector.
[in]jaThe column vector.
[in]aaThe data pointer vector.
[in]nrowsThe number of rows in the vector.
[in]ncolssThe number of cols in the vector.
[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.

◆ CsrMatrixPMatVecHost() [2/4]

int pargemslr::CsrMatrixPMatVecHost ( const int *  ia,
const int *  ja,
const complexs aa,
int  nrows,
int  ncols,
char  trans,
const complexs alpha,
const complexs x,
const complexs beta,
complexs y 
)

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

Parameters
[in]iaThe column pointer vector.
[in]jaThe column vector.
[in]aaThe data pointer vector.
[in]nrowsThe number of rows in the vector.
[in]ncolssThe number of cols in the vector.
[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.

◆ CsrMatrixPMatVecHost() [3/4]

int pargemslr::CsrMatrixPMatVecHost ( const int *  ia,
const int *  ja,
const double *  aa,
int  nrows,
int  ncols,
char  trans,
const double &  alpha,
const double *  x,
const double &  beta,
double *  y 
)

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

Parameters
[in]iaThe column pointer vector.
[in]jaThe column vector.
[in]aaThe data pointer vector.
[in]nrowsThe number of rows in the vector.
[in]ncolssThe number of cols in the vector.
[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.

◆ CsrMatrixPMatVecHost() [4/4]

int pargemslr::CsrMatrixPMatVecHost ( const int *  ia,
const int *  ja,
const float *  aa,
int  nrows,
int  ncols,
char  trans,
const float &  alpha,
const float *  x,
const float &  beta,
float *  y 
)

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

Parameters
[in]iaThe column pointer vector.
[in]jaThe column vector.
[in]aaThe data pointer vector.
[in]nrowsThe number of rows in the vector.
[in]ncolssThe number of cols in the vector.
[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.

◆ CsrMatrixPMatVecHostTemplate()

template<typename T >
int pargemslr::CsrMatrixPMatVecHostTemplate ( const int *  ia,
const int *  ja,
const T *  aa,
int  nrows,
int  ncols,
char  trans,
const T &  alpha,
const T *  x,
const T &  beta,
T *  y 
)

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

Parameters
[in]iaThe column pointer vector.
[in]jaThe column vector.
[in]aaThe data pointer vector.
[in]nrowsThe number of rows in the vector.
[in]ncolssThe number of cols in the vector.
[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.

◆ CsrMatrixRcmBfsHost()

template<typename T >
int pargemslr::CsrMatrixRcmBfsHost ( CsrMatrixClass< T > &  G,
int  root,
vector_int marker,
std::vector< std::vector< int > > &  level 
)

This function builds the level structure starting from node root with BFS.

Parameters
[in]GThe target graph in CSR format.
[in]rootThe node number of the root node.
[in]markerHelper array of length equals to the size of G.
[out]levelThe level structure
Returns
Return error message.

◆ CsrMatrixRcmClearLevelHost()

int pargemslr::CsrMatrixRcmClearLevelHost ( std::vector< std::vector< int > > &  level)

This function frees the level structure.

Parameters
[in,out]levelThe level structure
Returns
Return error message.

◆ CsrMatrixRcmHost()

template<typename T >
int pargemslr::CsrMatrixRcmHost ( CsrMatrixClass< T > &  A,
vector_int perm 
)

This function computes the RCM ordering of A+trans(A).
Connected component will be ordered by the node number of the first node in the original A in ascending order. So, if you have a block diagonal matrix, using this RCM function won't influence the block diagonal structure.

Parameters
[in]AThe target CSR matrix.
[out]permThe RCM permutation.
Returns
Return error message.

◆ CsrMatrixRcmNumberingHost()

template<typename T >
int pargemslr::CsrMatrixRcmNumberingHost ( CsrMatrixClass< T > &  G,
int  root,
vector_int marker,
vector_int perm,
int &  current_num 
)

This function computes the RCM numbering for a connect component start from root.

Parameters
[in]GThe target graph in CSR format.
[in]rootThe node number of the root node.
[in]markerHelper array of length equals to the size of G. On return, value inside will be change so that we know this connect component is visited.
[in,out]permThe permutation of this connect component will be in here.
[in]current_numStart number of the current connect component.
Returns
Return error message.

◆ CsrMatrixRcmPerphnHost()

template<typename T >
int pargemslr::CsrMatrixRcmPerphnHost ( CsrMatrixClass< T > &  G,
int &  root,
vector_int marker 
)

This function find an end of the pseudo-peripheral starting from the root.

Parameters
[in]GThe target graph in CSR format.
[in]rootThe node number of the root node.
[in]markerHelper array of length equals to the size of G.
Returns
Return error message.
Note
Reference: Gibbs, Norman E., William G. Poole, Jr, and Paul K. Stockmeyer. "An algorithm for reducing the bandwidth and profile of a sparse matrix." SIAM Journal on Numerical Analysis 13.2 (1976): 236-250.

◆ CsrMatrixRcmReverseHost()

int pargemslr::CsrMatrixRcmReverseHost ( vector_int perm,
int  start,
int  end 
)

This function reverse an array between [start,end].

Parameters
[in,out]permThe permutation vector
[in]startThe start index.
[in]endThe end index(includd).
Returns
Return error message.

◆ CsrMatrixRcmRootHost()

template<typename T >
int pargemslr::CsrMatrixRcmRootHost ( CsrMatrixClass< T > &  G,
vector_int marker,
int &  root 
)

This function find the next root for the RCM ordering of graph G.
It serachs for the lowest degree node in the next connect component.

Parameters
[in]GThe target graph in CSR format.
[in]markerHelper array of length equals to the size of G.
[out]rootThe node number of the root node.
Returns
Return error message.

◆ CsrMatrixRcmSwapHost()

int pargemslr::CsrMatrixRcmSwapHost ( vector_int perm,
int  a,
int  b 
)

This function swap two element in the perm array.

Parameters
[in,out]permThe permutation vector
[in]aOne of the index.
[in]bThe other index.
Returns
Return error message.

◆ CsrMatrixSortRow()

template<typename T >
int pargemslr::CsrMatrixSortRow ( CsrMatrixClass< T > &  A)

Sort columns inside each row of a real csr matrix by ascending order.

Parameters
[in,out]AThe target matrix.
Returns
Return error message.

◆ CsrMatrixTransposeHost()

template<typename T >
int pargemslr::CsrMatrixTransposeHost ( CsrMatrixClass< T > &  A,
CsrMatrixClass< T > &  AT 
)

Compute the transpose of a csr/csc matrix.

Parameters
[in]AThe input matrix.
[out]ATThe output transpose matrix.
Returns
Return error message.

◆ CsrSubMatrixAmdHost()

template<typename T >
int pargemslr::CsrSubMatrixAmdHost ( CsrMatrixClass< T > &  A,
vector_int rowscols,
vector_int perm 
)

This function computes the AMD ordering of a submatrix A[rowscols,rowscols] + trans(A[rowscols,rowscols]).

Parameters
[in]AThe target CSR matrix.
[in]rowscolsThe target rows/cols.
[out]permThe rcm order.
Returns
Return error message.

◆ CsrSubMatrixNdHost()

template<typename T >
int pargemslr::CsrSubMatrixNdHost ( CsrMatrixClass< T > &  A,
vector_int rowscols,
vector_int perm 
)

This function computes the ND ordering of a submatrix A[rowscols,rowscols].

Parameters
[in]AThe target CSR matrix.
[in]rowscolsThe target rows/cols.
[out]permThe rcm order.
Returns
Return error message.

◆ CsrSubMatrixRcmHost()

template<typename T >
int pargemslr::CsrSubMatrixRcmHost ( CsrMatrixClass< T > &  A,
vector_int rowscols,
vector_int perm 
)

This function computes the RCM ordering of a submatrix A[rowscols,rowscols] + trans(A[rowscols,rowscols]).

Parameters
[in]AThe target CSR matrix.
[in]rowscolsThe target rows/cols.
[out]permThe rcm order.
Returns
Return error message.

◆ DenseMatrixComplexHessEigenDecompositionHost() [1/2]

int pargemslr::DenseMatrixComplexHessEigenDecompositionHost ( DenseMatrixClass< complexd > &  A,
DenseMatrixClass< complexd > &  Q,
vector_seq_complexd w 
)

Compute the Eigen decomposition of double complex dense hessenberg matrix A, AQ = QD. A is not modified.

Parameters
[in]AThe input matrix A.
[out]Qthe Q matrix AQ = QD.
[in,out]wArray of real part of eigenvalues.
Returns
Return error message.

◆ DenseMatrixComplexHessEigenDecompositionHost() [2/2]

int pargemslr::DenseMatrixComplexHessEigenDecompositionHost ( DenseMatrixClass< complexs > &  A,
DenseMatrixClass< complexs > &  Q,
vector_seq_complexs w 
)

Compute the Eigen decomposition of single complex dense hessenberg matrix A, AQ = QD. A is not modified.

Parameters
[in]AThe input matrix A.
[out]Qthe Q matrix AQ = QD.
[in,out]wArray of real part of eigenvalues.
Returns
Return error message.

◆ DenseMatrixComplexHessSchurDecompositionHost() [1/4]

int pargemslr::DenseMatrixComplexHessSchurDecompositionHost ( DenseMatrixClass< complexd > &  A,
DenseMatrixClass< complexd > &  Q,
vector_seq_complexd w 
)

Compute the Schur decomposition of a double complex dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[out]Qthe Q matrix A = QUQ^H.
[out]wArray of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixComplexHessSchurDecompositionHost() [2/4]

int pargemslr::DenseMatrixComplexHessSchurDecompositionHost ( DenseMatrixClass< complexd > &  A,
int  start,
int  end,
DenseMatrixClass< complexd > &  Q,
vector_seq_complexd w 
)

Compute the Schur decomposition of a double complex dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[in]start
[in]endOutside [start, end) the matrix is already upper triangular.
[out]Qthe Q matrix A = QUQ^H.
[out]wArray of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixComplexHessSchurDecompositionHost() [3/4]

int pargemslr::DenseMatrixComplexHessSchurDecompositionHost ( DenseMatrixClass< complexs > &  A,
DenseMatrixClass< complexs > &  Q,
vector_seq_complexs w 
)

Compute the Schur decomposition of a single complex dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[out]Qthe Q matrix A = QUQ^H.
[out]wArray of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixComplexHessSchurDecompositionHost() [4/4]

int pargemslr::DenseMatrixComplexHessSchurDecompositionHost ( DenseMatrixClass< complexs > &  A,
int  start,
int  end,
DenseMatrixClass< complexs > &  Q,
vector_seq_complexs w 
)

Compute the Schur decomposition of a single complex dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[in]start
[in]endOutside [start, end) the matrix is already upper triangular.
[out]Qthe Q matrix A = QUQ^H.
[out]wArray of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixComplexOrderSchur() [1/2]

int pargemslr::DenseMatrixComplexOrderSchur ( DenseMatrixClass< complexd > &  A,
DenseMatrixClass< complexd > &  Q,
vector_seq_complexd w,
vector_int select 
)

Reorder the double complex schur decomposition, puch selected eigenvalues to the leading part.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wArray of eigenvalues.
[in]selectEigenvalues corresponting to positive entries in select vector will be put in the leading part.
Returns
Return error message.

◆ DenseMatrixComplexOrderSchur() [2/2]

int pargemslr::DenseMatrixComplexOrderSchur ( DenseMatrixClass< complexs > &  A,
DenseMatrixClass< complexs > &  Q,
vector_seq_complexs w,
vector_int select 
)

Reorder the single complex schur decomposition, puch selected eigenvalues to the leading part. A 2x2 block on the diagonal of Q need to be in the same cluster.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wArray of eigenvalues.
[in]selectEigenvalues corresponting to positive entries in select vector will be put in the leading part.
Returns
Return error message.

◆ DenseMatrixComplexOrderSchurClusters() [1/2]

int pargemslr::DenseMatrixComplexOrderSchurClusters ( DenseMatrixClass< complexd > &  A,
DenseMatrixClass< complexd > &  Q,
vector_seq_complexd w,
vector_int clusters 
)

Reorder the double complex schur decomposition, puch selected eigenvalues to the leading part.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wArray of eigenvalues.
[in]clustersEigenvalues corresponding to positive entries in the clusters vector will be put in the leading part.
The order of those eigenvalues are based on the values in the cluster vector in descending order.

◆ DenseMatrixComplexOrderSchurClusters() [2/2]

int pargemslr::DenseMatrixComplexOrderSchurClusters ( DenseMatrixClass< complexs > &  A,
DenseMatrixClass< complexs > &  Q,
vector_seq_complexs w,
vector_int clusters 
)

Reorder the single complex schur decomposition, puch selected eigenvalues to the leading part.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wArray of eigenvalues.
[in]clustersEigenvalues corresponding to positive entries in the clusters vector will be put in the leading part.
The order of those eigenvalues are based on the values in the cluster vector in descending order.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [1/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< complexd > &  A,
DenseMatrixClass< complexd > &  Q 
)

Transform a dense double complex matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [2/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< complexd > &  A,
int  start,
int  end,
DenseMatrixClass< complexd > &  Q 
)

Transform a dense double complex matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[in]start
[in]endOutside [start, end) the matrix is already hessenberg.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [3/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< complexs > &  A,
DenseMatrixClass< complexs > &  Q 
)

Transform a dense single complex matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [4/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< complexs > &  A,
int  start,
int  end,
DenseMatrixClass< complexs > &  Q 
)

Transform a dense single complex matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[in]start
[in]endOutside [start, end) the matrix is already hessenberg.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [5/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< double > &  A,
DenseMatrixClass< double > &  Q 
)

Transform a dense double matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [6/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< double > &  A,
int  start,
int  end,
DenseMatrixClass< double > &  Q 
)

Transform a dense double matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[in]start
[in]endOutside [start, end) the matrix is already hessenberg.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [7/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< float > &  A,
DenseMatrixClass< float > &  Q 
)

Transform a dense float matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixHessDecompositionHost() [8/8]

int pargemslr::DenseMatrixHessDecompositionHost ( DenseMatrixClass< float > &  A,
int  start,
int  end,
DenseMatrixClass< float > &  Q 
)

Transform a dense float matrix A into hessenberg matrix Q^HAQ = Hess.

Parameters
[in,out]AThe input matrix A, on return, the hessenberg matrix.
[in]start
[in]endOutside [start, end) the matrix is already hessenberg.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixInvertHost() [1/4]

int pargemslr::DenseMatrixInvertHost ( DenseMatrixClass< complexd > &  A)

Compute the inverese of a general double complex matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixInvertHost() [2/4]

int pargemslr::DenseMatrixInvertHost ( DenseMatrixClass< complexs > &  A)

Compute the inverese of a general single complex matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixInvertHost() [3/4]

int pargemslr::DenseMatrixInvertHost ( DenseMatrixClass< double > &  A)

Compute the inverese of a general double matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixInvertHost() [4/4]

int pargemslr::DenseMatrixInvertHost ( DenseMatrixClass< float > &  A)

Compute the inverese of a general float matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixInvertUpperTriangularHost() [1/4]

int pargemslr::DenseMatrixInvertUpperTriangularHost ( DenseMatrixClass< complexd > &  A)

Compute the inverese of a upper triangular double complex matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixInvertUpperTriangularHost() [2/4]

int pargemslr::DenseMatrixInvertUpperTriangularHost ( DenseMatrixClass< complexs > &  A)

Compute the inverese of a upper triangular single complex matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixInvertUpperTriangularHost() [3/4]

int pargemslr::DenseMatrixInvertUpperTriangularHost ( DenseMatrixClass< double > &  A)

Compute the inverese of a upper triangular double matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixInvertUpperTriangularHost() [4/4]

int pargemslr::DenseMatrixInvertUpperTriangularHost ( DenseMatrixClass< float > &  A)

Compute the inverese of a upper triangular float matrix.

Parameters
[in,out]AThe dense matrix, on exit the inverse of it.
Returns
Return error message.

◆ DenseMatrixMatMat() [1/4]

int pargemslr::DenseMatrixMatMat ( const complexd alpha,
const DenseMatrixClass< complexd > &  A,
char  transa,
const DenseMatrixClass< complexd > &  B,
char  transb,
const complexd beta,
DenseMatrixClass< complexd > &  C 
)

Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C.
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]BThe second csr matrix.
[in]transbWhether we transpose B or not.
[in]betaThe beta value.
[out]CThe output C=A*B. Need not to allocate memory in advance.
Returns
Return error message.

◆ DenseMatrixMatMat() [2/4]

int pargemslr::DenseMatrixMatMat ( const complexs alpha,
const DenseMatrixClass< complexs > &  A,
char  transa,
const DenseMatrixClass< complexs > &  B,
char  transb,
const complexs beta,
DenseMatrixClass< complexs > &  C 
)

Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C.
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]BThe second csr matrix.
[in]transbWhether we transpose B or not.
[in]betaThe beta value.
[out]CThe output C=A*B. Need not to allocate memory in advance.
Returns
Return error message.

◆ DenseMatrixMatMat() [3/4]

int pargemslr::DenseMatrixMatMat ( const double &  alpha,
const DenseMatrixClass< double > &  A,
char  transa,
const DenseMatrixClass< double > &  B,
char  transb,
const double &  beta,
DenseMatrixClass< double > &  C 
)

Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C.
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]BThe second csr matrix.
[in]transbWhether we transpose B or not.
[in]betaThe beta value.
[out]CThe output C=A*B. Need not to allocate memory in advance.
Returns
Return error message.

◆ DenseMatrixMatMat() [4/4]

int pargemslr::DenseMatrixMatMat ( const float &  alpha,
const DenseMatrixClass< float > &  A,
char  transa,
const DenseMatrixClass< float > &  B,
char  transb,
const float &  beta,
DenseMatrixClass< float > &  C 
)

Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C.
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]BThe second csr matrix.
[in]transbWhether we transpose B or not.
[in]betaThe beta value.
[out]CThe output C=A*B. Need not to allocate memory in advance.
Returns
Return error message.

◆ DenseMatrixMatMatTemplate()

template<typename T >
int pargemslr::DenseMatrixMatMatTemplate ( const T &  alpha,
const DenseMatrixClass< T > &  A,
char  transa,
const DenseMatrixClass< T > &  B,
char  transb,
const T &  beta,
DenseMatrixClass< T > &  C 
)

Dense float Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C.
C is assumed to be pre-allocated.

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.
[out]CThe output C=A*B. Need not to allocate memory in advance.
Returns
Return error message.

◆ DenseMatrixMatVec() [1/4]

int pargemslr::DenseMatrixMatVec ( const DenseMatrixClass< complexd > &  A,
char  trans,
const complexd alpha,
const VectorClass< complexd > &  x,
const complexd beta,
VectorClass< complexd > &  y 
)

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

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

◆ DenseMatrixMatVec() [2/4]

int pargemslr::DenseMatrixMatVec ( const DenseMatrixClass< complexs > &  A,
char  trans,
const complexs alpha,
const VectorClass< complexs > &  x,
const complexs beta,
VectorClass< complexs > &  y 
)

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

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

◆ DenseMatrixMatVec() [3/4]

int pargemslr::DenseMatrixMatVec ( const DenseMatrixClass< double > &  A,
char  trans,
const double &  alpha,
const VectorClass< double > &  x,
const double &  beta,
VectorClass< double > &  y 
)

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

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

◆ DenseMatrixMatVec() [4/4]

int pargemslr::DenseMatrixMatVec ( const DenseMatrixClass< float > &  A,
char  trans,
const float &  alpha,
const VectorClass< float > &  x,
const float &  beta,
VectorClass< float > &  y 
)

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

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

◆ DenseMatrixPlotHost()

template<typename T >
int pargemslr::DenseMatrixPlotHost ( DenseMatrixClass< T > &  A,
int  conditiona,
int  conditionb,
int  width 
)

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

Parameters
[in]AThe target matrix.
[in]conditionaFirst condition.
[in]conditionbSecend condition, only spy when conditiona == conditionb.
[in]widthThe plot width.
Returns
Return error message.

◆ DenseMatrixPMatVecTemplate()

template<typename T >
int pargemslr::DenseMatrixPMatVecTemplate ( char  trans,
int  nrows,
int  ncols,
const T &  alpha,
const T *  aa,
int  ldim,
const T *  x,
const T &  beta,
T *  y 
)

In place dense complex 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]nrowsNumber of rows.
[in]ncolsNumber of columns.
[in]alphaThe alpha value.
[in]aaPointer to the data.
[in]ldimThe leading dimension of A.
[in]xThe left vector.
[in]betaThe beta value.
[in,out]yThe result vector.
Returns
Return error message.

◆ DenseMatrixQRDecompositionHost() [1/4]

int pargemslr::DenseMatrixQRDecompositionHost ( DenseMatrixClass< complexd > &  A,
DenseMatrixClass< complexd > &  Q 
)

Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn.

Parameters
[in,out]AThe input matrix A. On return the R matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixQRDecompositionHost() [2/4]

int pargemslr::DenseMatrixQRDecompositionHost ( DenseMatrixClass< complexs > &  A,
DenseMatrixClass< complexs > &  Q 
)

Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn.

Parameters
[in,out]AThe input matrix A. On return the R matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixQRDecompositionHost() [3/4]

int pargemslr::DenseMatrixQRDecompositionHost ( DenseMatrixClass< double > &  A,
DenseMatrixClass< double > &  Q 
)

Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn.

Parameters
[in,out]AThe input matrix A. On return the R matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixQRDecompositionHost() [4/4]

int pargemslr::DenseMatrixQRDecompositionHost ( DenseMatrixClass< float > &  A,
DenseMatrixClass< float > &  Q 
)

Compute the thin QR decomposition of mxn matrix A = QR. When m >= n, Q is mxn and R is nxn. When n > m, Q is mxm and R is mxn.

Parameters
[in,out]AThe input matrix A. On return the R matrix.
[out]Qthe Q matrix.
Returns
Return error message.

◆ DenseMatrixRealHessEigenDecompositionHost() [1/2]

int pargemslr::DenseMatrixRealHessEigenDecompositionHost ( DenseMatrixClass< double > &  A,
DenseMatrixClass< double > &  Q,
vector_seq_double wr,
vector_seq_double wi 
)

Compute the Eigen decomposition of double dense hessenberg matrix A, AQ = QD. A is not modified.

Parameters
[in]AThe input matrix A.
[out]Qthe Q matrix AQ = QD.
[in,out]wrArray of real part of eigenvalues.
[in]wiArray of imag part of eigenvalues.
Returns
Return error message.

◆ DenseMatrixRealHessEigenDecompositionHost() [2/2]

int pargemslr::DenseMatrixRealHessEigenDecompositionHost ( DenseMatrixClass< float > &  A,
DenseMatrixClass< float > &  Q,
vector_seq_float wr,
vector_seq_float wi 
)

Compute the Eigen decomposition of float dense hessenberg matrix A, AQ = QD. A is not modified.

Parameters
[in]AThe input matrix A.
[out]Qthe Q matrix AQ = QD.
[in,out]wrArray of real part of eigenvalues.
[in]wiArray of imag part of eigenvalues.
Returns
Return error message.

◆ DenseMatrixRealHessSchurDecompositionHost() [1/4]

int pargemslr::DenseMatrixRealHessSchurDecompositionHost ( DenseMatrixClass< double > &  A,
DenseMatrixClass< double > &  Q,
vector_seq_double wr,
vector_seq_double wi 
)

Compute the Schur decomposition of a double dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[out]Qthe Q matrix A = QUQ^H.
[out]wrArray of real part of the eigenvalues.
[out]wiArray of imag part of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixRealHessSchurDecompositionHost() [2/4]

int pargemslr::DenseMatrixRealHessSchurDecompositionHost ( DenseMatrixClass< double > &  A,
int  start,
int  end,
DenseMatrixClass< double > &  Q,
vector_seq_double wr,
vector_seq_double wi 
)

Compute the Schur decomposition of a double dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[in]start
[in]endOutside [start, end) the matrix is already upper triangular.
[out]Qthe Q matrix A = QUQ^H.
[out]wrArray of real part of the eigenvalues.
[out]wiArray of imag part of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixRealHessSchurDecompositionHost() [3/4]

int pargemslr::DenseMatrixRealHessSchurDecompositionHost ( DenseMatrixClass< float > &  A,
DenseMatrixClass< float > &  Q,
vector_seq_float wr,
vector_seq_float wi 
)

Compute the Schur decomposition of a float dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[out]Qthe Q matrix A = QUQ^H.
[out]wrArray of real part of the eigenvalues.
[out]wiArray of imag part of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixRealHessSchurDecompositionHost() [4/4]

int pargemslr::DenseMatrixRealHessSchurDecompositionHost ( DenseMatrixClass< float > &  A,
int  start,
int  end,
DenseMatrixClass< float > &  Q,
vector_seq_float wr,
vector_seq_float wi 
)

Compute the Schur decomposition of a float dense hessenberg matrix A = QUQ^H.

Parameters
[in,out]AThe input matrix A, on return, the U matrix.
[in]start
[in]endOutside [start, end) the matrix is already upper triangular.
[out]Qthe Q matrix A = QUQ^H.
[out]wrArray of real part of the eigenvalues.
[out]wiArray of imag part of the eigenvalues.
Returns
Return error message.

◆ DenseMatrixRealOrderSchur() [1/2]

int pargemslr::DenseMatrixRealOrderSchur ( DenseMatrixClass< double > &  A,
DenseMatrixClass< double > &  Q,
vector_seq_double wr,
vector_seq_double wi,
vector_int select 
)

Reorder the double schur decomposition, puch selected eigenvalues to the leading part.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wrArray of real part of eigenvalues.
[out]wiArray of imag part of eigenvalues.
[in]selectEigenvalues corresponting to positive entries in select vector will be put in the leading part.
Returns
Return error message.

◆ DenseMatrixRealOrderSchur() [2/2]

int pargemslr::DenseMatrixRealOrderSchur ( DenseMatrixClass< float > &  A,
DenseMatrixClass< float > &  Q,
vector_seq_float wr,
vector_seq_float wi,
vector_int select 
)

Reorder the float schur decomposition, puch selected eigenvalues to the leading part.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wrArray of real part of eigenvalues.
[out]wiArray of imag part of eigenvalues.
[in]selectEigenvalues corresponting to positive entries in select vector will be put in the leading part.
Returns
Return error message.

◆ DenseMatrixRealOrderSchurClusters() [1/2]

int pargemslr::DenseMatrixRealOrderSchurClusters ( DenseMatrixClass< double > &  A,
DenseMatrixClass< double > &  Q,
vector_seq_double wr,
vector_seq_double wi,
vector_int clusters 
)

Reorder the double schur decomposition, puch selected eigenvalues to the leading part. A 2x2 block on the diagonal of Q need to be in the same cluster.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wrArray of real part of eigenvalues.
[out]wiArray of imag part of eigenvalues.
[in]clustersEigenvalues corresponding to positive entries in the clusters vector will be put in the leading part.
The order of those eigenvalues are based on the values in the cluster vector in descending order.
Returns
Return error message.

◆ DenseMatrixRealOrderSchurClusters() [2/2]

int pargemslr::DenseMatrixRealOrderSchurClusters ( DenseMatrixClass< float > &  A,
DenseMatrixClass< float > &  Q,
vector_seq_float wr,
vector_seq_float wi,
vector_int clusters 
)

Reorder the float schur decomposition, puch selected eigenvalues to the leading part. A 2x2 block on the diagonal of Q need to be in the same cluster.

Parameters
[in,out]AA and Q are the original schur decomposition, on return, the new decomposition.
[in,out]QA and Q are the original schur decomposition, on return, the new decomposition.
[out]wrArray of real part of eigenvalues.
[out]wiArray of imag part of eigenvalues.
[in]clustersEigenvalues corresponding to positive entries in the clusters vector will be put in the leading part.
The order of those eigenvalues are based on the values in the cluster vector in descending order.
Returns
Return error message.

◆ GetMatrixPPrecision() [1/6]

PrecisionEnum pargemslr::GetMatrixPPrecision ( const MatrixClass< complexd > *  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPPrecision() [2/6]

PrecisionEnum pargemslr::GetMatrixPPrecision ( const MatrixClass< complexs > *  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPPrecision() [3/6]

PrecisionEnum pargemslr::GetMatrixPPrecision ( const MatrixClass< double > *  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPPrecision() [4/6]

PrecisionEnum pargemslr::GetMatrixPPrecision ( const MatrixClass< float > *  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPPrecision() [5/6]

PrecisionEnum pargemslr::GetMatrixPPrecision ( const MatrixClass< int > *  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPPrecision() [6/6]

PrecisionEnum pargemslr::GetMatrixPPrecision ( const MatrixClass< long int > *  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPrecision() [1/6]

PrecisionEnum pargemslr::GetMatrixPrecision ( const MatrixClass< complexd > &  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPrecision() [2/6]

PrecisionEnum pargemslr::GetMatrixPrecision ( const MatrixClass< complexs > &  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPrecision() [3/6]

PrecisionEnum pargemslr::GetMatrixPrecision ( const MatrixClass< double > &  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPrecision() [4/6]

PrecisionEnum pargemslr::GetMatrixPrecision ( const MatrixClass< float > &  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPrecision() [5/6]

PrecisionEnum pargemslr::GetMatrixPrecision ( const MatrixClass< int > &  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ GetMatrixPrecision() [6/6]

PrecisionEnum pargemslr::GetMatrixPrecision ( const MatrixClass< long int > &  mat)

Get the precision of a matrix.

Parameters
[in]matthe matrix.
Returns
Return the precision in PrecisionEnum form.

◆ ParallelCsrMatrixAddHost()

template<typename T >
int pargemslr::ParallelCsrMatrixAddHost ( ParallelCsrMatrixClass< T > &  A,
ParallelCsrMatrixClass< T > &  B,
ParallelCsrMatrixClass< T > &  C 
)

Add sparse matrix in Parallel CSR format, C= A + B.

Parameters
[in]AMatrix A.
[in]BMatrix B.
[out]CMatrix C = A + B.
Returns
Return error message.

◆ ParallelCsrMatrixSetupIOOrder()

template<typename T >
int pargemslr::ParallelCsrMatrixSetupIOOrder ( ParallelCsrMatrixClass< T > &  parcsr_in,
vector_int local_perm,
int &  nI,
CsrMatrixClass< T > &  B_mat,
CsrMatrixClass< T > &  E_mat,
CsrMatrixClass< T > &  F_mat,
ParallelCsrMatrixClass< T > &  C_mat,
int  perm_option,
bool  perm_c 
)

Setup the local order splits the interior and exterior nodes.

Parameters
[in]parcsr_inThe input matrix.
[out]local_permThe local permutation.
[out]nIThe number of local interior nodes.
[out]B_matThe local B matrix.
[out]E_matThe local E matrix.
[out]F_matThe local F matrix.
[out]C_matThe global coupling matrix.
[in]perm_optionThe permutation option. kIluReorderingNo, kIluReorderingRCM, kIluReorderingAMD.
[in]boolperm_c Permute B only or B and C?
Returns
Return error message.

◆ ParallelCsrMatrixSetupPermutationParallelND()

template<typename T >
int pargemslr::ParallelCsrMatrixSetupPermutationParallelND ( ParallelCsrMatrixClass< T > &  A,
bool  vertexsep,
int &  nlev,
long int  minsep,
vector_int map_v,
vector_int mapptr_v 
)

The main ND reordering function using ND in the ParMETIS. The number of level is fixed once the np is given.

Parameters
[in]AThe target parallel CSR matrix.
[in]vertexsepUse edge seperator or vertex seperator (currently not in use).
[in,out]nlevThe target number of levels (currently not in use).
[in]minsepMinimal size of the seperator (currently not in use).
Returns
map_v The map from node number to domain number.
mapptr_v The vector holds the start domian number on each level.

◆ ParallelCsrMatrixSetupPermutationParallelRKway()

template<typename T >
int pargemslr::ParallelCsrMatrixSetupPermutationParallelRKway ( ParallelCsrMatrixClass< T > &  A,
bool  vertexsep,
int &  nlev,
long int  ncomp,
long int  minsep,
long int  kmin,
long int  kfactor,
vector_int map_v,
vector_int mapptr_v,
bool  bj_last 
)

The first part of the function conatins the following steps:
1) On each level, partition the current adjacency graph with kway partitioning.
2) Find the adjacency graph of the interface matrix.
3) Repeat until we reach the last level or the interface matrix is too small.
After that, we re-assign all nodes to their new processor, and apply communication to transfer data.

Parameters
[in]AThe target parallel CSR matrix.
[in]vertexsetUse edge seperator or vertex seperator.
[in,out]nlevThe target number of levels.
[in]ncompNumber of components.
[in]minsepMinimal size of the seperator.
[in]kminMinimal size of the number of subdomains on a level.
[in]kfactorAt each level, ncomp is deviced by kfactor, for example, when ncomp = 16 and kfactor = 2, size would be 16, 8, 4, ....
[out]map_vThe map from node number to domain number.
[out]mapptr_vThe vector holds the start domian number on each level.
[in]bj_lastShould we treat last level with block Jacobi?
Note
We assume that A_mat has symmetric pattern.

◆ ParallelCsrMatrixTransposeHost()

template<typename T >
int pargemslr::ParallelCsrMatrixTransposeHost ( ParallelCsrMatrixClass< T > &  A,
ParallelCsrMatrixClass< T > &  AT 
)

Compute the transpose of a parallel csr matrix.

Parameters
[in]AThe input matrix.
[out]ATThe output transpose matrix.
Returns
Return error message.

◆ ParallelNDGetSeparator()

int pargemslr::ParallelNDGetSeparator ( vector_long vtxdist,
vector_long xadj,
vector_long adjncy,
vector_long map,
long int &  ndom1,
long int &  ndom2,
long int &  edge_cut,
parallel_log parlog 
)

This function finds the vertex separator based on the map information. A very simple appraoch, find the edge separator and only keep one side.

Parameters
[in]vtxdistThe distribution of rows for parMETIS.
[in]xadjThe column ptr for parMETIS.
[in]adjncyThe column numbers for parMETIS.
[in,out]mapThe subdomain number of each node, on return, edge separator would be marked as 2.
[out]ndom1Global number of interior nodes in dom1.
[out]ndom2Global number of interior nodes in dom2.
[out]edge_cutGlobal number of edge_cut.
[in]parallel_logParallel info struct.
Returns
return PARGEMSLR_SUCCESS or error information.
Note
We assume that A_mat has symmetric pattern.

◆ ParallelRKwayGetSeparator()

int pargemslr::ParallelRKwayGetSeparator ( vector_long vtxdist,
vector_long xadj,
vector_long adjncy,
vector_long vtxdist_s,
vector_long xadj_s,
vector_long adjncy_s,
vector_long map,
int  num_dom,
vector_int vtxsep,
parallel_log parlog 
)

This function finds the edge separator based on the map information.

Parameters
[in]vtxdistThe distribution of rows for parMETIS.
[in]xadjThe column ptr for parMETIS.
[in]adjncyThe column numbers for parMETIS.
[in]vtxdist_sThe distribution of rows for the edge separator.
[in]xadj_sThe column ptr for the edge separator.
[in]adjncy_sThe column numbers for the edge separator.
[in]mapThe subdomain number of each node.
[in]num_domThe number ofsubdomain number of each node.
[out]vtxsepvtxsep[i] > 0 if node i is inside the edge separator.
[in]parallel_logParallel info struct.
Returns
return PARGEMSLR_SUCCESS or error information.
Note
We assume that A_mat has symmetric pattern.

◆ PargemslrArnoldiNoRestart()

template<class VectorType , class MatrixType , typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiNoRestart ( MatrixType &  A,
int  mstart,
int  msteps,
DenseMatrixClass< DataType > &  V,
DenseMatrixClass< DataType > &  H,
RealDataType  tol_orth,
RealDataType  tol_reorth 
)

The standard Arnoldi procedure on a matrix.
Can start from a certan vector instead of the first vector in V.

Note
Template function, matrix requires SetupVectorPtrStr, GetNumRowsLocal, MatVec, and GetDataLocation functions.
Parameters
[in]AThe target matrix with Matvec and GetNumRows function.
[in]mstartthe index of the starting vector (start from V_mstart).
[in]msteps(Maximum) Number of Arnoldi steps is msteps - mstart.
[in,out]VMatrix holding the Krylov basis, space should be reserved, V(:,0:mstart) should be set already.
[in,out]HThe upper Hessenberg matrix generated by Arnoldi, space should be reserved.
[in]tol_orththe tolorance for the Arnoldi.
[in]tol_reorththe tolorance for the reorthgonation, set to -1.0 to turn off reorthgonation.
Returns
# of steps we actually do, might break before reaching msteps

◆ PargemslrArnoldiThickRestartBuildThickRestartNewVector()

template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiThickRestartBuildThickRestartNewVector ( DenseMatrixClass< DataType > &  V,
DenseMatrixClass< DataType > &  H,
int  m,
RealDataType  tol_orth,
RealDataType  tol_reorth,
VectorType &  v 
)

The current V is a invarient subspace, pick a new vector to restart.

Parameters
[in]VThe V at the current step.
[in]HThe H at the current step.
[in]mThe size of Hm.
[in]tol_orthThe threshold for orthgonal.
[in]tol_reorthThe threshold for reorthgonal.
[in]vWorking vector pointer, set by SetupVectorPtrStr().
Returns
Return error message.

◆ PargemslrArnoldiThickRestartChooseEigenValuesComplex()

template<typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiThickRestartChooseEigenValuesComplex ( DenseMatrixClass< DataType > &  H,
DenseMatrixClass< DataType > &  Q,
DataType  h_last,
RealDataType(*)(DataType)  weight,
RealDataType  truncate,
int &  ncov,
int &  nicov,
int &  nsatis,
RealDataType  tol_eig,
RealDataType  eig_target_mag,
RealDataType  eig_truncate,
bool &  cut,
SequentialVectorClass< DataType > &  w,
vector_int icov,
vector_int iicov,
vector_int isatis,
SequentialVectorClass< RealDataType > &  dcov,
SequentialVectorClass< RealDataType > &  dicov,
SequentialVectorClass< RealDataType > &  dsatis 
)

Pick convergenced eigenvalues and unconvergenced eigenvalues.

Parameters
[in]HThe H matrix after Schur decomposition.
[in]QThe Q matrix after Schur decomposition.
[in]h_lastThe H(m+1.m) value.
[in]weightThe weight function.
[out]ncovThe number of convergenced eigenvalues.
[out]nicovThe number of unconvergenced eigenvalues.
[in]wThe eigenvalues.
[out]icovThe index of convergenced eigenvalues.
[out]iicovThe index of unconvergenced eigenvalues.
[out]dcovThe distance of convergenced eigenvalues.
[out]dicovThe distance of unconvergenced eigenvalues.
Returns
Return error message.

◆ PargemslrArnoldiThickRestartChooseEigenValuesReal()

template<typename T >
int pargemslr::PargemslrArnoldiThickRestartChooseEigenValuesReal ( DenseMatrixClass< T > &  H,
DenseMatrixClass< T > &  Q,
h_last,
T(*)(ComplexValueClass< T >)  weight,
truncate,
int &  ncov,
int &  nicov,
int &  nsatis,
tol_eig,
eig_target_mag,
eig_truncate,
bool &  cut,
SequentialVectorClass< T > &  wr,
SequentialVectorClass< T > &  wi,
vector_int icov,
vector_int iicov,
vector_int isatis,
SequentialVectorClass< T > &  dcov,
SequentialVectorClass< T > &  dicov,
SequentialVectorClass< T > &  dsatis 
)

Pick convergenced eigenvalues and unconvergenced eigenvalues.

Parameters
[in]HThe H matrix after Schur decomposition.
[in]QThe Q matrix after Schur decomposition.
[in]h_lastThe H(m+1.m) value.
[in]weightThe weight function.
[out]ncovThe number of convergenced eigenvalues.
[out]nicovThe number of unconvergenced eigenvalues.
[in]wrThe real part of eigenvalues.
[in]wiThe imag part of eigenvalues.
[out]icovThe index of convergenced eigenvalues.
[out]iicovThe index of unconvergenced eigenvalues.
[out]dcovThe distance of convergenced eigenvalues.
[out]dicovThe distance of unconvergenced eigenvalues.
Returns
Return error message.

◆ PargemslrArnoldiThickRestartNoLock()

template<class VectorType , class MatrixType , typename DataType , typename RealDataType >
int pargemslr::PargemslrArnoldiThickRestartNoLock ( MatrixType &  A,
int  msteps,
int  maxits,
int  rank,
int  rank2,
RealDataType  truncate,
RealDataType  tr_fact,
RealDataType  tol_eig,
RealDataType  eig_target_mag,
RealDataType  eig_truncate,
RealDataType(*)(ComplexValueClass< RealDataType >)  weight,
DenseMatrixClass< DataType > &  V,
DenseMatrixClass< DataType > &  H,
RealDataType  tol_orth,
RealDataType  tol_reorth 
)

The thick-restart Arnoldi without locking. This version should be used when convergence tolerance is high.

Note
Template function, matrix requires SetupVectorPtrStr, GetNumRowsLocal, MatVec, and GetDataLocation functions.
Parameters
[in]AThe target matrix with Matvec and GetNumRows function.
[in]mstepsInitial number of Arnoldi steps.
[in]maxitsMaximum number of restarts.
[in]rankThe number of convergenced eigenvalues to compute.
[in]rank2The number of convergenced eigenvalues to keep.
[in]tr_factReal value between (0,1), thick-restart factor. Control the size of thick-restart.
Each restart would restart with number of convergenced plus number of unconvergenced eigenvalue * tr_fact.
[in]tol_eigThe tolorance of eigenvalue convergence.
[in]weightCompute eigenvalue has larger weight(val) first.
[in,out]VMatrix holding the Krylov basis, space should be reserved, V(:,0:mstart) should be set already.
[in,out]HThe upper Hessenberg matrix generated by Arnoldi, space should be reserved.
[in]tol_orththe tolorance for the Arnoldi.
[in]tol_reorththe tolorance for the reorthgonation, set to -1.0 to turn off reorthgonation.
Returns
# of convergenced eigenvalues computed.

◆ PargemslrCgs2()

template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrCgs2 ( VectorType &  w,
DenseMatrixClass< DataType > &  V,
DenseMatrixClass< DataType > &  H,
RealDataType &  t,
int  k,
RealDataType  tol_orth 
)

The classical Gram-Schmidt with re-orthgonal (always one step of re-orth).

Parameters
[in]wParallel vector to project to I-VV'.
[in,out]VMatrix holding the Krylov basis.
[in,out]HThe upper Hessenberg matrix generated by Arnoldi.
[out]tScalar to hold ||w||.
[in]kNumber of columns used in V, exclude w.
[in]tol_orthtol used to check break.
Returns
return GEMSLR_SUCESS or error information.

◆ PargemslrMgs()

template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrMgs ( VectorType &  w,
DenseMatrixClass< DataType > &  V,
DenseMatrixClass< DataType > &  H,
RealDataType &  t,
int  k,
RealDataType  tol_orth,
RealDataType  tol_reorth 
)

The modified Gram-Schmidt with re-orthgonal.

Parameters
[in]wParallel vector to project to I-VV'.
[in,out]VMatrix holding the Krylov basis.
[in,out]HThe upper Hessenberg matrix generated by Arnoldi.
[out]tScalar to hold ||w||.
[in]kNumber of columns used in V, exclude w.
[in]tol_orthtol used to check break.
[in]tol_reorthtol used in reorth for higher accuracy.
Returns
return GEMSLR_SUCESS or error information.

◆ PargemslrOrthogonal()

template<class VectorType , typename DataType , typename RealDataType >
int pargemslr::PargemslrOrthogonal ( VectorType &  w,
DenseMatrixClass< DataType > &  V,
RealDataType &  t,
int  k,
RealDataType  tol_orth,
RealDataType  tol_reorth 
)

The modified Gram-Schmidt with re-orthgonal without H matrix.

Parameters
[in]wParallel vector to project to I-VV'.
[in,out]VMatrix holding the Krylov basis.
[out]tScalar to hold ||w||.
[in]kNumber of columns used in V, exclude w.
[in]tol_orthtol used to check break.
[in]tol_reorthtol used in reorth for higher accuracy.
Returns
return GEMSLR_SUCESS or error information.

◆ PargemslrSubSpaceIteration()

template<class VectorType , class MatrixType , typename DataType , typename RealDataType >
int pargemslr::PargemslrSubSpaceIteration ( MatrixType &  A,
int  k,
int  its,
DenseMatrixClass< DataType > &  V,
DenseMatrixClass< DataType > &  H,
RealDataType  tol 
)

The subspace iteration procedure on a matrix. Note that we do not check for convergence. Apply a fixed number of iterations.

Note
Template function, matrix requires SetupVectorPtrStr, GetNumRowsLocal, MatVec, and GetDataLocation functions.
Parameters
[in]AThe target matrix with Matvec and GetNumRows function.
[in]kthe dimension of the subspace.
[in]its(Maximum) Number of iterations.
[in,out]VMatrix holding the final Q.
[in,out]HMatrix holding the final R.
[in]tolNot yet in used.
Returns
Return error message.

◆ ParmetisKwayHost()

int pargemslr::ParmetisKwayHost ( vector_long vtxdist,
vector_long xadj,
vector_long adjncy,
long int &  num_dom,
vector_long map,
parallel_log parlog 
)

Call ParMETIS to partition A based on the adjancency graph of A.

Note
See ParMETIS documentation for more details.
Parameters
[in]vtxdistVector of the number of the first vertex on each processor.
[in]xadjLocal I vector in the CSR format.
[in]adjncyLocal J vector in the CSR format.
[in,out]num_domNumber of requested partitions. On exit the number of partitions we get.
[out]mapThe domain number of each vertex.
[in]parallel_logParallel info struct.
Returns
return GEMSLR_SUCESS or error information.

◆ ParmetisNodeND()

int pargemslr::ParmetisNodeND ( vector_long vtxdist,
vector_long xadj,
vector_long adjncy,
long int &  num_dom,
vector_long map,
parallel_log parlog 
)

Partition with log(p) levels nd, p is the number of processors.

Note
See ParMETIS documentation for more details.
Parameters
[in]vtxdistVector of the number of the first vertex on each processor.
[in]xadjLocal I vector in the CSR format.
[in]adjncyLocal J vector in the CSR format.
[in,out]num_domNumber of requested partitions. On exit the number of partitions we get.
[out]mapThe domain number of each vertex.
[in]parallel_logParallel info struct.
Returns
return GEMSLR_SUCESS or error information.

◆ SetupPermutationNDCombineLevels()

int pargemslr::SetupPermutationNDCombineLevels ( std::vector< vector_int > &  level_stri,
int &  ndom 
)

Compress a certain level of level_str from SetupPermutationNDRecursive into a give number of domains.

Parameters
[in]level_striThe level of a level_str.
[in,out]ndomThe target number of domains, on return the ndom we get.
Returns
Return error message.

◆ SetupPermutationNDRecursive()

template<typename T >
int pargemslr::SetupPermutationNDRecursive ( CsrMatrixClass< T > &  A,
bool  vertexsep,
int  clvl,
int &  tlvl,
int  minsep,
std::vector< std::vector< vector_int > > &  level_str 
)

Recursive ND partition call.

Parameters
[in]AThe target matrix.
[in]vertexsetUse edge seperator or vertex seperator.
[in]clvlThe current level.
[in,out]tlvlThe total number of levels.
[in]num_domThe target number of domains on this level.
[in]minsepThe minimal size of the edge seperator.
[in]level_strThe level structure. level_str[i][j] is the nodes in the i-th component, of the j-th level in this component.
Returns
Return error message.

◆ SetupPermutationParallelKwayVertexSep()

template<typename T >
int pargemslr::SetupPermutationParallelKwayVertexSep ( ParallelCsrMatrixClass< T > &  A,
long int &  ncomp,
vector_int map_v,
vector_long perm_sep,
parallel_log parlog 
)

This function build the k-way partition where k is the power of 2, and provide an approximate vertex saperator.

Parameters
[in]AThe target matrix.
[in,out]ncompThe k for the kway partition, on return, the actual number of subdomains.
[out]map_vThe subdomain number of each node, starting from 0 to ncomp-1. if map_v[i] == ncomp this is the separator.
[out]perm_sepThe permutation, A(perm_sep, perm_sep) is the global coupling matrix.
[in]parallel_logParallel info struct.
Returns
return PARGEMSLR_SUCCESS or error information.
Note
We assume that A has symmetric pattern.

◆ SetupPermutationParallelKwayVertexSepRecursive()

template<typename T >
int pargemslr::SetupPermutationParallelKwayVertexSepRecursive ( ParallelCsrMatrixClass< T > &  A,
int  clvl,
int  tlvl,
bool &  succeed,
vector_int map_v,
parallel_log parlog 
)

This recursie call function for building the k-way partition where k is the power of 2, and provide an approximate vertex saperator.

Parameters
[in]AThe target matrix.
[in]clvlThe current level starting from 0.
[in]tlvlThe total number of levels, map from level to domain: 2->2, 3->4, 4->8 ....
[in,out]succeedTell if the partition get desired result.
[out]map_vThe subdomain number of each node, starting from 0 to ncomp-1. if map_v[i] == ncomp this is the separator.
[in]parallel_logParallel info struct.
Returns
return PARGEMSLR_SUCCESS or error information.
Note
We assume that A has symmetric pattern.

◆ SetupPermutationParallelRKwayRecursive()

int pargemslr::SetupPermutationParallelRKwayRecursive ( vector_long vtxdist,
vector_long xadj,
vector_long adjncy,
int  clvl,
int &  tlvl,
long int  ncomp,
long int  minsep,
long int  kmin,
long int  kfactor,
vector_int map_v,
vector_int mapptr_v,
bool  bj_last,
parallel_log parlog 
)

The function conatins the following steps:
1) On each level, partition the current adjacency graph with kway partitioning.
2) Find the adjacency graph of the interface matrix.
3) Repeat until we reach the last level or the interface matrix is too small.

Parameters
[in]vtxdistThe distribution of rows for parMETIS.
[in]xadjThe column ptr for parMETIS.
[in]adjncyThe column numbers for parMETIS.
[in]clvlCurrent level number.
[in,out]tlvlThe max number of levels. If the algorithm stops at a level lower than desired level, reset tlvl to current number of levels.
[in]ncompThe k for the kway partition on the current level.
[in]minsepThe minimal size of the edge separator. The algorithm stops if the size of the current interface matrix is too small.
[in]kminThe minimal value of k for the kway partition.
[in]kfactorOn the next level, ncomp = ncomp/kfactor. We might want more subdomains on higher levels.
[out]map_vThe subdomain number of each node.
[out]mapptr_vThe vector holds the start domian number on each level.
[in]bj_lastShould we treat last level with block Jacobi?
[in]parallel_logParallel info struct.
Returns
return PARGEMSLR_SUCCESS or error information.
Note
We assume that A_mat has symmetric pattern.

◆ SetupPermutationParallelRKwayRecursive2()

template<typename T >
int pargemslr::SetupPermutationParallelRKwayRecursive2 ( ParallelCsrMatrixClass< T > &  A,
int  clvl,
int &  tlvl,
long int  ncomp,
long int  minsep,
long int  kmin,
long int  kfactor,
vector_int map_v,
vector_int mapptr_v,
bool  bj_last,
parallel_log parlog 
)

This is the recursive function to build the mapping information, using an approximate vetrex separator (not the minimal).

Parameters
[in]AThe target matrix.
[in]clvlCurrent level number.
[in,out]tlvlThe max number of levels. If the algorithm stops at a level lower than desired level, reset tlvl to current number of levels.
[in]ncompThe k for the kway partition on the current level.
[in]minsepThe minimal size of the edge separator. The algorithm stops if the size of the current interface matrix is too small.
[in]kminThe minimal value of k for the kway partition.
[in]kfactorOn the next level, ncomp = ncomp/kfactor. We might want more subdomains on higher levels.
[out]map_vThe subdomain number of each node.
[out]mapptr_vThe vector holds the start domian number on each level.
[in]bj_lastShould we treat last level with block Jacobi?
[in]parallel_logParallel info struct.
Returns
return PARGEMSLR_SUCCESS or error information.
Note
We assume that A_mat has symmetric pattern.

◆ SetupPermutationRKwayRecursive()

template<typename T >
int pargemslr::SetupPermutationRKwayRecursive ( CsrMatrixClass< T > &  A,
bool  vertexsep,
int  clvl,
int &  tlvl,
int  num_dom,
int  minsep,
int  kmin,
int  kfactor,
vector_int map_v,
vector_int mapptr_v 
)

Recursive KWay partition call.

Parameters
[in]AThe target matrix.
[in]vertexsetUse edge seperator or vertex seperator.
[in]clvlThe current level.
[in,out]tlvlThe total number of levels.
[in]num_domThe target number of domains on this level, must be no less than 2.
[in]minsepThe minimal size of the edge seperator.
[in]kminThe minimal number of domains on this level.
[in]kfactorThe reduce factor of the number of domains.
[in]map_vThe map from node number to domain number.
[in]mapptr_vThe vector holds the start domian number on each level.
Returns
Return error message.