OGS
|
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) | |
PETScMatrix & | operator= (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. | |
using MathLib::PETScMatrix::IndexType = PetscInt |
Definition at line 34 of file PETScMatrix.h.
|
inline |
Definition at line 37 of file PETScMatrix.h.
MathLib::PETScMatrix::PETScMatrix | ( | const PetscInt | nrows, |
const PETScMatrixOption & | mat_op = PETScMatrixOption() ) |
Constructor for a square matrix partitioning with more options.
nrows | The number of rows of the matrix or the local matrix. |
mat_op | The configuration information for creating a matrix. |
Definition at line 23 of file PETScMatrix.cpp.
References create(), MathLib::PETScMatrixOption::d_nz, MathLib::PETScMatrixOption::is_global_size, and MathLib::PETScMatrixOption::o_nz.
MathLib::PETScMatrix::PETScMatrix | ( | const PetscInt | nrows, |
const PetscInt | ncols, | ||
const PETScMatrixOption & | mat_op = PETScMatrixOption() ) |
Constructor for a rectangular matrix partitioning with more options.
nrows | The number of global or local rows. |
ncols | The number of global or local columns. |
mat_op | The configuration information for creating a matrix. |
Definition at line 40 of file PETScMatrix.cpp.
References create(), MathLib::PETScMatrixOption::d_nz, MathLib::PETScMatrixOption::is_global_size, and MathLib::PETScMatrixOption::o_nz.
|
inline |
MathLib::PETScMatrix::PETScMatrix | ( | PETScMatrix const & | A | ) |
Definition at line 58 of file PETScMatrix.cpp.
References A_.
|
inline |
Add value to a single entry.
i | The row index. |
j | The column index. |
value | The entry value. |
Definition at line 124 of file PETScMatrix.h.
References A_.
Referenced by add(), MathLib::addToMatrix(), and MathLib::setMatrix().
|
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.
References add(), MathLib::RowColumnIndices< IDX_TYPE >::columns, ncols_, and MathLib::RowColumnIndices< IDX_TYPE >::rows.
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.
row_pos | The global indices of the rows of the dense sub-matrix. |
col_pos | The global indices of the columns of the dense sub-matrix. |
sub_mat | A dense sub-matrix to be added. Its data of which must be row oriented stored. |
Definition at line 262 of file PETScMatrix.h.
References A_.
|
private |
Create the matrix, configure memory allocation and set the related member data.
d_nz | Number of nonzeros per row in the diagonal portion of local submatrix (same value is used for all local rows), |
o_nz | Number 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.
References A_, end_rank_, n_loc_cols_, n_loc_rows_, ncols_, nrows_, and start_rank_.
Referenced by PETScMatrix(), and PETScMatrix().
|
inlineprivate |
Definition at line 217 of file PETScMatrix.h.
References A_.
Referenced by ~PETScMatrix(), and operator=().
|
inline |
Perform MPI collection of assembled entries in buffer.
asm_type | Assembly type, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY |
Definition at line 67 of file PETScMatrix.h.
References A_.
Referenced by MathLib::applyKnownSolution(), MathLib::LinAlg::finalizeAssembly(), and viewer().
|
inline |
Get the number of columns.
Definition at line 76 of file PETScMatrix.h.
References ncols_.
Referenced by MathLib::addToMatrix(), and MathLib::setMatrix().
|
inline |
Get the number of local columns.
Definition at line 80 of file PETScMatrix.h.
References n_loc_cols_.
|
inline |
|
inline |
Get the number of rows.
Definition at line 74 of file PETScMatrix.h.
References nrows_.
Referenced by MathLib::addToMatrix(), and MathLib::setMatrix().
|
inline |
Get the start global index of the rows of the same rank.
Definition at line 82 of file PETScMatrix.h.
References start_rank_.
|
inline |
Get the end global index of the rows in the same rank.
Definition at line 84 of file PETScMatrix.h.
References end_rank_.
|
inline |
Get matrix reference.
Definition at line 86 of file PETScMatrix.h.
References A_.
Referenced by MathLib::LinAlg::axpy(), MathLib::LinAlg::aypx(), MathLib::LinAlg::matMult(), MathLib::LinAlg::matMultAdd(), MathLib::LinAlg::scale(), and MathLib::PETScLinearSolver::solve().
|
inline |
Get a matrix reference.
Definition at line 93 of file PETScMatrix.h.
References A_.
PETScMatrix & MathLib::PETScMatrix::operator= | ( | PETScMatrix const & | A | ) |
Definition at line 69 of file PETScMatrix.cpp.
References A_, destroy(), end_rank_, n_loc_cols_, n_loc_rows_, ncols_, nrows_, and start_rank_.
|
inline |
Set a single entry with a value.
i | The row index. |
j | The column index. |
value | The entry value. |
Definition at line 113 of file PETScMatrix.h.
References A_.
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.
row_pos | The row indices of the specified rows. |
Definition at line 92 of file PETScMatrix.cpp.
References A_.
Referenced by MathLib::applyKnownSolution().
|
inline |
Set all entries to zero.
Definition at line 95 of file PETScMatrix.h.
References A_.
Referenced by MathLib::setMatrix(), and MathLib::setMatrix().
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.
file_name | File name for output |
vw_format | File 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.
References A_, finalizeAssembly(), and viewer().
Referenced by viewer().
|
friend |
General interface for the matrix assembly.
mat | The matrix to be finalized. |
asm_type | Assembly type, either MAT_FLUSH_ASSEMBLY or MAT_FINAL_ASSEMBLY. |
Definition at line 158 of file PETScMatrix.cpp.
|
private |
PETSc matrix.
Definition at line 227 of file PETScMatrix.h.
Referenced by PETScMatrix(), add(), add(), create(), destroy(), finalizeAssembly(), getRawMatrix(), getRawMatrix(), operator=(), set(), setRowsColumnsZero(), setZero(), and viewer().
|
private |
Ending index in a rank.
Definition at line 245 of file PETScMatrix.h.
Referenced by create(), getRangeEnd(), and operator=().
|
private |
Number of the local columns.
Definition at line 239 of file PETScMatrix.h.
Referenced by create(), getNumberOfLocalColumns(), and operator=().
|
private |
Number of the local rows.
Definition at line 236 of file PETScMatrix.h.
Referenced by create(), getNumberOfLocalRows(), and operator=().
|
private |
Number of the global columns.
Definition at line 233 of file PETScMatrix.h.
Referenced by add(), create(), getNumberOfColumns(), and operator=().
|
private |
Number of the global rows.
Definition at line 230 of file PETScMatrix.h.
Referenced by create(), getNumberOfRows(), and operator=().
|
private |
Starting index in a rank.
Definition at line 242 of file PETScMatrix.h.
Referenced by create(), getRangeBegin(), and operator=().