Dense Matrices
    

In that page, methods and functions related to dense matrices are detailed.

Basic declaration :

Classes :

Matrix<T,Prop,RowMajor>
Matrix<T,Prop,RowSymPacked>
Matrix<T,Prop,RowHermPacked>
Matrix<T,Prop,RowUpTriangPacked>
Matrix<T,Prop,RowLoTriangPacked>

Example :

// dense matrix of doubles
Matrix<double> A;

// dense symmetric matrix
Matrix<float, Symmetric, RowSymPacked> B;

// dense hermitian matrix
Matrix<complex<double>, Hermitian, RowHermPacked> C;

// dense lower triangular matrix
Matrix<double, General, RowLoTriangPacked> D;

// dense upper triangular matrix
Matrix<double, General, RowUpTriangPacked> E;

Below we list the methods specific to dense matrices. All the matrices are inherited from the class VirtualMatrix

.

Methods :

Matrix constructors
Matrix operators
GetDataSize returns the number of elements effectively stored
GetData returns a pointer to the array containing the values
GetDataConst returns a pointer to the array containing the values
GetDataVoid returns a pointer to the array containing the values
GetDataConstVoid returns a pointer to the array containing the values
Resize changes the size of matrix (keeps previous elements)
SetData sets the pointer to the array containing the values
Nullify clears the matrix without releasing memory
Val access to a matrix element
Get access to a matrix element
Set modifies an entry of the matrix
Copy copies a matrix
Zero sets all elements to zero
SetIdentity sets matrix to identity matrix
Fill sets all elements to a given value
FillRand fills randomly the matrix
Print displays the matrix
Write writes the matrix in binary format
Read reads the matrix in binary format
WriteText writes the matrix in text format
ReadText reads the matrix in text format


VirtualMatrix

The virtual methods declared in the class VirtualMatrix are common to all matrices, and present if and only if SELDON_WITH_VIRTUAL has been defined before including Seldon. The methods Clear, Reallocate and Zero are becoming virtual with this flag. This flag SELDON_WITH_VIRTUAL is set by default in the file SeldonFlag.hxx. If you are including SeldonLib.hxx, you do not need to set this flag since it should be done in SeldonFlag.hxx. Below an example of use of polymorphism for matrices :

#include "SeldonLib.hxx"

using namespace Seldon;

int main()
{
  // A can be a sparse or dense matrix
  VirtualMatrix<double>* A;

  bool dense = false;
  
  // creation of the object depending on boolean dense
  if (dense)
    A = new Matrix<double>;
  else 
    A = new Matrix<double, General, ArrayRowSparse>;

  // allocation of the matrix
  int n = 100;
  A->Reallocate(n, n);
  A->Zero();
  
  // adding values to A
  A->AddInteraction(1, 3, 4.61);
  A->AddInteraction(0, 2, -0.72);

  // performing a matrix vector product
  Vector<double> x(n), y(n); x.FillRand();
  A->MltVector(x, y);

  return 0;
}

Below we list the methods of the class VirtualMatrix

GetM returns the number of rows in the matrix
GetN returns the number of columns in the matrix
GetSize returns the number of elements in the matrix
Clear removes all elements of the matrix
Reallocate changes the size of matrix (does not keep previous elements)
SetEntry changes an entry of the matrix
AddInteraction adds a given value to an entry of the matrix
AddInteraction adds a given value to an entry of the matrix
AddInteractionRow adds several values to entries in a given row of the matrix
AddInteractionColumn adds several values to entries in a given column of the matrix
AddDistantInteraction adds a given value to an entry of the matrix (the column belongs to another processor)
AddRowDistantInteraction adds a given value to an entry of the matrix (the row belongs to another processor)
GetMemorySize returns the memory used to store the matrix in bytes
ClearRow clears a row of the matrix
ApplySor performs a over-relaxation step
MltAddVector performs a matrix-vector product
MltVector performs a matrix-vector product
IsComplex returns true if the matrix is complex
IsSymmetric returns true if the matrix is symmetric


Functions :

Mlt multiplication by a scalar or matrix-vector product
MltAdd performs a matrix-vector or matrix-matrix product
Add adds two matrices
Copy copies a matrix into another one
Rank1Update Adds a contribution X.Y' to a matrix
Rank2Update Adds a contribution X.Y' + Y.X' to a symmetric matrix
Solve solves a triangular system
Transpose replaces a matrix by its transpose
TransposeConj replaces a matrix by its conjugate transpose
MaxAbs returns highest absolute value of A
Norm1 returns 1-norm of A
NormInf returns infinity-norm of A
GetRow returns a matrix row
SetRow changes a matrix row
GetCol returns a matrix column
SetCol changes a matrix column
ApplyPermutation permutes rows and columns of a matrix
ApplyInversePermutation permutes rows and columns of a matrix
ScaleMatrix multiplies rows and columns by coefficients
ScaleLeftMatrix multiplies rows by coefficients
ScaleRightMatrix multiplies columns by coefficients
SOR performs successive over-relaxation algorithm
GetLU performs a LU (or LDL^t) factorization
SolveLU solve linear system by using LU factorization
RefineSolutionLU improves solution computed by SolveLU
ReciprocalConditionNumber computes the inverse of matrix condition number
GetScalingFactors computes row and column scalings to equilibrate a matrix
GetInverse computes the matrix inverse
GetQR QR factorization of matrix
GetLQ LQ factorization of matrix
GetQ_FromQR Forms explicitly Q from QR factorization
MltQ_FromLQ multiplies vector by Q
GetQ_FromLQ Forms explicitly Q from LQ factorization
MltQ_FromQR multiplies vector by Q
SolveQR solves least-square problems by using QR factorization
SolveLQ solves least-square problems by using LQ factorization
GetEigenvalues computes eigenvalues
GetEigenvaluesEigenvectors computes eigenvalues and eigenvectors
GetSVD performs singular value decomposition (SVD)



Matrix constructors

Syntax :

  Matrix();
  Matrix(int, int );

Example :

// default constructor -> empty matrix
Matrix<int> V;
cout << "Number of elements "<< V.GetSize() << endl; // should return 0 
// then you can use Reallocate to set the number of rows and columns
V.Reallocate(3, 2);
V.Fill();

// we construct matrix with 4 rows and 3 columns
Matrix<double> W(4, 3);
// W is not initialized, you have to fill it
W.Fill(1.0);

Related topics :

Reallocate
Fill

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Matrix operators

Syntax :

  const T& operator (int i, int j) const;
  T& operator (int i, int j);
  T& operator [int i];
  const T& operator [int i] const;
  Matrix& operator =(const Matrix& )
  Vector& operator =(const T0& alpha)
  Vector& operator *=(const T0& alpha)

The operator [] gives a direct access to the 1-D array storing values of the matrix, its use should be considered carefully. The operator () can't be used to modify values of matrix for specific storages, e.g. RowSym, RowHerm, RowUpTriang, use Get or Val instead.

Example :

Matrix<double> V(3, 3);
// use of operator () to modify matrix
V(0, 0) = 2.0;
V(1, 0) = V(0, 0) + 1.0;

// operator [] should be used with caution
V[3] = V[0] + 1.4;

Matrix<double> W;
// use of operator = to copy contents of vector V
W = V;

// set all elements to a given value
W = 1;

// multiplication by a scalar
Matrix<double> A(3, 2);
A.Fill();
A *= 1.5;

Related topics :

Copy
Val

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

GetM

Syntax :

  int GetM() const;

This method returns the number of rows.

Example :

Matrix<float> V(3, 2);
// V.GetM() should return 3 
cout << "Number of rows of V " << V.GetM() << endl;

Location :

Class Matrix_Base
Matrix_Base.hxx
Matrix_Base.cxx

GetN

Syntax :

  int GetN() const;

This method returns the number of columns.

Example :

Matrix<float> V(3, 2);
// V.GetN() should return 2 
cout << "Number of columns of V " << V.GetN() << endl;

Location :

Class Matrix_Base
Matrix_Base.hxx
Matrix_Base.cxx

GetSize

Syntax :

  int GetSize() const;

This method returns the number of elements in the matrix.

Example :

Matrix<float, Symmetric, RowSymPacked> V(3, 3);
// V.GetSize() should return 9
cout << "Number of elements of V " << V.GetSize() << endl;

Location :

Class Matrix_Base
Matrix_Base.hxx
Matrix_Base.cxx

GetDataSize

Syntax :

  int GetDataSize() const;

This method returns the number of elements effectively stored in the matrix. This is different from GetSize for some storages, e.g. RowSymPacked, RowHermPacked, RowUpTriangPacked.

Example :

Matrix<float, Symmetric, RowSymPacked> V(3, 3);
// V.GetDataSize() should return 6
cout << "Number of elements of V " << V.GetDataSize() << endl;

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

GetData, GetDataConst, GetDataVoid, GetDataConstVoid

Syntax :

  T* GetData() const;
  const T* GetDataConst() const;
  void* GetDataVoid() const;
  const void* GetDataConstVoid() const;

These methods are useful to retrieve the pointer to the values. In practice, you can use those methods in order to interface with C/fortran subroutines or to perform some low level operations. But in this last case, you have to be careful, because debugging operations will be more tedious.

Example :

Matrix<double> V(3, 4); V.Fill();
double* data = V.GetData();
// you can use data as a normal C array
// here the sum of elements is computed
double sum = 0;
for (int i = 0; i < V.GetDataSize(); i++)
  sum += data[i];

// if you want to call a fortran subroutine daxpy
Matrix<double> X(3, 3); 
double coef = 2.0;
int m = X.GetM();
int n = X.GetN();
daxpy_(&coef, &m, &n, X.GetData(),);

// for complex numbers, conversion to void* is needed :
Matrix<complex<double> > Xc(4, 4);
complex<double> beta(1,1);
zaxpy(reinterpret_cast<const void*>(&beta), Xc.GetDataVoid());

Related topics :

SetData
Nullify

Location :

Class Matrix_Base
Matrix_Base.hxx
Matrix_Base.cxx

Clear

Syntax :

  void Clear();

This method removes all the elements of the matrix.

Example :

Matrix<double> A(3, 2);
A.Fill();
// clears matrix A
A.Clear();

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Reallocate

Syntax :

  void Reallocate(int, int);

This method changes the size of the matrix, but removes previous elements if NewAlloc is used.

Example :

Matrix<long int> A(5, 4);
V.Fill();
// resizes matrix A
A.Reallocate(4, 3);
// you need to initialize all elements of A
A.Zero();

Related topics :

Resize

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Resize

Syntax :

  void Resize(int, int);

This method changes the size of the matrix, and keeps previous elements.

Example :

Matrix<long double> A(3,3);
A.Fill();
// resizes matrix A
A.Resize(4,4);
// you need to initialize new elements if there are new
for (int i = 0; i < 4; i++)
  A(4,i) = A(i,4) = 0;

Related topics :

Reallocate

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

SetData

Syntax :

  void SetData(int, int, T*);

This method sets the pointer to the array containing elements. This method should be used carefully, and generally in conjunction with method Nullify.

Example :

// for example, you can define a function with a pointer as argument
void f(int m, int n, double* data)
{
  // and sets this array into a Matrix instance
  Matrix<double> A;
  // m : number of rows, n : number of columns
  A.SetData(m, n, data);
  // then you use a C++ method
  double coef = Norm1(A);
  // you don't release memory, because data is used after the function
  A.Nullify();
}

Related topics :

GetData
Nullify

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Nullify

Syntax :

  void Nullify();

This method clears the matrix without releasing memory. This method should be used carefully, and generally in conjunction with method Nullify. You can look at the example shown in the explanation of method SetData.

Related topics :

SetData
GetData

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Val

Syntax :

  T& Val(int, int);
  const T& Val(int, int) const;

This method is similar to operator (), except that it can be always used to modify the matrix, especially for storages like RowSym, RowHerm, RowUpTriang. It gives a direct access to the values of the matrix, so that it is up to the user to take care that only upper part of symmetric matrices is accessed. In the case of sparse matrices, this method modifies only entries that are belonging to the sparsity pattern. An exception is raised, if you try to modify an entry which does not exist. If you prefer a more flexible access, use Get instead.

Example :

Matrix<double, General, RowUpTriang> A(3,3);
// operator () does not work to change the value
// A(0, 0) = 1;  => Error during compilation
A.Val(0, 0) = 2.0;  // Okay it works

Related topics :

Operators

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Get

Syntax :

  T& Get(int, int);
  const T& Get(int, int) const;

This method is similar to operator (), except that it can be used to modify the matrix, especially for storages like RowSym, RowHerm, RowUpTriang. It gives access to A(i, j). The only restriction to that method occurs for hermitian matrix. In that case, only upper part can be modified, if you wish to modify lower part, the method Set should be considered.

Example :

Matrix<double, General, RowUpTriang> A(3,3);
// operator () does not work to change the value
// A(0,0) = 1;  => Error during compilation
A.Get(0, 0) = 2.0;  // Okay it works

Related topics :

Operators

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Set, SetEntry

Syntax :

  void Set(int, int, T);
  void SetEntry(int, int, T);

These methods allow the modification of A(i, j), there is no restriction on i or j. The method SetEntry is virtual (and belongs to the class VirtualMatrix).

Example :

Matrix<double, General, RowUpTriang> A(3,3);
// operator () does not work to change the value
// A(0,0) = 1;  => Error during compilation
A.Set(0, 0, 2.0);  // Okay it works

// if you want to use the polymorphism, you can use SetEntry
// it can be slower (virtual method)
VirtualMatrix<double>* B;
B = new Matrix<double, Symmetric, RowSymPacked>(3, 3);
B->SetEntry(0, 1, -2.5); // B(0, 1) = -2.5

Related topics :

Operators

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Copy

Syntax :

  void Copy(const Matrix&);

This method copies a matrix into the current matrix.

Example :

// copy of a matrix M
Matrix<double> M(3, 3), A;
M.FillRand();
A.Copy(M);
// this is equivalent to use operator =
A = M;

Related topics :

Matrix operators

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Zero

Syntax :

  void Zero();

This method fills memory of 0, is convenient for matrices made of doubles, integers, floats, but not for more complicated types. In that case, it is better to use the method Fill.

Example :

Matrix<double> V(5, 3);
// initialization
V.Fill();

Matrix<IVect> W(10, 2);
// W.Zero() is incorrect and would generate an error at the execution
// a good initialization is to use Fill
IVect zero(5); zero.Zero();
W.Fill(zero);

Related topics :

Fill

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

SetIdentity

Syntax :

  void SetIdentity();

This method sets all elements to 0, except on the diagonal set to 1. This forms the so-called identity matrix.

Example :

Matrix<double> V(5, 5);
// initialization
V.SetIdentity();

Related topics :

Fill

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Fill

Syntax :

  void Fill();
  template<class T0>
  void Fill(const T0& );

This method fills matrix with 0, 1, 2, etc or with a given value.

Example :

Matrix<int> A(2,2);
A.Fill();
// A should contain [0 1; 2 3]

A.Fill(2);
// A should contain [2 2; 2 2]

Related topics :

Zero

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

FillRand

Syntax :

  void FillRand();

This method fills the matrix with random values.

Example :

Matrix<double> A(5, 3);
A.FillRand();
// A should contain 15 random values

Related topics :

Fill

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Print

Syntax :

  void Print() const;
  void Print(int) const;
  void Print(int, int, int, int) const;

This method displays the matrix.

Example :

Matrix<string> A(2, 2);
A(0, 0) = string("hello");
A(0, 1) = string("world");
A(1, 0) = string("you");
A(1, 1) = string("welcome");
A.Print(); 
// should display :
// hello world
// you welcome

// you can also display a sub-matrix
A.Print(0, 0, 0, 1); 
// should display "hello world"

// A.Print(2); is equivalent to A.Print(0, 0, 2, 2);

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Write

Syntax :

  void Write(string) const;
  void Write(ofstream&) const;
  void Write(ofstream&, bool with_size) const;

This method writes the matrix on a file/stream in binary format. The file will contain the number of rows, columns, then the list of elements. You can also require that the number of rows and columns is not written by setting with_size to false.

Example :

Matrix<double> A(2, 2), B(2, 2); 
// you can write directly in a file
A.Fill();
A.Write("matrix.dat");

// or open a stream with other datas
ofstream file_out("matrix.dat");
int my_info = 3;
file_out.write(reinterpret_cast<char*>(>my_info), sizeof(int));
A.Write(file_out);
// you can write directly values of B without the number of rows and columns 
B.Write(file_out, false);
file_out.close();

Related topics :

Read
WriteText
ReadText

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

Read

Syntax :

  void Read(string);
  void Read(ifstream&);
  void Read(ifstream&, bool with_size);

This method sets the matrix from a file/stream in binary format. The file contains the number of rows, columns, then the list of elements. If with_size is false, the number of rows and columns won't be read.

Example :

Matrix<double> V, W(3, 4); 
// you can read directly on a file
V.Read("matrix.dat");

// or read from a stream
ifstream file_in("matrix.dat");
int my_info;
file_in.read(reinterpret_cast<char*<(>my_info), sizeof(int));
// V may be resized depending on the matrix stored on the file
V.Read(file_in);
// if the size is not written, W is considered to be allocated
// with the correct size (here 3x4)
W.Read(file_in, false);
file_in.close();

Related topics :

Write
WriteText
ReadText

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

WriteText

Syntax :

  void WriteText(string) const;
  void WriteText(ofstream&) const;

This method writes the matrix on a file/stream in text format. The file will contain the list of elements. The matrix is then directly readable in Matlab (load) or in Python (loadtxt).

Example :

Matrix<double> V(2); 
// you can write directly in a file
V.Fill();
// for more digits, cout.precision should be called
cout.precision(15);
V.WriteText("matrix.dat");

// or open a stream with other datas
ofstream file_out("matrix.dat");
int my_info = 3;
file_out << my_info << '\n';
V.WriteText(file_out);
file_out.close();

Related topics :

Write
Read
ReadText

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

ReadText

Syntax :

  void ReadText(string);
  void ReadText(ifstream&);

This method sets the matrix from a file/stream in text format. The file contains the list of elements.

Example :

Matrix<double> V; 
// you can read directly on a file
V.ReadText("matrix.dat");

// or read from a stream
ifstream file_in("matrix.dat");
int my_info;
file_in >> my_info;
V.ReadText(file_in);
file_in.close();

Related topics :

Write
Read
WriteText

Location :

Class Matrix_Pointers
Class Matrix_Symmetric
Class Matrix_SymPacked
Class Matrix_Hermitian
Class Matrix_HermPacked
Class Matrix_Triangular
Class Matrix_TriangPacked
Matrix_Pointers.hxx Matrix_Pointers.cxx
Matrix_Symmetric.hxx Matrix_Symmetric.cxx
Matrix_SymPacked.hxx Matrix_SymPacked.cxx
Matrix_Hermitian.hxx Matrix_Hermitian.cxx
Matrix_HermPacked.hxx Matrix_HermPacked.cxx
Matrix_Triangular.hxx Matrix_Triangular.cxx
Matrix_TriangPacked.hxx Matrix_TriangPacked.cxx

AddInteraction

Syntax :

  void AddInteraction(int, int, T);

This methods adds/inserts a value in the matrix. This is "almost" equivalent to use the method Get() in conjunction with operator +=. Indeed, this is truly equivalent for unsymmetric matrices. However, it is different for symmetric matrices, since in that case AddInteraction ignores values given on the lower part of the matrix, whereas Get() will modify the corresponding non-zero entry located in the upper part. The advantage here is that if you fill a matrix with AddInteraction/AddInteractionRow, you don't have to change your function if the matrix is symmetric or not. If you are using the methods Get or Set, you will have to be careful to not add twice the same value because A.Get(j, i) and A.Get(i, j) will point to the same address for symmetric matrices.

Example :

Matrix<double, General, ArrayRowSparse> A(5, 5);
// insertion
A.AddInteraction(1, 3, 2.6);
// we add
A.AddInteraction(1, 3, -1.0);
// now A should be equal to [1 3 1.6]

// for symmetric matrices, only upper part is considered
Matrix<double, General, ArrayRowSymSparse> B(10, 10);
B.AddInteraction(4, 7, -0.8);
B.AddInteraction(7, 4, -0.8); 
// B(4, 7) and B(7, 4) will be equal to -0.8
B.AddInteraction(3, 5, 2.2);
B.AddInteraction(5, 3, 0.9);
// B(3, 5) and B(5, 3) will be equal to 2.2

Related topics :

AddInteractionRow
AddInteractionColumn

Location :

Class Matrix_ArraySparse
Matrix_ArraySparse.hxx Matrix_ArraySparse.cxx

AddInteractionRow

Syntax :

  void AddInteractionRow(int, int, Vector<int>, Vector<T>);

This methods adds/inserts several values in a row of the matrix. This is more efficient to use that method rather than calling AddInteraction several times. For symmetric matrices, the interactions belonging to the lower part of the matrix will be ignored. The advantage of this behaviour is that you can use the same code for filling a symmetric and an unsymmetric matrix.

Example :

Matrix<double, General, ArrayRowSparse> A(5, 5);
// insertion, you don't need to sort column numbers
int irow = 2;
int nb_values = 3;
IVect col(nb_values);
Vector<double> values(nb_values);
col(0) = 0; values(0) = 0.1;
col(1) = 3; values(1) = 0.6;
col(2) = 2; values(2) = -1.4;

A.AddInteractionRow(irow, nb_values, col, values);

Related topics :

AddInteraction
AddInteractionColumn

Location :

Class Matrix_ArraySparse
Matrix_ArraySparse.hxx Matrix_ArraySparse.cxx

AddInteractionColumn

Syntax :

  void AddInteractionColumns(int, int, Vector<int>, Vector<T>);

This methods adds/inserts several values in a column of the matrix. For symmetric matrices, the interactions belonging to the lower part of the matrix will be ignored. The advantage of this behaviour is that you can use the same code for filling a symmetric and an unsymmetric matrix.

Example :

Matrix<double, General, ArrayRowSparse> A(5, 5);
// insertion, you don't need to sort row numbers
int icol = 2;
int nb_values = 3;
IVect row(nb_values);
Vector<double> values(nb_values);
row(0) = 0; values(0) = 0.1;
row(1) = 3; values(1) = 0.6;
row(2) = 2; values(2) = -1.4;

A.AddInteractionColumn(icol, nb_values, row, values);

Related topics :

AddInteraction
AddInteractionRow

Location :

Class Matrix_ArraySparse
Matrix_ArraySparse.hxx Matrix_ArraySparse.cxx

AddDistantInteraction

Syntax

 void AddDistantInteraction(int i, int jglob, int proc, const T& val);

This member function adds val for the local row i, and the global row jglob, proc being the processor that treats the global row jglob.

Example :

// default constructor
DistributedMatrix<double, General, ArrayRowSparse> A;

// then Reallocate is called with the local number of rows
A.Reallocate(n, n);

// after filling nglob, global_row, overlap_row and overlap_proc
// calling Init in order to provide these datas
A.Init(nglob, &global_row, &overlap_row, &overlap_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI::COMM_WORLD);

// for local interaction, you have to use AddInteraction
A.AddInteraction(i, j, val);

// when the column is distant (ie located on another processor), you have to use AddDistantInteraction
// jglob : global column number
// proc : distant processor
// val : value to add
A.AddDistantInteraction(i, jglob, proc, val);

Location :

matrix_sparse/DistributedMatrix.cxx

AddRowDistantInteraction

Syntax

 void AddRowDistantInteraction(int iglob, int j, int proc, const T& val);

This member function adds val for the global row iglob, and the local column j, proc being the processor that treats the global row iglob.

Example :

// default constructor
DistributedMatrix<double, General, ArrayRowSparse> A;

// then Reallocate is called with the local number of rows
A.Reallocate(n, n);

// after filling nglob, global_row, overlap_row and overlap_proc
// calling Init in order to provide these datas
A.Init(nglob, &global_row, &overlap_row, &overlap_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI::COMM_WORLD);

// for a local interaction, you have to use AddInteraction
A.AddInteraction(i, j, val);

// when the column is distant (ie located on another processor), you have to use AddDistantInteraction
// i : local row number
// jglob : global column number
// proc : distant processor
// val : value to add
// for the global matrix, it means that A_{global_row(i), jglob} is incremented with val
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you have to use AddRowDistantInteraction
// iglob : global row number
// j : local column number
// proc : distant processor
// val : value to add
// for the global matrix, it means that A_{iglob, global_row(j)} is incremented with val
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

matrix_sparse/DistributedMatrix.cxx

GetMemorySize

Syntax

 size_t GetMemorySize()

This member function returns the size used to store the matrix in bytes.

Example :

int n = 100;
Matrix<double> A(n, n);

// memory used to store A in bytes
size_t taille = A.GetMemorySize();

Location :

matrix/Matrix_Base.cxx

ClearRow

Syntax

 void ClearRow(int i)

This member function clears a row of the matrix. For dense matrices, the row will be filled of zeros. For sparse matrices, the row is cleared (deallocated).

Example :

int n = 100;
Matrix<double> A(n, n);
A.FillRand();

// if you want to clear the row 4
A.ClearRow(4);

Location :

matrix/Matrix_Base.cxx

ApplySor

Syntax

 void ApplySor(SeldonTrans, Vector& x, const Vector& b, omega, int nb_iter, int stage);

This member function applies nb_iter iterations of Successive Over Relaxation method to solve A x = b by starting with the vector given in x. x is overwritten with the updated solution. omega is the relaxation parameter. stage is the type of sweep (2 for a forward sweep, 3 for a backward sweep and 0 for a forward sweep followed by a backward sweep). This function should be implemented for most types of matrices (dense or sparse).

Example :

int n = 100;
Matrix<double> A(n, n);

Vector<double> x(n), b(n);
b.FillRand();
x.Zero();

int p = 10, stage = 0;
// p steps of S.S.O.R method to find an approximation solution of A x = b
A.ApplySor(SeldonNoTrans, x, b, 1.6, p, stage);

// if you want to solve A^T x = b :
x.Zero();
int p = 10, stage = 0;
// p steps of S.S.O.R method to find an approximation solution of A x = b
A.ApplySor(SeldonTrans, x, b, 1.6, p, stage);

Location :

matrix/Matrix_Base.cxx

MltAddVector

Syntax

 void MltAddVector(alpha, const Vector& x, beta, Vector& y);
 void MltAddVector(alpha, SeldonTranspose, const Vector& x, beta, Vector& y);

This member function adds to y the matrix-vector product A x or its transpose.

Example :

int n = 100;
Matrix<double> A(n, n);

Vector<double> x(n), y(n);
x.FillRand();
y.FillRand();

double alpha = 1.5, beta = 0.5;

// y = beta y + alpha A x
A.MltAddVector(alpha, x, beta, y);

// y = beta y + alpha A^T x
A.MltAddVector(alpha, SeldonTrans, x, beta, y);

Location :

matrix/Matrix_Base.cxx

MltVector

Syntax

 void MltVector(const Vector& x, Vector& y);
 void MltVector(SeldonTranspose, const Vector& x, Vector& y);

This member function puts in y the matrix-vector product A x or its transpose.

Example :

int n = 100;
Matrix<double> A(n, n);

Vector<double> x(n), y(n);
x.FillRand();

// y = A x
A.MltVector(x, y);

// y = A^T x
A.MltVector(SeldonTrans, x, y);

Location :

matrix/Matrix_Base.cxx

IsSymmetric

Syntax

bool IsSymmetric()

This member function returns true if the matrix is symmetric (by construction)..

Example :

int n = 100;
Matrix<double> A(n, n);

// sym should be false here
bool sym = A.IsSymmetric();

Matrix<double, Symmetric, RowSymPacked> B;

// sym should be true here
sym = B.IsSymmetric();

Location :

matrix/Matrix_Base.cxx

IsComplex

Syntax

bool IsComplex()

This member function returns true if the matrix is complex (by construction)..

Example :

int n = 100;
Matrix<double> A(n, n);

// cplx should be false here
bool cplx = A.IsComplex();

Matrix<complex<double> > B;

// cplx should be true here
sym = B.IsComplex();

Location :

matrix/Matrix_Base.cxx