Distributed Matrices
    

An implementation of distributed over several processors is proposed. It can be used to solve large linear systems with direct solvers (Mumps, SuperLU and Pastix are interfaced), or iterative solvers (all iteratives solvers provided in Seldon are compliant). Only square matrices (same number of rows and columns) are currently supported. Unitary tests for this class are present in the file test/unit/distributed_matrix_test.cc. This file must be compiled with the option USE_MPI := YES (in the Makefile), and executed with at least two processors (mpirun -np 2 ./test.x).

Basic use

// on each processor, you specify rows that are already treated by another processor
// for example if the processor 0 handles rows [0, 3, 5, 6] and
// processor 1 the rows [1 2, 4, 5, 7], the row 5 is already treated by
// processor 0. Of course, if each row is treated by a processor and only one
// this array should be left empty
IVect OverlapRow;
int n = 4;
if (MPI::COMM_WORLD.Get_rank() == 1)
  {
    n = 5;
    OverlapRow.Reallocate(1);
    // be careful because OverlapRow stores local numbers
    // here the global row 5 has a local number equal to 3 on processor 1
    OverlapRow(0) = 3;
  }

// in the constructor, you need to provide this array and also the communicator 
// (you can share a vector between the processors you want by
// constructing the appropriate communicator)
DistributedVector<double> U(OverlapRow, MPI_COMM_WORLD);

// another solution is to use the copy constructor
DistributedVector<double> V(U);

// you can use all the methods of a classical vector
U.Reallocate(n);
U.Fill(1.3);
V.FillRand();

// functions DotProd, DotProdConj and Norm2 have been overloaded
// so that iterative algorithms are working with distributed vectors
// other functions have not been overloaded
double scal = DotProd(U, V);

// for distributed matrices
// we also need to know the global row numbers
// and for each overlapped row, the original processor
IVect glob_number, original_proc;

// and also rows that are shared with other processors
// with the same numbering between SharedRows(i) of processor j
// and SharedRows(j) of processor i (to ensure an exact correspondence)
IVect ProcShared;
Vector<IVect> SharedRows;
if (MPI::COMM_WORLD.Get_rank() == 0)
  {
    glob_number.Reallocate(4);
    glob_number(0) = 0;
    glob_number(1) = 3;
    glob_number(2) = 5;
    glob_number(3) = 6;

    // global row 5 is shared with processor 1
    ProcShared.Reallocate(1);
    ProcShared(0) = 1;

    // local numbers are specified here (local row is 2, global row is 5)
    SharedRows.Reallocate(1);
    SharedRows(0).PushBack(2);
  }
else
  {
    glob_number.Reallocate(5);
    glob_number(0) = 1;
    glob_number(1) = 2;
    glob_number(2) = 4;
    glob_number(3) = 5;
    glob_number(4) = 7;

    // row 5 is already handled by processor 0
    original_proc.Reallocate(1);
    original_proc(0) = 0;

    // global row 5 is shared with processor 0
    ProcShared.Reallocate(1);
    ProcShared(0) = 0;

    // local numbers are specified here (local row is 3, global row is 5)
    SharedRows.Reallocate(1);
    SharedRows(0).PushBack(3);

  }

// initialisation of the matrix
// initial Seldon storages are used
// here n is the number of local rows (in our example 4 for proc 0, 5 for proc 1)
DistributedMatrix<double, General, ArrayRowSparse> A(n, n);

// you need to give also the number of global rows
// overlapped rows and processors, shared rows and MPI communicator
// shared rows can be given only for the first unknown, then Nvol is the number
// of rows for each unknown and nb_u the number of unknowns
int nglob = 8;
int Nvol = n, nb_u = 1;
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

// other syntaxes simpler for Init exist
// you can for example only provide global row numbers (see documentation of Init)

// you can use the methods of the base class Matrix
// for example if you use AddInteraction, it will concern local numbers
// this interaction will concern only the current processor
double val = 1.34;
// val is added to A(5, 6) (local number of global row 5 is 2, and 3 for global row 6)
if (MPI::COMM_WORLD.Get_rank() == 0)
  A.AddInteraction(2, 3, val);

// if you wish to add an interaction concerning several processors,
// you can use AddDistantInteraction if the row is local (and column number is global)
// or AddRowDistantInteraction when the column is local
// for example A(3, 7) = A(3, 7) + val gives 
// the third argument is the rank of the distant processor
if (MPI::COMM_WORLD.Get_rank() == 0)
  A.AddDistantInteraction(1, 7, 1, val);

// and A(7, 3) = A(7, 3) + val
if (MPI::COMM_WORLD.Get_rank() == 0)
  A.AddRowDistantInteraction(7, 1, 1, val);

// once the matrix is constructed, you can factorize it with Mumps/Pastix/SuperLU
SparseDistributedSolver<double> mat_lu;
mat_lu.Factorize(A);

// and solve the linear system
Vector<double> x(n);
// on input, x contains the right hand side, on output the solution
mat_lu.Solve(x);


Methods for distributed matrices :

Distributed matrix constructors
GetCommunicator returns the MPI communicator associated with the matrix
GetLocalM returns the local number of rows
GetLocalN returns the local number of columns
GetGlobalM returns the global number of rows
GetNodlScalar returns the number of rows for a scalar unknown
GetNbScalarUnknowns returns the number of scalar unknowns
AddDistantInteraction adds a value to A(i, j) where i is local and j global
AddRowDistantInteraction adds a value to A(i, j) where i is global and j local
GetMaxDataSizeDistantCol returns the number of values in distant columns stored in all the processors
GetMaxDataSizeDistantRow returns the number of values in distant rows stored in all the processors
IsReadyForMltAdd returns true if the structure is ready to perform a matrix vector without preparation
GetDistantColSize returns the number of distant non-zero entries for row i
IndexGlobalCol returns the global column number of distant non-zero entry j of local row i
ProcessorDistantCol returns the processor number of distant non-zero entry j of local row i
ValueDistantCol returns the value associated with distant non-zero entry j of local row i
GetDistantRowSize returns the number of distant non-zero entries for column i
IndexGlobalRow returns the global row number of distant non-zero entry j of local column i
ProcessorDistantRow returns the processor number of distant non-zero entry j of local column i
ValueDistantRow returns the value associated with distant non-zero entry j of local column i
Init sets pointers to the arrays containing global row numbers and overlapped rows
Reallocate changes the size of the matrix
Resize changes the size of the matrix and keeps previous entries
Clear erases the matrix
GetGlobalRowNumber returns local to global numbering
GetOverlapRowNumber returns the array containing the numbers of rows already handled by an another processor
GetOverlapProcNumber returns the array containing the processor numbers of rows already handled by an another processor
GetProcessorSharingRows returns the list of processors that share rows with the current processor
GetSharingRowNumbers returns the list of rows shared with each processor
GetMemorySize returns the size used by this object in bytes
GetNonZeros returns the number of non-zero entries stored in the matrix
GetDataSize returns the number of non-zero entries stored in the matrix
RemoveSmallEntry drops small non-zero entries present in the matrix
SetIdentity sets the matrix to the identity matrix
Fill fills non-zero entries with 0, 1, 2, etc
FillRand fills non-zero entries with random values
Write writes the distributed matrix on several files (one per processor) in binary format
WriteText writes the distributed matrix on several files (one per processor) in ascii format
Read reads the distributed matrix on several files (one per processor) in binary format
ReadText reads the distributed matrix on several files (one per processor) in ascii format


Functions for distributed matrices :

Mlt performs a matrix-vector product
MltAdd performs a matrix-vector product
Add adds a distributed matrix to another one
MaxAbs returns the maximum absolute value of entries in the matrix
GetRowSum computes the sum of absolute values for each row
GetColSum computes the sum of absolute values for each column
GetRowColSum computes the sum of absolute value for each column and row
Norm1 returns the 1-norm of the matrix
NormInf returns the infinite norm of the matrix
Transpose computes the transpose of a distributed matrix
Conjugate computes the conjugate of a distributed matrix
TransposeConj computes the transpose conjugate of a distributed matrix
GetRow extracts a row of a distributed matrix
GetCol extracts a column of a distributed matrix
SetRow sets a row of a distributed matrix
SetCol sets a column of a distributed matrix
ApplyPermutation applies a permutation (rows and columns) to a distributed matrix
ApplyInversePermutation applies the inverse of a permutation to a distributed matrix
SOR Applies Successive Over Relaxation (S.O.R) iterations to solve a linear system
Copy Copies/Converts a distributed matrix into another one
GetSubMatrix extracts a sub-matrix from a given distributed matrix
NormFro returns the Froebenius norm of the matrix
ScaleLeftMatrix scales rows of the matrix
ScaleRightMatrix scales columns of the matrix
ScaleMatrix scales rows and columns of the matrix
AssembleDistributed assembles a distributed matrix
EraseRow erases rows in the distributed matrix
EraseCol erases columns in the distributed matrix

DistributedMatrixIntegerArray

The class DistributedMatrixIntegerArray stores all the integers arrays needed for distributed matrices. You can use this class to encapsulate these arrays. Below, we list the global variables and methods of this class

OverlapRowNumbers row numbers already treated by another processor
OverlapProcNumbers number of the original processor for shared rows
GlobalRowNumbers global row numbers
ProcSharingRows list of processor numbers that share rows with the current processor
SharingRowNumbers for each processor, list of shared rows
nodl_scalar number of rows for a scalar unknown
nb_unknowns number of scalar unknowns
nglob global number of rows
nloc local number of rows
comm MPI communicator
Clear clears stored arrays
SetData copies pointers to fill arrays
Nullify nullifies stored arrays
ConstructArrays constructs stored arrays from global row numbers



GetCommunicator

Syntax

  MPI_Comm& GetCommunicator();

This method returns the MPI communicator used to distribute the vector or the matrix.

Example :

// for example if the processor 0 handles rows [0, 3, 5, 6] and
// processor 1 the rows [1 2, 4, 5, 7], the row 5 is already treated by
// processor 0 :
IVect OverlapRow;
int n = 4;
if (MPI::COMM_WORLD.Get_rank() == 1)
  {
    n = 5;
    OverlapRow.Reallocate(1);
    // be careful because OverlapRow stores local numbers
    // here the global row 5 has a local number equal to 3 on processor 1
    OverlapRow(0) = 3;
  }

DistributedVector<double> U(OverlapRow, MPI_COMM_WORLD);

// returns the communicator (here it should return MPI_COMM_WORLD)
MPI_Comm& comm = U.GetCommunicator();

Location :

DistributedMatrix.hxx
DistributedMatrixInline.cxx

Constructors of DistributedMatrix

Syntax

  DistributedMatrix();
  DistributedMatrix(int m, int n);
  DistributedMatrix(const DistributedMatrix&);

For a distributed matrix, the same constructors as for sequential matrices are available. The informations about communication with other processors are provided with the method Init.

Example :

// for the constructor, we specify the local number of rows and columns
DistributedMatrix<double, General, ArrayRowSparse> A(n, n);

// then you need to initialise the arrays used for parallel functions
// (e.g. matrix-vector product, factorization)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

// then the matrix can be modified with AddInteraction, AddDistantInteration and AddRowDistantInteraction
A.AddInteraction(0, 3, 2.5);

// if the copy constructor is used, there is no need to call Init since it will copy the pointers contained in A
DistributedMatrix<double, General, ArrayRowSparse> B(A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetLocalM, GetLocalN

Syntax

   int GetLocalM() const
   int GetLocalN() const

This method returns the local number of rows/columns of the distributed matrix. Only square matrices are supported, the matrix has the same number of rows and columns.

Example :

// for the constructor, you use the local number of rows
DistributedMatrix<double, General, ArrayRowSparse> A(n, n);

// then you need to initialise the arrays used for parallel functions
// (e.g. matrix-vector product)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

// GetLocalM, GetLocalN returns the n provided in the constructor or reallocate
cout << "Number of rows in this processor = " << A.GetLocalM()  << endl;

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetGlobalM

Syntax

 int GetGlobalM() const;

This function returns the global number of rows.

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, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

// you can access to nglob by calling GetGlobalM
nglob = A.GetGlobalM();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetNodlScalar

Syntax

 int GetNodlScalar() const;

This method returns the local number of rows for each scalar component. Several components can be specified when Init is called.

Example :

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

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

// total number of rows (with all components)
int nglob = 18;

// global numbers (for all components)
IVect glob_number;

// overlapped rows (with original processor)
IVect OverlapRow, original_proc;

// local number of rows (with all components)
int n;

// local number of rows for only one component
int Nvol;

// number of components
int nb_u = 2;

// corresponding rows that are shared between processors
// for only one component
IVect ProcShared;
Vector<IVect> SharedRows;

if (MPI::COMM_WORLD.Get_rank() == 0)
  {
    // for example, processor 0 owns rows 1, 3, 4, 6, 7, 8 and 10, 12, 13, 15, 16, 17 
    n = 12;
    Nvol = 6;

    glob_number.Reallocate(n);
    glob_number(0) = 1; glob_number(1) = 3; glob_number(2) = 4;
    glob_number(3) = 6; glob_number(4) = 7; glob_number(5) = 8;
    glob_number(6) = 10; glob_number(7) = 12; glob_number(8) = 13;
    glob_number(9) = 15; glob_number(10) = 16; glob_number(11) = 17;

    // rows 3, 6 and 12, 15 are shared with processor 1
    // local numbers are 1, 3 and 7, 9
    // for ProcShared/SharedRows, only one component is needed
    ProcShared.Reallocate(1); SharedRows.Reallocate(1);
    ProcShared(0) = 1;
    SharedRows(0).PushBack(1); SharedRows(0).PushBack(3);

  }
else
  {
    // for example, processor 1 owns rows 0, 2, 3, 5, 6 and 9, 11, 12, 14, 15
    n = 10;
    Nvol = 5;

    glob_number.Reallocate(n);
    glob_number(0) = 0; glob_number(1) = 2; glob_number(2) = 3;
    glob_number(3) = 5; glob_number(4) = 6;
    glob_number(5) = 9; glob_number(6) = 11; glob_number(7) = 12;
    glob_number(8) = 14; glob_number(9) = 15;

    // rows 3, 6 and 12, 15 are originally in processor 0
    OverlapRow.Reallocate(4);
    // local numbers are 2, 4 and 7, 9
    OverlapRow(0) = 2; OverlapRow(1) = 4;
    OverlapRow(2) = 7; OverlapRow(3) = 9;
    original_proc.Reallocate(4); original_proc.Fill(0);

    // for ProcShared/SharedRows, only one component is needed
    ProcShared.Reallocate(1); SharedRows.Reallocate(1);
    ProcShared(0) = 0;
    SharedRows(0).PushBack(2); SharedRows(0).PushBack(4);
  }

// after filling nglob, global_row, overlap_row and overlap_proc
// calling Init in order to provide these datas
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

// you can access to Nvol by calling GetNodlScalar
Nvol = A.GetNodlScalar();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetNbScalarUnknowns

Syntax

 int GetNbScalarUnknowns() const;

This method returns the local number of scalar components. Several components can be specified when Init is called.

Example :

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

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

// total number of rows (with all components)
int nglob = 18;

// global numbers (for all components)
IVect glob_number;

// overlapped rows (with original processor)
IVect OverlapRow, original_proc;

// local number of rows (with all components)
int n;

// local number of rows for only one component
int Nvol;

// number of components
int nb_u = 2;

// corresponding rows that are shared between processors
// for only one component
IVect ProcShared;
Vector<IVect> SharedRows;

if (MPI::COMM_WORLD.Get_rank() == 0)
  {
    // for example, processor 0 owns rows 1, 3, 4, 6, 7, 8 and 10, 12, 13, 15, 16, 17 
    n = 12;
    Nvol = 6;

    glob_number.Reallocate(n);
    glob_number(0) = 1; glob_number(1) = 3; glob_number(2) = 4;
    glob_number(3) = 6; glob_number(4) = 7; glob_number(5) = 8;
    glob_number(6) = 10; glob_number(7) = 12; glob_number(8) = 13;
    glob_number(9) = 15; glob_number(10) = 16; glob_number(11) = 17;

    // rows 3, 6 and 12, 15 are shared with processor 1
    // local numbers are 1, 3 and 7, 9
    // for ProcShared/SharedRows, only one component is needed
    ProcShared.Reallocate(1); SharedRows.Reallocate(1);
    ProcShared(0) = 1;
    SharedRows(0).PushBack(1); SharedRows(0).PushBack(3);

  }
else
  {
    // for example, processor 1 owns rows 0, 2, 3, 5, 6 and 9, 11, 12, 14, 15
    n = 10;
    Nvol = 5;

    glob_number.Reallocate(n);
    glob_number(0) = 0; glob_number(1) = 2; glob_number(2) = 3;
    glob_number(3) = 5; glob_number(4) = 6;
    glob_number(5) = 9; glob_number(6) = 11; glob_number(7) = 12;
    glob_number(8) = 14; glob_number(9) = 15;

    // rows 3, 6 and 12, 15 are originally in processor 0
    OverlapRow.Reallocate(4);
    // local numbers are 2, 4 and 7, 9
    OverlapRow(0) = 2; OverlapRow(1) = 4;
    OverlapRow(2) = 7; OverlapRow(3) = 9;
    original_proc.Reallocate(4); original_proc.Fill(0);

    // for ProcShared/SharedRows, only one component is needed
    ProcShared.Reallocate(1); SharedRows.Reallocate(1);
    ProcShared(0) = 0;
    SharedRows(0).PushBack(2); SharedRows(0).PushBack(4);
  }

// after filling nglob, global_row, overlap_row and overlap_proc
// calling Init in order to provide these datas
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

// you can access to nb_u by calling GetNbScalarUnknowns
Nvol = A.GetNbScalarUnknowns();

Location :

DistributedMatrix.hxx
DistributedMatrix.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 :

DistributedMatrix.hxx
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 :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetMaxDataSizeDistantCol

Syntax

 int GetMaxDataSizeDistantCol() const

This member function returns the number of non-zero entries in distant columns for all the processors (it is obtained with a MPI_Allreduce instruction).

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);

// in order to obtain the number of non-zero entries in distant columns
// (values that have been created with AddDistantInteraction)
// a matrix vector-product must be performed 
Mlt(A, x, y);

// number of non-zero entries stored for all the processors (only distant columns)
// the result is the same for each processor
int nnz_dist_col = A.GetMaxDataSizeDistantCol();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetMaxDataSizeDistantRow

Syntax

 int GetMaxDataSizeDistantRow() const

This member function returns the number of non-zero entries in distant rows for all the processors (it is obtained with a MPI_Allreduce instruction).

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, &amp;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);

// in order to obtain the number of non-zero entries in distant rows
// (values that have been created with AddRowDistantInteraction)
// a matrix vector-product must be performed
Mlt(A, x, y);

// number of non-zero entries stored for all the processors (only distant rows)
// the result is the same for each processor
int nnz_dist_row = A.GetMaxDataSizeDistantRow();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

IsReadyForMltAdd

Syntax

 bool IsReadyForMltAdd() const

This member function returns true if the structure is ready to perform a matrix-vector product. This function is not really useful since if the user performs a matrix-vector product with a structure that is not ready, the program will create all the arrays needed to perform the matrix-vector product. The function IsReadyForMltAdd can be used for informative purposes, for example to avoid taking into account the cost of constructing all these arrays (they induce MPI communications) if the structure was not ready.

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);

// the structure will be said "ready" for a matrix-vector operation
// if a first matrix-vector product has been performed
bool struct_ready = A.IsReadyForMltAdd();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetDistantColSize

Syntax

 int GetDistantColSize(int i) const

This member function returns the number of non-zero entries with distant column numbers for the local row i.

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);

// if you want to know all the non-zero entries stored in local row i (created wih AddDistantInteraction(i, ...) )
int nnz_row_i = A.GetDistantColSize(i);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

IndexGlobalCol

Syntax

 int IndexGlobalCol(int i, int j) const

This member function returns the global column number of the j-th non-zero entry with a distant column number for the local row i.

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);

// if you want to know all the non-zero entries stored in local row i (created wih AddDistantInteraction(i, ...) )
int nnz_row_i = A.GetDistantColSize(i);

// loop over these distant non-zero entries
for (int j = 0; j < nnz_row_i; j++)
  {
    // for the j-th non-zero entry
    // global column number 
    jglob = A.IndexGlobalCol(i, j);
    // distant processor
    proc = A.ProcessorDistantCol(i, j);
    // value
    val = A.ValueDistantCol(i, j);
  }

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ProcessorDistantCol

Syntax

 int ProcessorDistantCol(int i, int j) const

This member function returns the distant processor number of the j-th non-zero entry with a distant column number for the local row i.

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);

// if you want to know all the non-zero entries stored in local row i (created wih AddDistantInteraction(i, ...) )
int nnz_row_i = A.GetDistantColSize(i);

// loop over these distant non-zero entries
for (int j = 0; j < nnz_row_i; j++)
  {
    // for the j-th non-zero entry
    // global column number 
    jglob = A.IndexGlobalCol(i, j);
    // distant processor
    proc = A.ProcessorDistantCol(i, j);
    // value
    val = A.ValueDistantCol(i, j);
  }

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ValueDistantCol

Syntax

 T ValueDistantCol(int i, int j) const

This member function returns the value of the j-th non-zero entry with a distant column number for the local row i.

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);

// if you want to know all the non-zero entries stored in local row i (created wih AddDistantInteraction(i, ...) )
int nnz_row_i = A.GetDistantColSize(i);

// loop over these distant non-zero entries
for (int j = 0; j < nnz_row_i; j++)
  {
    // for the j-th non-zero entry
    // global column number 
    jglob = A.IndexGlobalCol(i, j);
    // distant processor
    proc = A.ProcessorDistantCol(i, j);
    // value
    val = A.ValueDistantCol(i, j);
  }

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetDistantRowSize

Syntax

 int GetDistantRowSize(int i) const

This member function returns the number of non-zero entries with distant row numbers for the local column i.

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 row is distant (ie located on another processor), 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(i)} is incremented with val
A.AddRowDistantInteraction(iglob, j, proc, val);

// if you want to know all the non-zero entries stored in local column j (created wih AddRowDistantInteraction(iglob, j, ...) )
int nnz_col_j = A.GetDistantRowSize(j);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

IndexGlobalRow

Syntax

 int IndexGlobalRow(int i, int j) const

This member function returns the global row number of the j-th non-zero entry with a distant row number for the local column i.

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 row is distant (ie located on another processor), you have to use AddRowDistantInteraction
// iglob : global row number
// i : 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.AddDistantInteraction(iglob, i, proc, val);

// if you want to know all the non-zero entries stored in local column i (created wih AddDistantInteraction(iglob, i, ...) )
int nnz_col_i = A.GetDistantRowSize(i);

// loop over these distant non-zero entries
for (int j = 0; j < nnz_col_i; j++)
  {
    // for the j-th non-zero entry
    // global row number 
    iglob = A.IndexGlobalRow(i, j);
    // distant processor
    proc = A.ProcessorDistantRow(i, j);
    // value
    val = A.ValueDistantRow(i, j);
  }

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ProcessorDistantRow

Syntax

 int ProcessorDistantRow(int i, int j) const

This member function returns the distant processor number of the j-th non-zero entry with a distant row number for the local column i.

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 row is distant (ie located on another processor), you have to use AddRowDistantInteraction
// iglob : global row number
// i : 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.AddDistantInteraction(iglob, i, proc, val);

// if you want to know all the non-zero entries stored in local column i (created wih AddDistantInteraction(iglob, i, ...) )
int nnz_col_i = A.GetDistantRowSize(i);

// loop over these distant non-zero entries
for (int j = 0; j < nnz_col_i; j++)
  {
    // for the j-th non-zero entry
    // global row number 
    iglob = A.IndexGlobalRow(i, j);
    // distant processor
    proc = A.ProcessorDistantRow(i, j);
    // value
    val = A.ValueDistantRow(i, j);
  }

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ValueDistantRow

Syntax

 T ValueDistantRow(int i, int j) const

This member function returns the value of the j-th non-zero entry with a distant row number for the local column i.

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 row is distant (ie located on another processor), you have to use AddRowDistantInteraction
// iglob : global row number
// i : 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.AddDistantInteraction(iglob, i, proc, val);

// if you want to know all the non-zero entries stored in local column i (created wih AddDistantInteraction(iglob, i, ...) )
int nnz_col_i = A.GetDistantRowSize(i);

// loop over these distant non-zero entries
for (int j = 0; j < nnz_col_i; j++)
  {
    // for the j-th non-zero entry
    // global row number 
    iglob = A.IndexGlobalRow(i, j);
    // distant processor
    proc = A.ProcessorDistantRow(i, j);
    // value
    val = A.ValueDistantRow(i, j);
  }

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Init

Syntax

 void Init(int nglob, IVect* global_row, IVect* overlap_row, IVect* overlap_proc,
           int Nvol, int nb_u, IVect* SharingProc, Vector<IVect>* SharingRows, MPI_Comm& comm);
 void Init(const DistributedMatrix_Base<T>&);
 void Init(const DistributedMatrixIntegerArray&);
 void Init(const IVect& global_row, const MPI_Comm& comm, const DistributedMatrix_Base<T>&);
 

The function Init must be called after the constructor in order to provide the global row numbers (array global_row), the rows already treated by another processor (array overlap_row), the processors that treat originally these rows (array overlap_proc) and the MPI communicator. Nvol is the number of rows for a scalar component, nb_u the number of scalar components, SharingProc stores the processors that share rows with current processor, and SharingRows the local row numbers that are shared with processor SharingProc(i). These numbers are assumed to be sorted such that shared row numbers with processor j in processor i correspond to row numbers shared with processor i in processor j. If the user does not know how to construct these arrays, they can be constructed automatically by providing the global row numbers only (last syntax of Init).

Example :

// on each processor, you specify rows that are already treated by another processor
// for example if the processor 0 handles rows [0, 3, 5, 6] and
// processor 1 the rows [1 2, 4, 5, 7], the row 5 is already treated by
// processor 0. Of course, if each row is treated by a processor and only one
// this array should be left empty
IVect OverlapRow;
int nglob = 8;
int n = 4;
if (MPI::COMM_WORLD.Get_rank() == 1)
  {
    n = 5;
    OverlapRow.Reallocate(1);
    // be careful because OverlapRow stores local numbers
    // here the global row 5 has a local number equal to 3 on processor 1
    OverlapRow(0) = 3;
  }

// for distributed matrices
// we also need to know the global row numbers
// and for each overlapped row, the original processor
IVect glob_number, original_proc;

// for shared rows, all the row numbers shared with processors
IVect SharingProc(1); Vector<IVect> SharingRows(1);
if (MPI::COMM_WORLD.Get_rank() == 0)
  {
    glob_number.Reallocate(4);
    glob_number(0) = 0;
    glob_number(1) = 3;
    glob_number(2) = 5;
    glob_number(3) = 6;

    // global row 5 is shared with processor 1 (local number = 2)
    SharingProc(0) = 1;
    SharingRows(0).PushBack(2);
  }
else
  {
    glob_number.Reallocate(5);
    glob_number(0) = 1;
    glob_number(1) = 2;
    glob_number(2) = 4;
    glob_number(3) = 5;
    glob_number(4) = 7;

    // row 5 is already handled by processor 0
    original_proc.Reallocate(1);
    original_proc(0) = 0;

    // here SharingRows is similar to OverlapRow because all shared dofs of this processor are already treated by processor 0
    // global row 5 is shared with processor 0 (local number = 3)
    SharingProc(0) = 0;
    SharingRows(0).PushBack(3);
  }

// for the constructor, you use the local number of rows
DistributedMatrix<double, General, ArrayRowSparse> A(n, n);

// then you need to initialise the arrays used for parallel functions (e.g. matrix-vector product)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       n, 1, &SharingProc, &SharingRows, MPI_COMM_WORLD);

// other solution : initialize from arrays already contained in another matrix
DistributedMatrix<double, General, ArrayRowSparse> B(n, n);

B.Init(A);

// last solution : encapsulate arrays in a DistributedMatrixIntegerArray
// and inits the matrix with this object
DistributedMatrixIntegerArray info;
info.GlobalRowNumbers = glob_number; info.nloc = n; info.nglob = nglob;
info.OverlapRowNumbers = OverlapRow; info.OverlapProcNumbers = original_proc;
info.ProcSharingRows = SharingProc; info.SharingRows ; SharingRows;
info.nodl_scalar = n; info.comm = MPI_COMM_WORLD; info.nb_unknowns_scal = 1;

DistributedMatrix<double, General, ArrayRowSparse> C(n, n);

C.Init(info);

Another example does not construct these arrays but calls Init to construct them automatically.

// local number of rows 
int n = 4;
if (MPI::COMM_WORLD.Get_rank() == 1)
  n = 5;

// each processor constructs its global numbers
// rows [0, 3, 5, 6] for processor 0
// rows [1 2, 4, 5, 7] for processor 1
// you notice that global row 5 is here shared by the two processors
IVect glob_number(n);
if (MPI::COMM_WORLD.Get_rank() == 0)
  {
    glob_number(0) = 0;
    glob_number(1) = 3;
    glob_number(2) = 5;
    glob_number(3) = 6;
  }
else
  {
    glob_number(0) = 1;
    glob_number(1) = 2;
    glob_number(2) = 4;
    glob_number(3) = 5;
    glob_number(4) = 7;
  }

// for the constructor, you use the local number of rows
DistributedMatrix<double, General, ArrayRowSparse> A(n, n);

// constructing object DistributedMatrixIntegerArray only from glob_number
DistributedMatrixIntegerArray info;
A.Init(glob_number, MPI_COMM_WORLD, info); // info is constructed here

// then you can reuse info for another matrix
DistributedMatrix<double, General, ArrayRowSparse> B(n, n);
B.Init(info); 

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Reallocate

Syntax

 void Reallocate(int m, int n);

This function reallocates the matrix, previous entries are lost.

Example :

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

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

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Resize

Syntax

 void Resize(int m, int n);

This function resizes the matrix, previous entries are kept.

Example :

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

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

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Clear

Syntax

 void Clear();

This function clears the matrix.

Example :

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

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

// the matrix is erased
A.Clear();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetGlobalRowNumber

Syntax

 IVect& GetGlobalRowNumber() const;

This function returns the global row numbers (local to global numbering).

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, &glob_number, &OverlapRow, &original_proc,
       n, 1, &SharingProc, &SharingRows, MPI_COMM_WORLD);

// you can access to global_row by calling GetGlobalRowNumber
IVect&amp; global_row = A.GetGlobalRowNumber();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetOverlapRowNumber

Syntax

 IVect& GetOverlapRowNumber() const;

This function returns the overlapped row numbers, i.e. rows that are already treated by another processor.

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, &glob_number, &OverlapRow, &original_proc,
       n, 1, &SharingProc, &SharingRows, MPI_COMM_WORLD);

// you can access to overlap_row by calling GetOverlapRowNumber
IVect&amp; overlap_row = A.GetOverlapRowNumber();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetOverlapProcNumber

Syntax

 IVect& GetOverlapProcNumber() const;

This function returns the processor numbers associated with overlapped row numbers, i.e. rows that are already treated by another processor.

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, &glob_number, &OverlapRow, &original_proc,
       n, 1, &SharingProc, &SharingRows, MPI_COMM_WORLD);

// you can access to original_proc by calling GetOverlapProcNumber
IVect&amp; original_proc = A.GetOverlaProcNumber();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetProcessorSharingRows

Syntax

 IVect& GetProcessorSharingRows() const;

This function returns the processor numbers that share rows with the current processor.

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, &glob_number, &OverlapRow, &original_proc,
       n, 1, &SharingProc, &SharingRows, MPI_COMM_WORLD);

// you can access to SharingProc by calling GetProcessorSharingRows
IVect& SharingProc = A.GetProcessorSharingRows();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetSharingRowNumbers

Syntax

 Vector<IVect>& GetSharingRowNumbers() const;

This function returns the row numbers that are shared with other processors.

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, &glob_number, &OverlapRow, &original_proc,
       n, 1, &SharingProc, &SharingRows, MPI_COMM_WORLD);

// you can access to SharingRows by calling GetSharingRowNumbers
IVect& SharingProc = A.GetSharingRowNumbers();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetMemorySize

Syntax

 size_t GetMemorySize() const;

This function returns the memory used (in bytes) to store the distributed matrix for the current matrix. The user should perform a MPI reduce operation to obtain the global size of the matrix.

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// once the matrix is filled, you can retrieve the size used for each processor
size_t taille = A.GetMemorySize();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetNonZeros

Syntax

 int GetNonZeros() const;

This function returns the number of non-zero entries stored in the distributed matrix for the current processor. The user should perform a MPI reduce operation to obtain the global number of non-zero entries.

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// once the matrix is filled, you can retrieve the number of non-zero entries
int nnz = A.GetNonZeros();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetDataSize

Syntax

 int GetDataSize() const;

This function returns the number of "equivalent" non-zero entries stored in the distributed matrix for the current processor. The user should perform a MPI reduce operation to obtain the global number of non-zero entries.

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// once the matrix is filled, you can retrieve the equivalent number of non-zero entries
int nnz = A.GetDataSize();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

RemoveSmallEntry

Syntax

 void RemoveSmallEntry(T epsilon) const;

This member function drops small values contained in the matrix. The result may depend on the number of processors, since non-zero entries are shared between processors. Here, each value stored is dropped if below epsilon, no assembling is performed.

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// you can drop small values if you want
A.RemoveSmallEntry(1e-12);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

SetIdentity

Syntax

 void SetIdentity();

This member function reinitializes the matrix as the identity matrix.

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// you can reinitialize A as the identity matrix
A.SetIdentity();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Fill

Syntax

  void Fill();
  void Fill(T x);

This member function fills the matrix with 0, 1, 2, etc or with a given value. The result may depend on the number of processors, since non-zero entries are shared between processors. Here each non-zero entry is modified without taking into account duplicate values due to shared row numbers.

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// you can initialize all non-zero entries to a given value
A.Fill(1.0);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

FillRand

Syntax

  void FillRand();

This member function fills the matrix with random values. The result may depend on the number of processors, since non-zero entries are shared between processors. Here each non-zero entry is modified without taking into account duplicate values due to shared row numbers.

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// you can initialize all non-zero entries to random values
A.FillRand();

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Write

Syntax

  void Write(string) const;
  void Write(ostream&) const;

This member function is not implemented for distributed matrices

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// you can write the matrix in files
A.Write("mat.dat");

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

WriteText

Syntax

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

This member function writes the matrix on files or output streams in text format. Each processor writes its part of the matrix on a different file (ending with _P0.dat, _P1.dat, etc).

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// you can write the matrix in files
// the matrix here will be written in files mat_P0.dat, mat_P1.dat, etc
A.WriteText("mat.dat");

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Read

Syntax

  void Read(string);
  void Read(istream&);

This member function is not implemented for distributed matrices

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ReadText

Syntax

  void ReadText(string);
  void ReadText(istream&);

This member function is not implemented for distributed matrices

Example :

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

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Mlt

Syntax

  void Mlt(const DistributedMatrix<T>& A, const Vector<T>& x, Vector<T>& y, bool assemble); 
  void Mlt(const SeldonTranspose& trans, const DistributedMatrix<T>& A, const Vector<T>& x, Vector<T>& y, bool assemble); 
  void Mlt(const T& alpha, DistributedMatrix<T>& A); 

This function can perform matrix-vector product (with the original matrix or with its transpose) and can be used to multiply the matrix by a scalar. The matrix-matrix product is currently not implemented. If the communicator contains only one processor, the sequential function Mlt (with the class Matrix) will be called, otherwise an error will be displayed during the execution.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// once the matrix is constructed, matrix-vector products can be performed
// each vector is stored as a usual Seldon vector (class Vector)
// each processor is assumed to store the values of the vector for the
// global rows given during the construction of the matrix (here, the array glob_number)
Vector<double> x(n), y(n), z;

// Mlt will assume that the values of x for rows that are shared between processors are the same
// which is here not the case when calling FillRand
x.FillRand();
y.FillRand();
z = y;

// a solution is to call AssembleVector to ensure this property
AssembleVector(x, MPI::SUM, ProcShared, SharedRows, MPI_COMM_WORLD, Nvol, nb_u, 23);

// classical matrix-vector product y = A x
// the default value of the last argument (assemble) is true
// the result y is assembled and will contain the same result as a sequential matrix-vector product,
// the values being distributed in the different processors
Mlt(A, x, y);

// assemble = true, induces a final call to the function AssembleVector (to add values of shared rows)
// if the user wants to perform a matrix-vector product without performing this final assembling step, he puts false in the last optional argument
Mlt(A, x, z, false);

// and perform the assembling phase later
AssembleVector(z, MPI::SUM, ProcShared, SharedRows, MPI_COMM_WORLD, Nvol, nb_u, 23);

// Mlt can be used to perform a matrix-vector with the transpose matrix or the transpose conjugate
Mlt(SeldonTrans, A, x, y);

// Finally Mlt can be used to multiply a distributed matrix by a scalar coefficient
Mlt(-2.1, A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

MltAdd

Syntax

  void MltAdd(const T& alpha, const DistributedMatrix<T>& A, const Vector<T>& x, const T& beta, Vector<T>& y, bool assemble); 
  void MltAdd(const T& alpha, const SeldonTranspose&, const DistributedMatrix<T>& A, const Vector<T>& x, const T& beta, Vector<T>& y, bool assemble); 

This function can perform matrix-vector product (with the original matrix or with its transpose). The matrix-matrix product is currently not implemented. If the communicator contains only one processor, the sequential function MltAdd (with the class Matrix) will be called, otherwise an error will be displayed during the execution.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// once the matrix is constructed, matrix-vector products can be performed
// each vector is stored as a usual Seldon vector (class Vector)
// each processor is assumed to store the values of the vector for the
// global rows given during the construction of the matrix (here, the array glob_number)
Vector<double> x(n), y(n), z;

// MltAdd will assume that the values of x for rows that are shared between processors are the same
// which is here not the case when calling FillRand
x.FillRand();
y.FillRand();
z = y;

// a solution is to call AssembleVector to ensure this property
AssembleVector(x, MPI::SUM, ProcShared, SharedRows, MPI_COMM_WORLD, Nvol, nb_u, 23);

double alpha = 2.5, beta = 0.4;

// matrix-vector product y = beta y + alpha A x
// the default value of the last argument (assemble) is true
// the result y is assembled and will contain the same result as a sequential matrix-vector product,
// the values being distributed in the different processors
MltAdd(alpha, A, x, beta, y);

// assemble = true, induces a final call to the function AssembleVector (to add values of shared rows)
// if the user wants to perform a matrix-vector product without performing this final assembling step, he puts false in the last optional argument
MltAdd(alpha, A, x, beta, z, false);

// and perform the assembling phase later
AssembleVector(z, MPI::SUM, ProcShared, SharedRows, MPI_COMM_WORLD, Nvol, nb_u, 23);

// MltAdd can be used to perform a matrix-vector with the transpose matrix or the transpose conjugate
MltAdd(alpha, SeldonTrans, A, x, beta, y);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Add

Syntax

  void Add(const T& alpha, const DistributedMatrix<T>& A, DistributedMatrix<T>& B); 

This function adds a distributed matrix multiplied by a scalar to another distributed matrix. The two matrices are assumed to have the same global row numbers.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// a way to have the same global row numbers
// is to use Init with the matrix A
DistributedMatrix<double, General, ArrayRowSparse> B;

B.Init(A);

double alpha = 2.1;
// we compute B = B + alpha A
Add(alpha, A, B);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

MaxAbs

Syntax

  Treal MaxAbs(const DistributedMatrix<T>& A);

This function returns the maximum absolute value of non-zero entries of a distributed matrix A. Currently the implemented function is not equivalent to the sequential function because of shared rows. It can be used to have an estimation of the actual maximum absolute value, but it will not return the exact value.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we estimate the maximum non-zero entry (modulus)
double maxA = MaxAbs(A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetRowSum

Syntax

  void GetRowSum(Vector<Treal>& sumA, const DistributedMatrix<T>& A);

This function computes in sumA the sum of absolute values of a_ij for each row i. Currently the implemented function is not equivalent to the sequential function because of shared rows. It can be used to have an estimation of the exact sum.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we estimate \sum_j |a_{i,j}| for each row i
Vector<double> sumA;
GetRowSum(sumA, A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetColSum

Syntax

  void GetColSum(Vector<Treal>& sumA, const DistributedMatrix<T>& A);

This function computes in sumA the sum of absolute values of a_ij for each column j. Currently the implemented function is not equivalent to the sequential function because of shared rows. It can be used to have an estimation of the exact sum.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we estimate \sum_i |a_{i,j}| for each column j
Vector<double> sumA;
GetColSum(sumA, A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetRowColSum

Syntax

  void GetRowColSum(Vector<Treal>& rowA, Vector<Treal>& colA, const DistributedMatrix<T>& A);

This function computes in rowA the sum of absolute values of a_ij for each row i and in colA the sum of absolute values for each column j. Currently the implemented function is not equivalent to the sequential function because of shared rows. It can be used to have an estimation of the exact sum.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we estimate \sum |a_{i,j}| for each row i and column j
Vector<double> rowA, colA;
GetRowColSum(rowA, colA, A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Norm1

Syntax

  Treal Norm1(const DistributedMatrix<T>& A);

This function returns an estimation of the 1-norm of a distributed matrix A. Currently the implemented function is not equivalent to the sequential function because of shared rows. It can be used to have an estimation of the actual 1-norm, but it will not return the exact value.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we estimate the 1-norm
double maxA = Norm1(A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

NormInf

Syntax

  Treal NormInf(const DistributedMatrix<T>& A);

This function returns an estimation of the infinite norm of a distributed matrix A. Currently the implemented function is not equivalent to the sequential function because of shared rows. It can be used to have an estimation of the actual infinite norm, but it will not return the exact value.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we estimate the infinite norm
double maxA = NormInf(A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Transpose

Syntax

  void Transpose(DistributedMatrix<T>& A);
  void Transpose(const DistributedMatrix<T>& A, DistributedMatrix<T>& B);

This function transposes a distributed matrix.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// A can be replaced by its transpose :
Transpose(A);

// the transpose can also be computed into another matrix
DistributedMatrix<double, General, ArrayRowSparse> B;

Transpose(A, B);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Conjugate

Syntax

  void Conjugate(DistributedMatrix<T>& A);

This function conjugates a distributed matrix.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// A can be replaced by its conjugate:
Conjugate(A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

TransposeConj

Syntax

  void TransposeConj(DistributedMatrix<T>& A);
  void TransposeConj(const DistributedMatrix<T>& A, DistributedMatrix<T>& B);

This function transposes and conjugates a distributed matrix.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// A can be replaced by its transpose conjugate :
TransposeConj(A);

// the transpose conjugate can also be computed into another matrix
DistributedMatrix<double, General, ArrayRowSparse> B;

TransposeConj(A, B);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetRow

Syntax

  void GetRow(const DistributedMatrix<T>& A, int, Vector<T, VectSparse>&);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetCol

Syntax

  void GetCol(const DistributedMatrix<T>& A, int, Vector<T, VectSparse>&);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

SetRow

Syntax

  void SetRow(const Vector<T, VectSparse>&, int, DistributedMatrix<T>& A);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

SetCol

Syntax

  void SetCol(const Vector<T, VectSparse>&, int, DistributedMatrix<T>& A);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ApplyPermutation

Syntax

  void ApplyPermutation(DistributedMatrix<T>& A, const IVect& row_perm, const IVect& col_perm);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ApplyInversePermutation

Syntax

  void ApplyInversePermutation(DistributedMatrix<T>& A, const IVect& row_perm, const IVect& col_perm);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

SOR

Syntax

  void SOR(const DistributedMatrix<T>& A, Vector<T>& X, const Vector<T>& B, const T& omega, int iter, int type_ssor);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

Copy

Syntax

  void Copy(const DistributedMatrix<T>& A, DistributedMatrix<T>&);

This function converts a distributed sparse matrix into another one. The format can be different (RowSparse, ColSparse, ArrayRowSparse, etc), the local part of the matrix is converted by calling the appropriate Copy function (located in the file Matrix_Conversions.cxx), the distributed part is always stored with the same format.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// The target matrix just needs to be declared (Init is not needed since it will take the arrays provided by the original matrix A)
DistributedMatrix<double, General, ColSparse> B;

// it should be emphasized that B has a storage ColSparse and A ArrayRowSparse
// Copy will perform a conversion
Copy(A, B);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

GetSubMatrix

Syntax

  void GetSubMatrix(const DistributedMatrix<T>& A, int m, int n, DistributedMatrix<T>& B);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

NormFro

Syntax

  Treal NormFro(const DistributedMatrix<T>& A);

This function is not implemented for distributed matrices. If the communicator has only one processor, the sequential function is called, otherwise an error is displayed.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ScaleLeftMatrix

Syntax

  void ScaleLeftMatrix(DistributedMatrix<T>& A, const Vector<T>& scale);

This function scales each row of the matrix by a given coefficient.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we get the sum of absolute values for each row
Vector<double> scale(n);
GetRowSum(scale, A);

// here the values of scale for shared rows should be the same
// this is case for the result of GetRowSum
for (int i = 0; i < n; i++)
  scale(i) = 1.0/scale(i);

// the row i of the matrix is multiplied by scale(i)
ScaleLeftMatrix(A, scale);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ScaleRightMatrix

Syntax

  void ScaleRightMatrix(DistributedMatrix<T>& A, const Vector<T>& scale);

This function scales each column of the matrix by a given coefficient.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we get the sum of absolute values for each column
Vector<double> scale(n);
GetColSum(scale, A);

// here the values of scale for shared rows should be the same
// this is case for the result of GetColSum
for (int i = 0; i < n; i++)
  scale(i) = 1.0/scale(i);

// the column i of the matrix is multiplied by scale(i)
ScaleRightMatrix(A, scale);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

ScaleMatrix

Syntax

  void ScaleMatrix(DistributedMatrix<T>& A, const Vector<T>& scale_row, const Vector<T>& scale_col);

This function scales each row and column of the matrix by a given coefficient.

Example :

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

// global row numbers are provided with Init (glob_number)
A.Init(nglob, &glob_number, &OverlapRow, &original_proc,
       Nvol, nb_u, &ProcShared, &SharedRows, MPI_COMM_WORLD);

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

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// we get the sum of absolute values for each row and column
Vector<double> scale_row(n), scale_col(n);
GetRowColSum(scale_row, scale_col, A);

// here the values of scale for shared rows should be the same
// this is case for the result of GetRowColSum
for (int i = 0; i < n; i++)
  {
    scale_row(i) = 1.0/scale_row(i);
    scale_col(i) = 1.0/scale_col(i);
  }

// a_ij is multiplied by scale_row(i) scale_col(j)
ScaleMatrix(A, scale_row, scale_col);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

AssembleDistributed

Syntax

 void AssembleDistributed(DistributedMatrix& A, Prop& prop, const MPI_Comm& comm, 
                          IVect& col_numbers, IVect& local_col_numbers,
                          IVect& PtrA, IVect& IndA, Vector& ValA, bool sym_pattern, bool reorder);

This function assembles a distributed matrix A such that each processor has only a row (or column). col_number contains the global row numbers, local_col_numbers contains the local row numbers of the original matrix A that are kept. Values associated with overlapped rows and distant columns are sent to the original processor. The columns of the resulting distributed matrices are given in CSC format, PtrA contains the beginning indexes for each column, IndA the global row numbers, ValA the values. These arrays can be used directly by direct solvers. This function is called in the class DistributedSolver which interfaces distributed direct solvers in Seldon (such as Pastix or Mumps). If sym_pattern is true, the pattern of the resulting matrix is symmetrized. If reorder is true, the row numbers are reordered such that row numbers of processor 0 are 0, 1, ..., N0-1, row numbers of processor 1 are N0, N0+1, ..., N0+N1-1. By default, reorder is set to false. reorder are false.

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, MPI_COMM_WORLD);

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// you can assemble the matrix and retrieve columns of the matrix
General prop;
// if sym_pattern is true, the pattern is symmetrized (by adding 0-values if necessary), the values being unchanged
bool sym_pattern = false;

IVect col_numbers, local_col_numbers;
IVect PtrA, IndA; Vector<double> ValA;
// the function fills the following arrays :
// col_numbers : global numbers of columns stored by the current processor (each column is owned by only one processor)
// local_col_numbers : local numbers of global columns stored
// PtrA, IndA, ValA : columns stored in CSC format (beginning of indexes, global row numbers and values) 
AssembleDistributed(A, prop, A.GetCommunicator(),
                    col_numbers, local_col_numbers,
                    PtrA, IndA, ValA, sym_pattern);

// for symmetric matrices, the upper part of the matrix is assembled, by rows :
DistributedMatrix<double, Symmetric, ArrayRowSymSparse> Asym;
Symmetric prop_sym;
// in these case, col_numbers, local_col_numbers contains row numbers (and not columns)
// and matrix PtrA, IndA, ValA are rows stored in CSR format
AssembleDistributed(Asym, prop_sym, Asym.GetCommunicator(),
                    col_numbers, local_col_numbers,
                    PtrA, IndA, ValA, sym_pattern);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

EraseRow

Syntax

 void EraseRow(IVect& num, DistributedMatrix& A);

This function erases some rows (local row numbers are provided in num) of the matrix A.

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, MPI_COMM_WORLD);

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// then you can erase some rows
IVect num(2);
EraseRow(num, A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

EraseCol

Syntax

 void EraseCol(IVect& num, DistributedMatrix& A);

This function erases some columns (local column numbers are provided in num) of the matrix A.

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, MPI_COMM_WORLD);

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

// when the column is distant, you can use AddDistantInteraction
A.AddDistantInteraction(i, jglob, proc, val);

// when the row is distant, you can use AddRowDistantInteraction
A.AddRowDistantInteraction(iglob, j, proc, val);

// then you can erase some columns
IVect num(2);
EraseCol(num, A);

Location :

DistributedMatrix.hxx
DistributedMatrix.cxx

OverlapRowNumbers, OverlapProcNumbers

Syntax

 Vector<int> OverlapRowNumbers;
 Vector<int> OverlapProcNumbers;

The attribute OverlapRowNumbers contains local row numbers that are already treated by another processor. The processor number treating the row is stored in OverlapProcNumbers.

For instance, you can have the following distribution of rows:

Proc 0 Proc 1 Proc 2
Row numbers 0
2
4
1
2
5
6
0
3
4
6

We obtain the following values for OverlapRowNumbers:

Proc 0 Proc 1 Proc 2
OverlapRowNumbers None 1 0
2
3

and for OverlapProcNumbers:

Proc 0 Proc 1 Proc 2
OverlapProcNumbers None 0 0
0
1

Location :

DistributedMatrix.hxx

GlobalRowNumbers

Syntax

 Vector<int> GlobalRowNumbers;

The attribute GlobalRowNumbers contains global row numbers. Some rows can be shared by different processors. This array contains sorted numbers.

For instance, you can have the following distribution of rows:

Proc 0 Proc 1 Proc 2
GlobalRowNumbers 0
2
4
1
2
5
6
0
3
4
6

Location :

DistributedMatrix.hxx

SharingRowNumbers, ProcSharingRows

Syntax

 Vector<int> ProcSharingRows;
 Vector<Vector<int> > SharingRowNumbers;

The attribute ProcSharingsRows contains the numbers of processors that shares rows with the current processor. The array SharingRowNumbers(i) will contain the list of local rows shared with processor ProcSharingRows(i). The arrays are sorted such that rows correspond between processors.

For instance, you can have the following distribution of rows:

Proc 0 Proc 1 Proc 2
Row numbers 0
2
4
1
2
5
6
0
3
4
6

We obtain the following values for ProcSharingRows:

Proc 0 Proc 1 Proc 2
ProcSharingRows 1
2
0
2
0
1

and for SharingRowNumbers:

Proc 0 Proc 1 Proc 2
SharingRowNumbers (1)
(0, 2)
(1)
(3)
(0, 2)
(3)

Location :

DistributedMatrix.hxx

nodl_scalar, nb_unknowns

Syntax

  int nodl_scalar;
  int nb_unknowns;

The attribute nodl_scalar contains the number of rows for a scalar unknown. The attribute nb_unknowns is the number of scalar unknowns. By default, it is equal to 1, and nodl_scalar is equal to nloc. Actually this number is only relevant for the array SharingRowNumbers. You can reduce the memory used by this array if the rows can be numbered as

k = j*(nodl_scalar) + i

where j is an unknown number and i the dof (degree of freedom) number. nodl_scalar is the local number of degrees of freedom for the current processor. In that case, you can provide a nb_unknowns equal to the number of unknowns (j will go from 0 to nb_unknowns-1), and SharingRowNumbers will contain only the dof numbers (i). It is assumed that each unknown has the same distribution pattern over processors.

Location :

DistributedMatrix.hxx

nglob, nloc

Syntax

  int nloc;
  int nglob;

The attribute nloc is the local number of rows. The attribute nglob is the global number of rows. For instance, if you have the following distribution of rows:

Proc 0 Proc 1 Proc 2
Row numbers 0
2
4
1
2
5
6
7
0
3
4
6

nglob will be equal to 8 and nloc will be equal to 3, 5 and 4 for the processor 0, 1 and 2 respectively.

Location :

DistributedMatrix.hxx

comm

Syntax

  MPI_Comm comm;

The attribute comm is the MPI communicator.

Location :

DistributedMatrix.hxx

Clear

Syntax

  void Clear();

This method clears all the arrays stored in the object.

Location :

DistributedMatrix.hxx

SetData

Syntax

  void SetData(int nloc, int nglob, int nodl, int nb_unknowns, MPI_Comm comm,
               IVect& overlap_row, IVect& overlap_proc, IVect& global_rows,
               IVect& proc_rows, Vector<IVect>& sharing_rows);

This method fills the array of the object with provided arguments. It is up to the user to nullify these objects (or not) to avoid segmentation faults. This method is a low level method and should be used carefully.

Example:

int nloc =2, nglob = 10, nodl = nloc, nb_unknowns = 1;
MPI_Comm comm = MPI_COMM_WORLD;
// you fill arrays overlap_row, overlap_proc, etc as you wish
IVect overlap_row, overlap_proc, global_rows;
IVect proc_rows; Vector<IVect> sharing_rows;

// then you can provide these arrays to the object
DistributedMatrixIntegerArray info;
info.SetData(nloc, nglob, nodl, nb_unknowns, comm, overlap_row, overlap_proc,
            global_rows, proc_rows, sharing_rows);

// and then nullify to avoid segmentation faults
overlap_row.Nullify(); overlap_proc.Nullify(); global_rows.Nullify();
proc_rows.Nullify(); sharing_rows.Nullify();

Location :

DistributedMatrix.hxx

Nullify

Syntax

  void Nullify();

This method nullifies the array of the object. This method is a low level method and should be used carefully.

Example:

int nloc =2, nglob = 10, nodl = nloc, nb_unknowns = 1;
MPI_Comm comm = MPI_COMM_WORLD;
// you fill arrays overlap_row, overlap_proc, etc as you wish
IVect overlap_row, overlap_proc, global_rows;
IVect proc_rows; Vector<IVect> sharing_rows;

// then you can provide these arrays to the object
DistributedMatrixIntegerArray info;
info.SetData(nloc, nglob, nodl, nb_unknowns, comm, overlap_row, overlap_proc,
            global_rows, proc_rows, sharing_rows);

// when you are done with the object, you can nullify it to avoir segmentation fault
info.Nullify();

Location :

DistributedMatrix.hxx

ConstructArrays

Syntax

  void ConstructArrays(const IVect& global_row_numbers, IVect& overlap_num,
                   IVect& proc_num, IVect& sharing_proc, Vector<IVect>& sharing_rows, MPI_Comm comm);
  void ConstructArrays(const Vector<IVect>& all_rows, IVect& global_row_numbers, IVect& overlap_num,
                   IVect& proc_num, IVect& sharing_proc, Vector<IVect>& sharing_rows, MPI_Comm comm, bool distribute = true);

These static methods construct the arrays given in argument from the global row numbers. The second syntax can be used if the root processor (rank = 0) has already access to the global row numbers for all the processors (input array all_rows). The last optional argument is true if you want to distribute these numbers to the local arrays global_row_numbers.

Example:

// assuming that proc 0 has rows (0, 2, 5, 8), proc 1 has rows (1, 2, 3, 5, 7, 9), and proc 2 has rows (0, 4, 5, 6, 8)
DistributedMatrixIntegerArray info;

// first solution : you fill these global row numbers for each proc
Vector<int> row_num;
if (MPI::COMM_WORLD.Get_rank() == 0)
{
  row_num.Reallocate(4);
  row_num(0) = 0; row_num(1) = 2; row_num(2) = 5; row_num(3) = 8;
}
else if (MPI::COMM_WORLD.Get_rank() == 1)
{
  row_num.Reallocate(6);
  row_num(0) = 1; row_num(1) = 2; row_num(2) = 3; row_num(3) = 5; row_num(4) = 7; row_num(5) = 9;
}
else
{
  row_num.Reallocate(5);
  row_num(0) = 0; row_num(1) = 4; row_num(2) = 5; row_num(3) = 6; row_num(4) = 8;
}

// and compute other arrays by calling ConstructArrays :
info.GlobalRowNumbers = row_num; info.nglob = 10; info.comm = MPI_COMM_WORLD;
info.nloc = row_num.GetM(); info.nodl_scalar = info.nloc; info.nb_unknowns_scal = 1;
info.ConstructArrays(row_num, info.OverlapRowNumbers, info.OverlapProcNumbers,
                     info.ProcSharingRows, info.SharingRowNumbers, info.comm);


// second solution : the global numbers are filled only on root
Vector<IVect> all_rows(3);
if (MPI::COMM_WORLD.Get_rank() == 0)
{
  all_rows(0).Reallocate(4);
  all_rows(0)(0) = 0; all_rows(0)(1) = 2; all_rows(0)(2) = 5; all_rows(0)(3) = 8;

  all_rows(1).Reallocate(6);
  all_rows(1)(0) = 1; all_rows(1)(1) = 2; all_rows(1)(2) = 3; all_rows(1)(3) = 5; all_rows(1)(4) = 7; all_rows(1)(5) = 9;

  all_rows(2).Reallocate(5);
  all_rows(2)(0) = 0; all_rows(2)(1) = 4; all_rows(2)(2) = 5; all_rows(2)(3) = 6; all_rows(2)(4) = 8;
}

// info.GlobalRowNumbers is computed by distributing all_rows, and other arrays are constructed :
info.nglob = 10; info.comm = MPI_COMM_WORLD;
info.ConstructArrays(all_rows, info.GlobalRowNumbers, info.OverlapRowNumbers, info.OverlapProcNumbers,
                     info.ProcSharingRows, info.SharingRowNumbers, info.comm);
                     
info.nloc = info.GlobalRowNumbers.GetM(); info.nodl_scalar = info.nloc; info.nb_unknowns_scal = 1;

Location :

DistributedMatrix.hxx