Warning: this documentation for the development version is under construction.

Functions for Matrices

In that page, we detail functions that are not related to Lapack

Transpose replaces a matrix by its transpose
TransposeConj replaces a matrix by its conjugate transpose
SetRow modifies a row of the matrix
GetRow extracts a row from the matrix
SetCol modifies a column of the matrix
GetCol extracts a column from the matrix
MaxAbs returns highest absolute value of A
Norm1 returns 1-norm of A
NormInf returns infinity-norm of A
ApplyInversePermutation permutes row and column numbers of a matrix
ScaleMatrix multiplies rows and columns of a matrix by coefficients
ScaleLeftMatrix multiplies rows of a matrix by coefficients
Copy copies a sparse matrix into another one (conversion of format if needed)
ConvertMatrix_to_Coordinates conversion of a sparse matrix into coordinates format
ConvertMatrix_from_Coordinates conversion of a matrix given as a triplet (i, j, val) to a sparse matrix
ConvertToCSC converts a sparse matrix to CSC (Compressed Sparse Column) format
ConvertToCSR converts a sparse matrix to CSR (Compressed Sparse Row) format
ConvertToSparse converts dense matrices to sparse matrices by specifying a threshold.
SOR applies successive over-relaxations to matrix
Cg, Gmres, BiCgSTAB, etc solves iteratively a linear system
IsComplexMatrix returns true if the matrix is complex
IsSparseMatrix returns true if the matrix is sparse
IsSymmetricMatrix returns true if the matrix is symmetric

Transpose, TransposeConj

Syntax :

  void Transpose(Matrix&);
  void TransposeConj(Matrix&);

Transpose overwrites a matrix by its transpose, while TransposeConj overwrites a matrix by its conjugate transpose. These functions are only available for dense matrices.

Example :

  
Matrix<double> A(5, 5);
A.Fill();

Transpose(A);

Matrix<complex<double> > B(5, 5);
// you fill B as you want, then overwrites it by conj(transpose(B))
TransposeConj(B);

Location :

Functions_Matrix.cxx

SetRow

Syntax :

  void SetRow(const Vector&, int, Matrix&);

This function modifies a row in the provided matrix. This function is available only for dense matrices.

Example :

  
Matrix<double> A(5, 5);
A.Fill();

// now you construct a row
Vector<double> row(5);
row.FillRand();

// and you put it in A
int irow = 1;
SetRow(row, irow, A);

Location :

Functions.cxx

GetRow

Syntax :

  void GetRow(const Matrix&, int, Vector&);

This function extracts a row from the provided matrix. This function is available only for dense matrices.

Example :

  
Matrix<double> A(5, 5);
A.Fill();

// now you extract the first row
Vector<double> row;
GetRow(A, 0, row);

Location :

Functions.cxx

SetCol

Syntax :

  void SetCol(const Vector&, int, Matrix&);

This function modifies a column in the provided matrix. This function is available only for dense matrices.

Example :

  
Matrix<double> A(5, 5);
A.Fill();

// now you construct a column
Vector<double> col(5);
col.FillRand();

// and you put it in A
int icol = 1;
SetCol(col, icol, A);

Location :

Functions.cxx

GetCol

Syntax :

  void GetCol(const Matrix&, int, Vector&);

This function extracts a column from the provided matrix. This function is available only for dense matrices.

Example :

  
Matrix<double> A(5, 5);
A.Fill();

// now you extract the first column
Vector<double> col;
GetCol(A, 0, col);

Location :

Functions.cxx

MaxAbs

Syntax :

  T MaxAbs(const Matrix<T>&);

This function returns the highest absolute value of a matrix.

Example :

  
Matrix<complex<double> > A(5, 5);
A.Fill();

double module_max = MaxAbs(A);

Location :

Functions_Matrix.cxx

Norm1

Syntax :

  T Norm1(const Matrix<T>&);

This function returns the 1-norm of a matrix.

Example :

  
Matrix<complex<double> > A(5, 5);
A.Fill();

double anorm_one = Norm1(A);

Location :

Functions_Matrix.cxx

NormInf

Syntax :

  T NormInf(const Matrix<T>&);

This function returns the infinite norm of a matrix.

Example :

  
Matrix<complex<double> > A(5, 5);
A.Fill();

double anorm_inf = NormInf(A);

Location :

Functions_Matrix.cxx

ApplyInversePermutation

Syntax :

  void ApplyInversePermutation(Matrix&, const Vector<int>&, const Vector<int>&);

This function permutes a given matrix with the provided new row numbers and column numbers. ApplyInversePermutation(A, row, col) does the same operation as A(row, col) = A in Matlab. This function is only available for RowMajor, ColMajor, ArrayRowSparse, ArrayRowSymSparse, ArrayRowComplexSparse and ArrayRowSymComplexSparse.

Example :

  
// you fill A as you wish
Matrix<double, Symmetric, ArrayRowSymSparse> A(5, 5);
// then new row and column numbers
IVect row(5);
// for symmetric matrix, row and column permutation must be the same
// if you want to keep symmetry
ApplyInversePermutation(A, row, row);

// for unsymmetric matrices, you can specify different permutations
IVect col(5);
Matrix<double, General, ArrayRowSparse> B(5, 5);
ApplyInversePermutation(B, row, col);

Location :

Permutation_ScalingMatrix.cxx

ScaleMatrix

Syntax :

  void ScaleMatrix(Matrix&, const Vector&, const Vector&);

This function multiplies each row and column with a coefficient, i.e. A is replaced by L*A*R where L and R are diagonal matrices and you provide the diagonal when you call ScaleMatrix. This function is only available for ArrayRowSparse, ArrayRowSymSparse, ArrayRowComplexSparse and ArrayRowSymComplexSparse.

Example :

  
// you fill A as you wish
Matrix<double, Symmetric, ArrayRowSymSparse> A(5, 5);
// then scaling vectors
Vector<double> scale(5);
// for symmetric matrix, row and column scaling must be the same
// if you want to keep symmetry
ScaleMatrix(A, scale, scale);

// for unsymmetric matrices, you can specify different scalings
IVect scale_right(5);
Matrix<double, General, ArrayRowSparse> B(5, 5);
ApplyInversePermutation(B, scale, scale_right);

Location :

Permutation_ScalingMatrix.cxx

ScaleLeftMatrix

Syntax :

  void ScaleLeftMatrix(Matrix&, const Vector&);

This function multiplies each row with a coefficient, i.e. A is replaced by L*A where L is diagonal and you provide the diagonal when you call ScaleLeftMatrix. This function is only available for ArrayRowSparse, and ArrayRowComplexSparse.

Example :

  
// you fill A as you wish
Matrix<double, General, ArrayRowSparse> A(5, 5);
// then scaling vector
Vector<double> scale(5);
ScaleLeftMatrix(A, scale);

Location :

Permutation_ScalingMatrix.cxx

Copy

Syntax :

  void Copy(const Matrix&, Matrix2&);

This function copies a sparse matrix into another one. If the types of the matrices differ, a conversion is performed. However, not all the conversion have been implemented, so you may have a compilation error when copying some matrices.

Example :

  
// you fill A as you wish
Matrix<double, General, ArrayRowSparse> A(5, 5);

// then you can copy it to another form
Matrix<double, General, RowSparse> B;

// B does not need to be allocated
Copy(A, B);

// For other types that don't compile
// you can use ConvertMatrix_to_Coordinates :
Matrix<double, Symmetric, ColSymSparse> C;
// conversion to triplet form (i, j, value)
IVect IndRow, IndCol;
Vector<double> Val;
ConvertMatrix_to_Coordinates(A, IndRow, IndCol, Val, 0, true);
// then C is filled
ConvertMatrix_from_Coordinates(IndRow, IndCol, Val, C, 0);

Location :

Matrix_Conversions.cxx

ConvertMatrix_to_Coordinates

Syntax :

  void ConvertMatrix_to_Coordinates(const Matrix& A, Vector& IndRow, Vector& IndCol, Vector& Val, int index, bool sym);

This function converts a sparse matrix to the triplet form (i, j, val) (coordinate format). The row and column numbers will start at index, therefore you can switch between 1-based indices and 0-based indices.

Example :

  
// you fill A as you wish
Matrix<double, General, ArrayRowSparse> A(5, 5);

// conversion to triplet form (i, j, value)
IVect IndRow, IndCol;
Vector<double> Val;
ConvertMatrix_to_Coordinates(A, IndRow, IndCol, Val, 0, true);
// number of non-zero entries :
int nnz = IndRow.GetM();

for (int i = 0; i < nnz; i++)
  {
    cout << "Row index : " << IndRow(i) << endl;
    cout << "Column index : " << IndCol(i) << endl;
    cout << "value : " << Val(i) << endl;
  }

Location :

Matrix_Conversions.cxx

ConvertMatrix_from_Coordinates

Syntax :

  void ConvertMatrix_from_Coordinates(Vector& IndRow, Vector& IndCol, Vector& Val, Matrix& A, int index);

This function converts a triplet form (i, j, val) (coordinate format) to a sparse matrix. The row and column numbers are assumed to start at index, therefore you can switch between 1-based indices and 0-based indices.

Example :

  
// creating a sparse matrix
// A = | 1    0   0   2|
//     | -1   2   3   0|
//     |  0  1.2  0  2.5|
int nnz = 7;
IVect IndRow(nnz), IndCol(nnz);
Vector<double> Val(nnz);
IndRow(0) = 0; IndCol(0) = 0; Val(0) = 1.0;
IndRow(1) = 0; IndCol(1) = 3; Val(1) = 2.0;
IndRow(2) = 1; IndCol(2) = 0; Val(2) = -1.0;
IndRow(3) = 1; IndCol(3) = 1; Val(3) = 2.0;
IndRow(4) = 1; IndCol(4) = 2; Val(4) = 3.0;
IndRow(5) = 2; IndCol(5) = 1; Val(5) = 1.2;
IndRow(6) = 2; IndCol(6) = 3; Val(6) = 2.5;

// conversion to a Seldon structure
Matrix<double, General, RowSparse>
ConvertMatrix_from_Coordinates(IndRow, IndCol, Val, A, 0);

Location :

Matrix_Conversions.cxx

ConvertToCSC

Syntax :

  void ConvertToCSC(const Matrix& A, Property& sym, Vector& Ptr, Vector& Ind, Vector& Val, bool sym_pat);

This function converts a matrix to Compressed Sparse Column (CSC) format. Val stores the values of non-zero entries, Ind the row indexes of the non-zero entries, and Ptr the locations in Val of non-zero entries starting a column. This is the storage represented by ColSparse in Seldon. If Property is Symmetric, only upper part of the matrix will be converted (ColSymSparse storage). If sym_pat is true, the sparsity pattern is symmetrized, that is to say that if a(i, j) exists, then a(j, i) also exists. Default value of sym_pat is false. This feature is used for exemple in the interface of Pastix solver, since this solver requires that the sparsity pattern is symmetric (values may be not symmetric).

Example :

  

// you fill A as you wish
Matrix<double, General, ArrayRowSparse> A(5, 5);

// then you can retrieve Ptr, Ind, Val arrays of CSC format
General prop;
ConvertToCSC(A, prop, Ptr, Ind, Val);

ConvertToCSR

Syntax :

  void ConvertToCSR(const Matrix& A, Property& sym, Vector& Ptr, Vector& Ind, Vector& Val);

This function converts a matrix to Compressed Sparse Row (CSR) format. Val stores the values of non-zero entries, Ind the column indexes of the non-zero entries, and Ptr the locations in Val of non-zero entries starting a row. This is the storage represented by RowSparse in Seldon. If Property is Symmetric, only upper part of the matrix will be converted (RowSymSparse storage).

Example :

  

// you fill A as you wish
Matrix<double, General, ArrayColSparse> A(5, 5);

// then you can retrieve Ptr, Ind, Val arrays of CSR format
General prop;
ConvertToCSR(A, prop, Ptr, Ind, Val);

Location :

Matrix_Conversions.cxx

ConvertToSparse

Syntax :

  void ConvertToSparse(const Matrix& A, Matrix& B, const T& threshold);

This function converts a dense matrix to a sparse matrix. All values whose modulus is below or equal to threshold are skipped.

Example :

  

// you fill A as you wish
Matrix<double, General, RowMajor> A(5, 5);

// then you can convert it to a sparse matrix
Matrix<double, General, RowSparse> B;
ConvertToSparse(A, B, 1e-12);

// and retrieve the number of non-zero entries
int nnz = B.GetDataSize();

Location :

Matrix_Conversions.cxx

IsComplexMatrix

Syntax :

  void IsComplexMatrix(const Matrix&)

This function returns true if the matrix is complex. It does not check if all the values are real, but merely returns true if the value type is complex (e.g. T = complex<double>).

Example :

  

// complex matrix
Matrix<complex<double>, General, RowMajor> A;

// IsComplexMatrix should return true
if (IsComplexMatrix(A))
  { 
    cout << "A is complex" << endl;
  }

Location :

Functions_Matrix.cxx

IsSymmetricMatrix

Syntax :

  void IsSymmetricMatrix(const Matrix&)

This function returns true if the matrix is symmetric. It does not check if a(i,j) = a(j,i) for all i and j, but merely returns true if the property of the matrix is set to symmetric.

Example :

  

// complex matrix
Matrix<complex<double>, Symmetric, RowSymPacked> A;

// IsSymmetricMatrix should return true
if (IsSymmetricMatrix(A))
  { 
    cout << "A is symmetric" << endl;
  }

Location :

Functions_Matrix.cxx

IsSparseMatrix

Syntax :

  void IsSparseMatrix(const Matrix&)

This function returns true if the matrix is sparse. It does not check if all the values are different from 0, but merely returns true if the storage associated with the matrix is a sparse storage.

Example :

  

// complex matrix
Matrix<complex<double>, General, ArrayRowSparse> A;

// IsSparseMatrix should return true
if (IsSparseMatrix(A))
  { 
    cout << "A is sparse" << endl;
  }

Location :

Functions_Matrix.cxx