DSMat

This class is used to represent Distributed Sparse Matrix. More…

#include "DSMat.h"

Inherits from MatrixProd

Public Functions

Name
void copy_DSMat(DSMat & other)
Copy members.
DSMat & operator=(DSMat & other)
Copy operator.
DSMat()
Creates an empty object.
~DSMat()
Deletes the object.
void CreateCompactCSR(CSRMAT< iReg > & mat_out, iReg & shift_diag) const
Creates a compacted CSR matrix from a stripe of DSMat.
void Load_CompactCSR(const CSRMAT< iReg > & mat_in)
This function loads a compacted CSR matrix into a stripe of DSMat.
void PartKway(const type_MPI_iReg rank, const type_MPI_iReg size, VEC< iGlo > & vec_part, VEC< iGlo > & vec_iperm)
Redistributes the entries currently stored locally into this process to all the other process according to a partition generated by metis.
void Redist(const type_MPI_iReg rank, const type_MPI_iReg size, const VEC< iGlo > & vec_part, const VEC< iGlo > & vec_iperm, DSMat & mat_AP)
Creates through metis a partition of the DSMat instance.
void LocalRCM(VEC< iReg > & vec_perm, VEC< iReg > & vec_iperm)
Computes a local permutation (i.e. each subdomain has a local independent ordering) on an input DSMat using the reverse Cuthill-McKee (RCM) algorithm.
void LocalPermute(const type_MPI_iReg rank, const type_MPI_iReg size, const bool Inverse_Flag, const VEC< iReg > vec_perm_in, const VEC< iReg > vec_iperm_in, DSMat & mat_AP)
Performs a local permutation (i.e. each subdomain has a local independent ordering) on an input DSMat. An input flag is provided to decide whether to produce a direct or an inverse ordering.
void mk_OutptDom()
Creates OutptDom from the Outgoing communication graph.
void evaluate_and_store_maxsize_comm()
Evaluate maximum number of incoming unknowns of the local csr matrices and store the result in maxsize_comm.
void set_sym_flag(bool input_sym_flag)
Set the symmetric flag.
void set_Prepared_flag(bool input_Prepared_flag)
Set flag for MxV prepare.
void set_nrows(iReg input_nrows, iReg input_nrows_i)
Store nrows and nrows_i.
matInfo check()
This function checks some basic properties of the DSMAT and returns this information in the ad hoc structure info.
bool get_sym_flag() const
Get symmetric flag.
virtual iReg get_nrows() const
Get nrows.
iReg get_nrows_i() const
Get nrows_i.
iReg get_nterm() const
Get local nterm (total number of terms in the process stripe).
iReg get_maxsize_comm() const
Get maximum number of unknowns among all the local CSR matrix (maxsize_comm).
iReg * get_ptr_ptDom_data()
Get pointer to ptDom.
const iReg * get_ptr_ptDom_data() const
iReg * get_ptr_OutptDom_data()
Get OutptDom.
const iExt * get_ptr_comm_index_L_RcvLst_data() const
Get pointer to comm_index_L_RcvLst.
const iExt * get_ptr_comm_index_L_SndLst_data() const
Get pointer to comm_index_L_SndLst.
const iExt * get_ptr_comm_index_R_RcvLst_data() const
Get pointer to comm_index_R_RcvLst.
const iExt * get_ptr_comm_index_R_SndLst_data() const
Get pointer to comm_index_R_SndLst.
rExt * get_ptr_WL_rcv_data()
Get pointer to WL_rcv.
rExt * get_ptr_WL_snd_data()
Get pointer to WL_snd.
rExt * get_ptr_WR_rcv_data()
Get pointer to WR_rcv.
rExt * get_ptr_WR_snd_data()
Get pointer to WR_snd.
MPI_Request * get_ptr_Recv_req_data()
Get pointer to Recv_req.
MPI_Request * get_ptr_Send_req_data()
Get pointer to Send_req.
const iReg * get_ptr_L_jc_RcvLst_data() const
Get pointer to L_jc_RcvLst.
const iReg * get_ptr_L_jc_SndLst_data() const
Get pointer to L_jc_SndLst.
iReg get_L_jc_SndLst_size() const
Get size of L_jc_SndLst.
const iReg * get_ptr_R_jc_RcvLst_data() const
Get pointer to R_jc_RcvLst.
const iReg * get_ptr_R_jc_SndLst_data() const
Get pointer to R_jc_SndLst.
iReg get_R_jc_SndLst_size() const
Get size of R_jc_SndLst.
iReg * get_ptr_in_list_com_scr_data()
Get global incoming comunication list.
iReg * get_ptr_out_list_com_scr_data()
Get global outgoing comunication list.
iGlo Get_GlobalNNZ() const
Evaluates the total number of nnz of the distributed csr matrices.
virtual iGlo Get_Globalnrows() const
rExt Get_GlobalMaxRowNorm() const
Evaluates the maximum value of the norm of a row of the matrix.
struct LocMax Get_GlobalLocMaxRowNorm() const
Evaluates the maximum value and location of the norm of a row of the matrix.
rExt Get_GlobalAvgRowNorm() const
Evaluates the average value of the norm of a row of the matrix.
CSRMAT< iExt > Convert_arrayCSR_To_ScalarLRGCSR(iReg * new_order =nullptr) const
Creates a CSRMAT starting from the array of local csr matrices. Useful in printing (used by fprintAll_PARDDMAT).
void cpt_matBalance(const type_MPI_iReg rank, const type_MPI_iReg size, struct bal_info & info)
Function to compute the load balance of current matrix. The result is stored in process with rank 0.
void print_matBalance(const type_MPI_iReg rank, const struct bal_info info)
Function to print the load balance of current matrix.
void RowStats(const iReg nb, const rExt * UpBounds, rExt * Stats) const
Function to compute row statistics of the DSMat matrix. The norm of the rows will be divided in the user provided classes and written in the Stats array. The array UpBound defines the upper bound of each rownorm category, i.e. if the row norm is below to UpBound[i], it belongs to the class i.
void fPrintAll_PARDDMAT_ASCII(const string & filename) const
Prints the distributed dense matrix in ASCII filename (overwrite if it exists). It uses Convert_arrayCSR_To_ScalarLRGCSR and CreateGlobalOffset.
void fPrintAll_PARDDMAT_BIN(const string & filename) const
Prints the distributed dense matrix in BINARY filename (overwrite if it exists). It uses Convert_arrayCSR_To_ScalarLRGCSR and CreateGlobalOffset.
virtual void Prepare_MxV(const iReg nRHS)
Prepares the data structure necessary for MxV. L/R_SndRcvLst, Send/Recv_req,Recv_flag,WL_rcv,WR_snd,WR_rcv,WR_snd,WL_snd.
virtual void UndoPrepare_MxV()
Undo the data structure necessary for MxV.
virtual void MxV(DDMat &restrict x, DDMat &restrict b, bool Barrier_flag)
Computes the Matrix by vector Product.
void MxM(DSMat &restrict B, DSMat &restrict C)
Computes the Matrix by Matrix Product A x B = C (with A equal to this).
void TreatBC()
void Accumulate(const rExt fa, const rExt fb, DSMat &restrict B)
Accumulate a matrix B over A (with ‘this’ == A) such that A = fa x A + fb x B.
void Transpose(DSMat &restrict MatT)
Transpose the matrix.
void SparseToDense(DDMat &restrict DenseMat)
Compacts the DSMat into a DDMat.
void RemoveEmptyCSRs()
Remove empty CSRs.

Public Attributes

Name
CommGraph * MatGraph
Incoming (for MxV and MxM) Communication graph.
CommGraph * OutMatGraph
Outgoing (for MxV and MxM) Communication graph If the matrix is symmetric OutMatGraph just points to the same MatGraph data.
VEC< iReg > * ptDom
vector of pointers to the compacted column indices of each incoming block. That is, ptDom[ib+1] - ptDom[ib] gives the number of columns of the incoming block ib (ib block on this row of blocks).
VEC< iReg > * OutptDom
transposed local compact vector to the unknowns
VEC< CSRMAT< iReg > > csr
vector of CSR.
CSRMAT< iReg > * matR
array to the Right CSR: each of this will be a reference to the previous vector of CSR (no memory associated)
CSRMAT< iReg > * matL
array to the left CSR: each of this will be a reference to a CSR mat (no memory associated)
CSRMAT< iReg > * matC
pointer to the diagonal entry of the vector csr (no memory associated)

Additional inherited members

Public Functions inherited from MatrixProd

Name
void copy_MatrixProd(MatrixProd & other)
Copy members.
bool get_Prepared_flag() const
Get flag for MxV prepare.
iReg get_nRHS_prep() const
Get the number of rhs the matrix is prepared for.
virtual ~MatrixProd() =0
Deletes the object.

Protected Attributes inherited from MatrixProd

Name
bool Prepared_flag
flag for MxV prepare.
iReg nRHS_prep
number of rhs the matrix is prepared for

Detailed Description

class DSMat;

This class is used to represent Distributed Sparse Matrix.

class DSMat.

Public Functions Documentation

function copy_DSMat

inline void copy_DSMat(
    DSMat & other
)

Copy members.

function operator=

inline DSMat & operator=(
    DSMat & other
)

Copy operator.

function DSMat

DSMat()

Creates an empty object.

function ~DSMat

~DSMat()

Deletes the object.

function CreateCompactCSR

void CreateCompactCSR(
    CSRMAT< iReg > & mat_out,
    iReg & shift_diag
) const

Creates a compacted CSR matrix from a stripe of DSMat.

function Load_CompactCSR

void Load_CompactCSR(
    const CSRMAT< iReg > & mat_in
)

This function loads a compacted CSR matrix into a stripe of DSMat.

function PartKway

void PartKway(
    const type_MPI_iReg rank,
    const type_MPI_iReg size,
    VEC< iGlo > & vec_part,
    VEC< iGlo > & vec_iperm
)

Redistributes the entries currently stored locally into this process to all the other process according to a partition generated by metis.

function Redist

void Redist(
    const type_MPI_iReg rank,
    const type_MPI_iReg size,
    const VEC< iGlo > & vec_part,
    const VEC< iGlo > & vec_iperm,
    DSMat & mat_AP
)

Creates through metis a partition of the DSMat instance.

function LocalRCM

void LocalRCM(
    VEC< iReg > & vec_perm,
    VEC< iReg > & vec_iperm
)

Computes a local permutation (i.e. each subdomain has a local independent ordering) on an input DSMat using the reverse Cuthill-McKee (RCM) algorithm.

function LocalPermute

void LocalPermute(
    const type_MPI_iReg rank,
    const type_MPI_iReg size,
    const bool Inverse_Flag,
    const VEC< iReg > vec_perm_in,
    const VEC< iReg > vec_iperm_in,
    DSMat & mat_AP
)

Performs a local permutation (i.e. each subdomain has a local independent ordering) on an input DSMat. An input flag is provided to decide whether to produce a direct or an inverse ordering.

function mk_OutptDom

void mk_OutptDom()

Creates OutptDom from the Outgoing communication graph.

function evaluate_and_store_maxsize_comm

void evaluate_and_store_maxsize_comm()

Evaluate maximum number of incoming unknowns of the local csr matrices and store the result in maxsize_comm.

function set_sym_flag

inline void set_sym_flag(
    bool input_sym_flag
)

Set the symmetric flag.

function set_Prepared_flag

inline void set_Prepared_flag(
    bool input_Prepared_flag
)

Set flag for MxV prepare.

function set_nrows

inline void set_nrows(
    iReg input_nrows,
    iReg input_nrows_i
)

Store nrows and nrows_i.

function check

matInfo check()

This function checks some basic properties of the DSMAT and returns this information in the ad hoc structure info.

function get_sym_flag

inline bool get_sym_flag() const

Get symmetric flag.

function get_nrows

inline virtual iReg get_nrows() const

Get nrows.

Reimplements: MatrixProd::get_nrows

function get_nrows_i

inline iReg get_nrows_i() const

Get nrows_i.

function get_nterm

inline iReg get_nterm() const

Get local nterm (total number of terms in the process stripe).

function get_maxsize_comm

inline iReg get_maxsize_comm() const

Get maximum number of unknowns among all the local CSR matrix (maxsize_comm).

function get_ptr_ptDom_data

inline iReg * get_ptr_ptDom_data()

Get pointer to ptDom.

function get_ptr_ptDom_data

inline const iReg * get_ptr_ptDom_data() const

function get_ptr_OutptDom_data

inline iReg * get_ptr_OutptDom_data()

Get OutptDom.

function get_ptr_comm_index_L_RcvLst_data

inline const iExt * get_ptr_comm_index_L_RcvLst_data() const

Get pointer to comm_index_L_RcvLst.

function get_ptr_comm_index_L_SndLst_data

inline const iExt * get_ptr_comm_index_L_SndLst_data() const

Get pointer to comm_index_L_SndLst.

function get_ptr_comm_index_R_RcvLst_data

inline const iExt * get_ptr_comm_index_R_RcvLst_data() const

Get pointer to comm_index_R_RcvLst.

function get_ptr_comm_index_R_SndLst_data

inline const iExt * get_ptr_comm_index_R_SndLst_data() const

Get pointer to comm_index_R_SndLst.

function get_ptr_WL_rcv_data

inline rExt * get_ptr_WL_rcv_data()

Get pointer to WL_rcv.

function get_ptr_WL_snd_data

inline rExt * get_ptr_WL_snd_data()

Get pointer to WL_snd.

function get_ptr_WR_rcv_data

inline rExt * get_ptr_WR_rcv_data()

Get pointer to WR_rcv.

function get_ptr_WR_snd_data

inline rExt * get_ptr_WR_snd_data()

Get pointer to WR_snd.

function get_ptr_Recv_req_data

inline MPI_Request * get_ptr_Recv_req_data()

Get pointer to Recv_req.

function get_ptr_Send_req_data

inline MPI_Request * get_ptr_Send_req_data()

Get pointer to Send_req.

function get_ptr_L_jc_RcvLst_data

inline const iReg * get_ptr_L_jc_RcvLst_data() const

Get pointer to L_jc_RcvLst.

function get_ptr_L_jc_SndLst_data

inline const iReg * get_ptr_L_jc_SndLst_data() const

Get pointer to L_jc_SndLst.

function get_L_jc_SndLst_size

inline iReg get_L_jc_SndLst_size() const

Get size of L_jc_SndLst.

function get_ptr_R_jc_RcvLst_data

inline const iReg * get_ptr_R_jc_RcvLst_data() const

Get pointer to R_jc_RcvLst.

function get_ptr_R_jc_SndLst_data

inline const iReg * get_ptr_R_jc_SndLst_data() const

Get pointer to R_jc_SndLst.

function get_R_jc_SndLst_size

inline iReg get_R_jc_SndLst_size() const

Get size of R_jc_SndLst.

function get_ptr_in_list_com_scr_data

inline iReg * get_ptr_in_list_com_scr_data()

Get global incoming comunication list.

function get_ptr_out_list_com_scr_data

inline iReg * get_ptr_out_list_com_scr_data()

Get global outgoing comunication list.

function Get_GlobalNNZ

iGlo Get_GlobalNNZ() const

Evaluates the total number of nnz of the distributed csr matrices.

function Get_Globalnrows

virtual iGlo Get_Globalnrows() const

Reimplements: MatrixProd::Get_Globalnrows

Evaluates the total number of rows of the distributed matrix.

function Get_GlobalMaxRowNorm

rExt Get_GlobalMaxRowNorm() const

Evaluates the maximum value of the norm of a row of the matrix.

function Get_GlobalLocMaxRowNorm

struct LocMax Get_GlobalLocMaxRowNorm() const

Evaluates the maximum value and location of the norm of a row of the matrix.

function Get_GlobalAvgRowNorm

rExt Get_GlobalAvgRowNorm() const

Evaluates the average value of the norm of a row of the matrix.

function Convert_arrayCSR_To_ScalarLRGCSR

CSRMAT< iExt > Convert_arrayCSR_To_ScalarLRGCSR(
    iReg * new_order =nullptr
) const

Creates a CSRMAT starting from the array of local csr matrices. Useful in printing (used by fprintAll_PARDDMAT).

Return: CSRMAT

function cpt_matBalance

void cpt_matBalance(
    const type_MPI_iReg rank,
    const type_MPI_iReg size,
    struct bal_info & info
)

Function to compute the load balance of current matrix. The result is stored in process with rank 0.

function print_matBalance

void print_matBalance(
    const type_MPI_iReg rank,
    const struct bal_info info
)

Function to print the load balance of current matrix.

function RowStats

void RowStats(
    const iReg nb,
    const rExt * UpBounds,
    rExt * Stats
) const

Function to compute row statistics of the DSMat matrix. The norm of the rows will be divided in the user provided classes and written in the Stats array. The array UpBound defines the upper bound of each rownorm category, i.e. if the row norm is below to UpBound[i], it belongs to the class i.

function fPrintAll_PARDDMAT_ASCII

void fPrintAll_PARDDMAT_ASCII(
    const string & filename
) const

Prints the distributed dense matrix in ASCII filename (overwrite if it exists). It uses Convert_arrayCSR_To_ScalarLRGCSR and CreateGlobalOffset.

function fPrintAll_PARDDMAT_BIN

void fPrintAll_PARDDMAT_BIN(
    const string & filename
) const

Prints the distributed dense matrix in BINARY filename (overwrite if it exists). It uses Convert_arrayCSR_To_ScalarLRGCSR and CreateGlobalOffset.

function Prepare_MxV

virtual void Prepare_MxV(
    const iReg nRHS
)

Prepares the data structure necessary for MxV. L/R_SndRcvLst, Send/Recv_req,Recv_flag,WL_rcv,WR_snd,WR_rcv,WR_snd,WL_snd.

Reimplements: MatrixProd::Prepare_MxV

function UndoPrepare_MxV

virtual void UndoPrepare_MxV()

Undo the data structure necessary for MxV.

Reimplements: MatrixProd::UndoPrepare_MxV

function MxV

virtual void MxV(
    DDMat &__restrict__ x,
    DDMat &__restrict__ b,
    bool Barrier_flag
)

Computes the Matrix by vector Product.

Parameters:

  • x vector that is multiplied by the matrix.
  • b vector where result is stored.
  • Barrier_flag if true global Barrier is called before return.

Reimplements: MatrixProd::MxV

function MxM

void MxM(
    DSMat &__restrict__ B,
    DSMat &__restrict__ C
)

Computes the Matrix by Matrix Product A x B = C (with A equal to this).

Parameters:

  • B matrix that multiplies the matrix.
  • C matrix where result is stored.

function TreatBC

void TreatBC()

This function modifies a DSMat instance by detecting rows/columns associated to Dirichlet boundary conditions and setting to those rows/columns a diagonal entry equal to the maximum diagonal entry of the whole matrix. In such a way the minimum eigenvalue of the matrix is not erroneously associated to the arbitrary value given to the diagonal entries of boundary dofs. Detection is carried out by performing one Jacobi smoothing step on the constant vector and setting as BC those dofs whose value is zero

function Accumulate

void Accumulate(
    const rExt fa,
    const rExt fb,
    DSMat &__restrict__ B
)

Accumulate a matrix B over A (with ‘this’ == A) such that A = fa x A + fb x B.

Parameters:

  • fa multiplication factor of A.
  • fb multiplication factor of B.
  • B matrix B to sum.

NOTE: the pattern of A must contain the pattern of B.

function Transpose

void Transpose(
    DSMat &__restrict__ MatT
)

Transpose the matrix.

Parameters:

  • MatT Tranposed Matrix.

function SparseToDense

void SparseToDense(
    DDMat &__restrict__ DenseMat
)

Compacts the DSMat into a DDMat.

Parameters:

  • DenseMat dense matrix

function RemoveEmptyCSRs

void RemoveEmptyCSRs()

Remove empty CSRs.

Public Attributes Documentation

variable MatGraph

CommGraph * MatGraph = nullptr;

Incoming (for MxV and MxM) Communication graph.

variable OutMatGraph

CommGraph * OutMatGraph = nullptr;

Outgoing (for MxV and MxM) Communication graph If the matrix is symmetric OutMatGraph just points to the same MatGraph data.

variable ptDom

VEC< iReg > * ptDom = nullptr;

vector of pointers to the compacted column indices of each incoming block. That is, ptDom[ib+1] - ptDom[ib] gives the number of columns of the incoming block ib (ib block on this row of blocks).

variable OutptDom

VEC< iReg > * OutptDom = nullptr;

transposed local compact vector to the unknowns

vector of pointers to the compacted row indices of each outgoing block. That is, OutptDom[ib+1] - OutptDom[ib] gives the number of rows of the outgoing block ib (ib block on the column of blocks including this rank). If the matrix is symmetric OutptDom just points to the same ptDom data

variable csr

VEC< CSRMAT< iReg > > csr;

vector of CSR.

The block C (diagonal) must be allocated also if it is null because matC must point to something

variable matR

CSRMAT< iReg > * matR = nullptr;

array to the Right CSR: each of this will be a reference to the previous vector of CSR (no memory associated)

variable matL

CSRMAT< iReg > * matL = nullptr;

array to the left CSR: each of this will be a reference to a CSR mat (no memory associated)

variable matC

CSRMAT< iReg > * matC = nullptr;

pointer to the diagonal entry of the vector csr (no memory associated)


Updated on 12 February 2021 at 11:59:50 CET