OGS
MathLib::PETScMatrix Class Reference

Detailed Description

Wrapper class for PETSc matrix routines for matrix.

Definition at line 31 of file PETScMatrix.h.

#include <PETScMatrix.h>

Public Types

using IndexType = PetscInt
 

Public Member Functions

 PETScMatrix ()
 
 PETScMatrix (const PetscInt nrows, const PETScMatrixOption &mat_op=PETScMatrixOption())
 Constructor for a square matrix partitioning with more options.
 
 PETScMatrix (const PetscInt nrows, const PetscInt ncols, const PETScMatrixOption &mat_op=PETScMatrixOption())
 Constructor for a rectangular matrix partitioning with more options.
 
 ~PETScMatrix ()
 
 PETScMatrix (PETScMatrix const &A)
 
PETScMatrixoperator= (PETScMatrix const &A)
 
void finalizeAssembly (const MatAssemblyType asm_type=MAT_FINAL_ASSEMBLY)
 Perform MPI collection of assembled entries in buffer.
 
PetscInt getNumberOfRows () const
 Get the number of rows.
 
PetscInt getNumberOfColumns () const
 Get the number of columns.
 
PetscInt getNumberOfLocalRows () const
 Get the number of local rows.
 
PetscInt getNumberOfLocalColumns () const
 Get the number of local columns.
 
PetscInt getRangeBegin () const
 Get the start global index of the rows of the same rank.
 
PetscInt getRangeEnd () const
 Get the end global index of the rows in the same rank.
 
Mat & getRawMatrix ()
 Get matrix reference.
 
Mat const & getRawMatrix () const
 
void setZero ()
 Set all entries to zero.
 
void setRowsColumnsZero (std::vector< PetscInt > const &row_pos)
 Set the specified rows to zero except diagonal entries, i.e. \(A(k, j) = \begin{cases} 0.0, &j\not=k, j=1,2,\dots,k-1, k+1, \dots, n \\ 1.0, &j = k \end{cases}\), where \(k \in \mbox{row\_pos}\) This function must be called by all ranks.
 
void set (const PetscInt i, const PetscInt j, const PetscScalar value)
 Set a single entry with a value.
 
void add (const PetscInt i, const PetscInt j, const PetscScalar value)
 Add value to a single entry.
 
template<class T_DENSE_MATRIX >
void add (RowColumnIndices< PetscInt > const &indices, const T_DENSE_MATRIX &sub_matrix)
 Add sub-matrix at positions given by global indices, in which negative index indicates ghost entry.
 
template<class T_DENSE_MATRIX >
void add (std::vector< PetscInt > const &row_pos, std::vector< PetscInt > const &col_pos, const T_DENSE_MATRIX &sub_mat)
 Add a dense sub-matrix to a PETSc matrix.
 
void viewer (const std::string &file_name, const PetscViewerFormat vw_format=PETSC_VIEWER_ASCII_MATLAB)
 

Private Member Functions

void destroy ()
 
void create (const PetscInt d_nz, const PetscInt o_nz)
 Create the matrix, configure memory allocation and set the related member data.
 

Private Attributes

Mat A_ = nullptr
 PETSc matrix.
 
PetscInt nrows_
 Number of the global rows.
 
PetscInt ncols_
 Number of the global columns.
 
PetscInt n_loc_rows_
 Number of the local rows.
 
PetscInt n_loc_cols_
 Number of the local columns.
 
PetscInt start_rank_
 Starting index in a rank.
 
PetscInt end_rank_
 Ending index in a rank.
 

Friends

bool finalizeMatrixAssembly (PETScMatrix &mat, const MatAssemblyType asm_type=MAT_FINAL_ASSEMBLY)
 General interface for the matrix assembly.
 

Member Typedef Documentation

◆ IndexType

Definition at line 34 of file PETScMatrix.h.

Constructor & Destructor Documentation

◆ PETScMatrix() [1/4]

MathLib::PETScMatrix::PETScMatrix ( )
inline

Definition at line 37 of file PETScMatrix.h.

37{}

◆ PETScMatrix() [2/4]

MathLib::PETScMatrix::PETScMatrix ( const PetscInt nrows,
const PETScMatrixOption & mat_op = PETScMatrixOption() )

Constructor for a square matrix partitioning with more options.

Parameters
nrowsThe number of rows of the matrix or the local matrix.
mat_opThe configuration information for creating a matrix.

Definition at line 23 of file PETScMatrix.cpp.

24 : nrows_(nrows),
25 ncols_(nrows),
26 n_loc_rows_(PETSC_DECIDE),
27 n_loc_cols_(mat_opt.n_local_cols)
28{
29 if (!mat_opt.is_global_size)
30 {
31 n_loc_rows_ = nrows;
32 n_loc_cols_ = nrows;
33 nrows_ = PETSC_DECIDE;
34 ncols_ = PETSC_DECIDE;
35 }
36
37 create(mat_opt.d_nz, mat_opt.o_nz);
38}
PetscInt n_loc_cols_
Number of the local columns.
void create(const PetscInt d_nz, const PetscInt o_nz)
Create the matrix, configure memory allocation and set the related member data.
PetscInt n_loc_rows_
Number of the local rows.
PetscInt nrows_
Number of the global rows.
PetscInt ncols_
Number of the global columns.

References create(), MathLib::PETScMatrixOption::d_nz, MathLib::PETScMatrixOption::is_global_size, and MathLib::PETScMatrixOption::o_nz.

◆ PETScMatrix() [3/4]

MathLib::PETScMatrix::PETScMatrix ( const PetscInt nrows,
const PetscInt ncols,
const PETScMatrixOption & mat_op = PETScMatrixOption() )

Constructor for a rectangular matrix partitioning with more options.

Parameters
nrowsThe number of global or local rows.
ncolsThe number of global or local columns.
mat_opThe configuration information for creating a matrix.

Definition at line 40 of file PETScMatrix.cpp.

42 : nrows_(nrows),
43 ncols_(ncols),
44 n_loc_rows_(PETSC_DECIDE),
45 n_loc_cols_(mat_opt.n_local_cols)
46{
47 if (!mat_opt.is_global_size)
48 {
49 nrows_ = PETSC_DECIDE;
50 ncols_ = PETSC_DECIDE;
51 n_loc_rows_ = nrows;
52 n_loc_cols_ = ncols;
53 }
54
55 create(mat_opt.d_nz, mat_opt.o_nz);
56}

References create(), MathLib::PETScMatrixOption::d_nz, MathLib::PETScMatrixOption::is_global_size, and MathLib::PETScMatrixOption::o_nz.

◆ ~PETScMatrix()

MathLib::PETScMatrix::~PETScMatrix ( )
inline

Definition at line 57 of file PETScMatrix.h.

57{ destroy(); }

References destroy().

◆ PETScMatrix() [4/4]

MathLib::PETScMatrix::PETScMatrix ( PETScMatrix const & A)

Definition at line 58 of file PETScMatrix.cpp.

59 : nrows_(A.nrows_),
60 ncols_(A.ncols_),
61 n_loc_rows_(A.n_loc_rows_),
62 n_loc_cols_(A.n_loc_cols_),
63 start_rank_(A.start_rank_),
64 end_rank_(A.end_rank_)
65{
66 MatConvert(A.A_, MATSAME, MAT_INITIAL_MATRIX, &A_);
67}
PetscInt end_rank_
Ending index in a rank.
Mat A_
PETSc matrix.
PetscInt start_rank_
Starting index in a rank.

References A_.

Member Function Documentation

◆ add() [1/3]

void MathLib::PETScMatrix::add ( const PetscInt i,
const PetscInt j,
const PetscScalar value )
inline

Add value to a single entry.

Parameters
iThe row index.
jThe column index.
valueThe entry value.

Definition at line 124 of file PETScMatrix.h.

125 {
126 MatSetValue(A_, i, j, value, ADD_VALUES);
127 }

References A_.

Referenced by add(), MathLib::addToMatrix(), and MathLib::setMatrix().

◆ add() [2/3]

template<class T_DENSE_MATRIX >
void MathLib::PETScMatrix::add ( RowColumnIndices< PetscInt > const & indices,
const T_DENSE_MATRIX & sub_matrix )
inline

Add sub-matrix at positions given by global indices, in which negative index indicates ghost entry.

In order to use MatZeroRows to apply Dirichlet boundary condition, entries in the rows with the negative global indices are skipped to added to the global matrix, meanwhile entries in the columns with the negative global indices are added the global matrix. By using MatZeroRows to apply Dirichlet boundary condition, the off diagonal entries in ghost rows of the global matrix are set to zero, while the off diagonal entries in ghost rows of the global matrix are assembled and kept for linear solver.

For the setting of Dirichlet boundary condition in PETSc, please refer to the PETSc:Documentation:FAQ.

Definition at line 146 of file PETScMatrix.h.

148 {
149 // Set global column indices to positive to allow all entries of columns
150 // to be added to the global matrix. For the ghost columns, only the
151 // off diagonal entries are added due to the negative indices of the
152 // corresponding rows.
153 std::vector<PetscInt> cols;
154 cols.reserve(indices.columns.size());
155 for (auto col : indices.columns)
156 {
157 // Ghost entries, and its original index is 0.
158 if (col == -ncols_)
159 {
160 cols.push_back(0);
161 }
162 else
163 {
164 cols.push_back(std::abs(col));
165 }
166 }
167
168 add(indices.rows, cols, sub_matrix);
169 }
void add(const PetscInt i, const PetscInt j, const PetscScalar value)
Add value to a single entry.

References add(), MathLib::RowColumnIndices< IDX_TYPE >::columns, ncols_, and MathLib::RowColumnIndices< IDX_TYPE >::rows.

◆ add() [3/3]

template<class T_DENSE_MATRIX >
void MathLib::PETScMatrix::add ( std::vector< PetscInt > const & row_pos,
std::vector< PetscInt > const & col_pos,
const T_DENSE_MATRIX & sub_mat )

Add a dense sub-matrix to a PETSc matrix.

Parameters
row_posThe global indices of the rows of the dense sub-matrix.
col_posThe global indices of the columns of the dense sub-matrix.
sub_matA dense sub-matrix to be added. Its data of which must be row oriented stored.

Definition at line 262 of file PETScMatrix.h.

265{
266 const PetscInt nrows = static_cast<PetscInt>(row_pos.size());
267 const PetscInt ncols = static_cast<PetscInt>(col_pos.size());
268
269 MatSetValues(A_, nrows, &row_pos[0], ncols, &col_pos[0], &sub_mat(0, 0),
270 ADD_VALUES);
271};

References A_.

◆ create()

void MathLib::PETScMatrix::create ( const PetscInt d_nz,
const PetscInt o_nz )
private

Create the matrix, configure memory allocation and set the related member data.

Parameters
d_nzNumber of nonzeros per row in the diagonal portion of local submatrix (same value is used for all local rows),
o_nzNumber of nonzeros per row in the off-diagonal portion of local submatrix (same value is used for all local rows)

Definition at line 140 of file PETScMatrix.cpp.

141{
142 MatCreate(PETSC_COMM_WORLD, &A_);
143 MatSetSizes(A_, n_loc_rows_, n_loc_cols_, nrows_, ncols_);
144
145 MatSetType(A_, MATAIJ);
146 MatSetFromOptions(A_);
147
148 MatSeqAIJSetPreallocation(A_, d_nz, PETSC_NULLPTR);
149 MatMPIAIJSetPreallocation(A_, d_nz, PETSC_NULLPTR, o_nz, PETSC_NULLPTR);
150 // If pre-allocation does not work one can use MatSetUp(A_), which is much
151 // slower.
152
153 MatGetOwnershipRange(A_, &start_rank_, &end_rank_);
154 MatGetSize(A_, &nrows_, &ncols_);
155 MatGetLocalSize(A_, &n_loc_rows_, &n_loc_cols_);
156}

References A_, end_rank_, n_loc_cols_, n_loc_rows_, ncols_, nrows_, and start_rank_.

Referenced by PETScMatrix(), and PETScMatrix().

◆ destroy()

void MathLib::PETScMatrix::destroy ( )
inlineprivate

Definition at line 217 of file PETScMatrix.h.

218 {
219 if (A_ != nullptr)
220 {
221 MatDestroy(&A_);
222 }
223 A_ = nullptr;
224 }

References A_.

Referenced by ~PETScMatrix(), and operator=().

◆ finalizeAssembly()

void MathLib::PETScMatrix::finalizeAssembly ( const MatAssemblyType asm_type = MAT_FINAL_ASSEMBLY)
inline

Perform MPI collection of assembled entries in buffer.

Parameters
asm_typeAssembly type, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY

Definition at line 67 of file PETScMatrix.h.

68 {
69 MatAssemblyBegin(A_, asm_type);
70 MatAssemblyEnd(A_, asm_type);
71 }

References A_.

Referenced by MathLib::applyKnownSolution(), MathLib::LinAlg::finalizeAssembly(), and viewer().

◆ getNumberOfColumns()

PetscInt MathLib::PETScMatrix::getNumberOfColumns ( ) const
inline

Get the number of columns.

Definition at line 76 of file PETScMatrix.h.

76{ return ncols_; }

References ncols_.

Referenced by MathLib::addToMatrix(), and MathLib::setMatrix().

◆ getNumberOfLocalColumns()

PetscInt MathLib::PETScMatrix::getNumberOfLocalColumns ( ) const
inline

Get the number of local columns.

Definition at line 80 of file PETScMatrix.h.

80{ return n_loc_cols_; }

References n_loc_cols_.

◆ getNumberOfLocalRows()

PetscInt MathLib::PETScMatrix::getNumberOfLocalRows ( ) const
inline

Get the number of local rows.

Definition at line 78 of file PETScMatrix.h.

78{ return n_loc_rows_; }

References n_loc_rows_.

◆ getNumberOfRows()

PetscInt MathLib::PETScMatrix::getNumberOfRows ( ) const
inline

Get the number of rows.

Definition at line 74 of file PETScMatrix.h.

74{ return nrows_; }

References nrows_.

Referenced by MathLib::addToMatrix(), and MathLib::setMatrix().

◆ getRangeBegin()

PetscInt MathLib::PETScMatrix::getRangeBegin ( ) const
inline

Get the start global index of the rows of the same rank.

Definition at line 82 of file PETScMatrix.h.

82{ return start_rank_; }

References start_rank_.

◆ getRangeEnd()

PetscInt MathLib::PETScMatrix::getRangeEnd ( ) const
inline

Get the end global index of the rows in the same rank.

Definition at line 84 of file PETScMatrix.h.

84{ return end_rank_; }

References end_rank_.

◆ getRawMatrix() [1/2]

Mat & MathLib::PETScMatrix::getRawMatrix ( )
inline

◆ getRawMatrix() [2/2]

Mat const & MathLib::PETScMatrix::getRawMatrix ( ) const
inline

Get a matrix reference.

Warning
This method is dangerous insofar as you can do arbitrary things also with a const PETSc matrix.

Definition at line 93 of file PETScMatrix.h.

93{ return A_; }

References A_.

◆ operator=()

PETScMatrix & MathLib::PETScMatrix::operator= ( PETScMatrix const & A)

Definition at line 69 of file PETScMatrix.cpp.

70{
71 nrows_ = A.nrows_;
72 ncols_ = A.ncols_;
73 n_loc_rows_ = A.n_loc_rows_;
74 n_loc_cols_ = A.n_loc_cols_;
75 start_rank_ = A.start_rank_;
76 end_rank_ = A.end_rank_;
77
78 if (A_ != nullptr)
79 {
80 // TODO this is the slowest option for copying
81 MatCopy(A.A_, A_, DIFFERENT_NONZERO_PATTERN);
82 }
83 else
84 {
85 destroy();
86 MatConvert(A.A_, MATSAME, MAT_INITIAL_MATRIX, &A_);
87 }
88
89 return *this;
90}

References A_, destroy(), end_rank_, n_loc_cols_, n_loc_rows_, ncols_, nrows_, and start_rank_.

◆ set()

void MathLib::PETScMatrix::set ( const PetscInt i,
const PetscInt j,
const PetscScalar value )
inline

Set a single entry with a value.

Parameters
iThe row index.
jThe column index.
valueThe entry value.

Definition at line 113 of file PETScMatrix.h.

114 {
115 MatSetValue(A_, i, j, value, INSERT_VALUES);
116 }

References A_.

◆ setRowsColumnsZero()

void MathLib::PETScMatrix::setRowsColumnsZero ( std::vector< PetscInt > const & row_pos)

Set the specified rows to zero except diagonal entries, i.e. \(A(k, j) = \begin{cases} 0.0, &j\not=k, j=1,2,\dots,k-1, k+1, \dots, n \\ 1.0, &j = k \end{cases}\), where \(k \in \mbox{row\_pos}\) This function must be called by all ranks.

Parameters
row_posThe row indices of the specified rows.

Definition at line 92 of file PETScMatrix.cpp.

93{
94 // Each rank (compute core) processes only the rows that belong to the rank
95 // itself.
96 const PetscScalar one = 1.0;
97 const PetscInt nrows = static_cast<PetscInt>(row_pos.size());
98
99 // Each process will only zero its own rows.
100 // This avoids all reductions in the zero row routines
101 // and thus improves performance for very large process counts.
102 // See PETSc doc about MAT_NO_OFF_PROC_ZERO_ROWS.
103 MatSetOption(A_, MAT_NO_OFF_PROC_ZERO_ROWS, PETSC_TRUE);
104
105 // Keep the non-zero pattern for the assignment operator.
106 MatSetOption(A_, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE);
107
108 if (nrows > 0)
109 {
110 MatZeroRows(A_, nrows, &row_pos[0], one, PETSC_NULLPTR, PETSC_NULLPTR);
111 }
112 else
113 {
114 MatZeroRows(A_, 0, PETSC_NULLPTR, one, PETSC_NULLPTR, PETSC_NULLPTR);
115 }
116}

References A_.

Referenced by MathLib::applyKnownSolution().

◆ setZero()

void MathLib::PETScMatrix::setZero ( )
inline

Set all entries to zero.

Definition at line 95 of file PETScMatrix.h.

95{ MatZeroEntries(A_); }

References A_.

Referenced by MathLib::setMatrix(), and MathLib::setMatrix().

◆ viewer()

void MathLib::PETScMatrix::viewer ( const std::string & file_name,
const PetscViewerFormat vw_format = PETSC_VIEWER_ASCII_MATLAB )

View the global vector for test purpose. Do not use it for output a big vector.

Parameters
file_nameFile name for output
vw_formatFile format listed as: PETSC_VIEWER_DEFAULT Default format PETSC_VIEWER_ASCII_MATLAB MATLAB format PETSC_VIEWER_ASCII_DENSE Print matrix as dense PETSC_VIEWER_ASCII_IMPL Implementation-specific format (which is in many cases the same as the default) PETSC_VIEWER_ASCII_INFO Basic information about object PETSC_VIEWER_ASCII_INFO_DETAIL More detailed info about object PETSC_VIEWER_ASCII_COMMON Identical output format for all objects of a particular type PETSC_VIEWER_ASCII_INDEX (for vectors) Prints the vector element number next to each vector entry PETSC_VIEWER_ASCII_SYMMODU Print parallel vectors without indicating the processor ranges PETSC_VIEWER_ASCII_VTK Outputs the object to a VTK file PETSC_VIEWER_NATIVE Store the object to the binary file in its native format (for example, dense matrices are stored as dense), DMDA vectors are dumped directly to the file instead of being first put in the natural ordering PETSC_VIEWER_DRAW_BASIC Views the vector with a simple 1d plot PETSC_VIEWER_DRAW_LG Views the vector with a line graph PETSC_VIEWER_DRAW_CONTOUR Views the vector with a contour plot

Definition at line 118 of file PETScMatrix.cpp.

120{
121 PetscViewer viewer;
122 PetscViewerASCIIOpen(PETSC_COMM_WORLD, file_name.c_str(), &viewer);
123 PetscViewerPushFormat(viewer, vw_format);
124
126
127 PetscObjectSetName((PetscObject)A_, "Stiffness_matrix");
128 MatView(A_, viewer);
129
130// This preprocessor is only for debugging, e.g. dump the matrix and exit the
131// program.
132// #define EXIT_TEST
133#ifdef EXIT_TEST
134 MatDestroy(A_);
135 PetscFinalize();
136 exit(0);
137#endif
138}
void finalizeAssembly(const MatAssemblyType asm_type=MAT_FINAL_ASSEMBLY)
Perform MPI collection of assembled entries in buffer.
Definition PETScMatrix.h:67
void viewer(const std::string &file_name, const PetscViewerFormat vw_format=PETSC_VIEWER_ASCII_MATLAB)

References A_, finalizeAssembly(), and viewer().

Referenced by viewer().

Friends And Related Symbol Documentation

◆ finalizeMatrixAssembly

bool finalizeMatrixAssembly ( PETScMatrix & mat,
const MatAssemblyType asm_type = MAT_FINAL_ASSEMBLY )
friend

General interface for the matrix assembly.

Parameters
matThe matrix to be finalized.
asm_typeAssembly type, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY.

Definition at line 158 of file PETScMatrix.cpp.

159{
160 mat.finalizeAssembly(asm_type);
161 return true;
162}

Member Data Documentation

◆ A_

Mat MathLib::PETScMatrix::A_ = nullptr
private

◆ end_rank_

PetscInt MathLib::PETScMatrix::end_rank_
private

Ending index in a rank.

Definition at line 245 of file PETScMatrix.h.

Referenced by create(), getRangeEnd(), and operator=().

◆ n_loc_cols_

PetscInt MathLib::PETScMatrix::n_loc_cols_
private

Number of the local columns.

Definition at line 239 of file PETScMatrix.h.

Referenced by create(), getNumberOfLocalColumns(), and operator=().

◆ n_loc_rows_

PetscInt MathLib::PETScMatrix::n_loc_rows_
private

Number of the local rows.

Definition at line 236 of file PETScMatrix.h.

Referenced by create(), getNumberOfLocalRows(), and operator=().

◆ ncols_

PetscInt MathLib::PETScMatrix::ncols_
private

Number of the global columns.

Definition at line 233 of file PETScMatrix.h.

Referenced by add(), create(), getNumberOfColumns(), and operator=().

◆ nrows_

PetscInt MathLib::PETScMatrix::nrows_
private

Number of the global rows.

Definition at line 230 of file PETScMatrix.h.

Referenced by create(), getNumberOfRows(), and operator=().

◆ start_rank_

PetscInt MathLib::PETScMatrix::start_rank_
private

Starting index in a rank.

Definition at line 242 of file PETScMatrix.h.

Referenced by create(), getRangeBegin(), and operator=().


The documentation for this class was generated from the following files: