Class ae108::cpppetsc::Matrix
template <class Policy>
ClassList > ae108 > cpppetsc > Matrix
Classes
Type | Name |
---|---|
class | AssemblyView |
Public Types
Type | Name |
---|---|
typedef TaggedEntity< size_type, struct GlobalColsTag > | GlobalCols |
typedef TaggedEntity< size_type, struct GlobalRowsTag > | GlobalRows |
typedef TaggedEntity< size_type, struct LocalColsTag > | LocalCols |
typedef TaggedEntity< size_type, struct LocalRowsTag > | LocalRows |
typedef PetscReal | real_type |
typedef PetscInt | size_type |
typedef PetscScalar | value_type |
Public Functions
Type | Name |
---|---|
Matrix (size_type global_rows, size_type global_columns) Allocates a matrix with global_rows rows and global_columns columns. |
|
Matrix (LocalRows local_rows, LocalCols local_columns, GlobalRows global_rows, GlobalCols global_columns) Allocates a matrix with the provided layout. |
|
Matrix (UniqueEntity< Mat > mat) Wraps the given UniqueEntity in a Matrix . |
|
void | addAlphaX (const value_type alpha, const Matrix & X) Adds alpha * X to the matrix. |
AssemblyView | assemblyView () Returns a view that permits assembly of the matrix (e.g. adding entries to the current entries). If the number of nonzero entries per row is known in advance, then it is more efficient to use the method preallocatedAssemblyView(...). |
size_type | blockSize () const Returns the block size of the matrix. |
Mat | data () const Returns the internal matrix. |
void | finalize () Perform final assembly of the matrix. |
std::pair< size_type, size_type > | localColumnRange () const Returns the local column indices in the format: (first column index, last column index + 1). |
std::pair< size_type, size_type > | localRowRange () const Returns the local row indices in the format: (first row index, last row index + 1). |
std::pair< size_type, size_type > | localSize () const Returns a pair of (number of rows, number of columns). |
real_type | norm () const Computes the Frobenius norm of the matrix. |
value_type | operator() (size_type row, size_type column) const Returns the entry at the provided coordinates. |
AssemblyView | preallocatedAssemblyView (const size_type nonZeroesPerRow) Returns a view that permits assembly of the matrix (e.g. adding entries to the current entries). |
void | print () const Print the matrix to world stdout. |
void | replaceRowsByEye (const std::vector< size_type > & rows) Replaces row indices by the corresponding rows of the identity matrix. Does not change the nonzero pattern. |
void | scale (const value_type factor) Scale the matrix by factor. |
void | setZero () Replaces all entries by zero, but keeps the structure for sparse matrices. |
std::pair< size_type, size_type > | size () const Returns a pair of (number of rows, number of columns). |
Public Static Functions
Type | Name |
---|---|
Matrix | clone (const Matrix & matrix) Returns a deep copy of the matrix. |
Matrix | fromAtB (const Matrix & A, const Matrix & B) creates a matrix A^t * B. |
Matrix | fromInvertedBlockDiagonalOf (const Matrix & matrix) Creates a matrix from the inverted block diagonal. |
Matrix | fromLayoutOf (const Matrix & matrix) Creates a matrix with the same layout and nonzero pattern as the provided matrix. |
Matrix | fromList (const std::initializer_list< std::initializer_list< value_type > > list) Creates a Matrix from a list of rows. |
Matrix | fromMesh (const Mesh< Policy > & mesh) Creates a matrix with a structure defined by the given mesh. |
Matrix | fromPAPt (const Matrix & P, const Matrix & A) Creates a matrix P * A * P^t. |
Matrix | fromProduct (const Matrix & A, const Matrix & B, const Matrix & C) Creates a matrix A * B * C. |
Matrix | fromProduct (const Matrix & A, const Matrix & B) Creates a matrix A * B. |
Matrix | fromPtAP (const Matrix & P, const Matrix & A) Creates a matrix P^t * A * P. |
Public Types Documentation
typedef GlobalCols
using ae108::cpppetsc::Matrix< Policy >::GlobalCols = TaggedEntity<size_type, struct GlobalColsTag>;
typedef GlobalRows
using ae108::cpppetsc::Matrix< Policy >::GlobalRows = TaggedEntity<size_type, struct GlobalRowsTag>;
typedef LocalCols
using ae108::cpppetsc::Matrix< Policy >::LocalCols = TaggedEntity<size_type, struct LocalColsTag>;
typedef LocalRows
using ae108::cpppetsc::Matrix< Policy >::LocalRows = TaggedEntity<size_type, struct LocalRowsTag>;
typedef real_type
using ae108::cpppetsc::Matrix< Policy >::real_type = PetscReal;
typedef size_type
using ae108::cpppetsc::Matrix< Policy >::size_type = PetscInt;
typedef value_type
using ae108::cpppetsc::Matrix< Policy >::value_type = PetscScalar;
Public Functions Documentation
function Matrix [1/3]
Allocates a matrix with global_rows rows and global_columns columns.
explicit ae108::cpppetsc::Matrix::Matrix (
size_type global_rows,
size_type global_columns
)
function Matrix [2/3]
Allocates a matrix with the provided layout.
explicit ae108::cpppetsc::Matrix::Matrix (
LocalRows local_rows,
LocalCols local_columns,
GlobalRows global_rows,
GlobalCols global_columns
)
Remark:
Use PETSC_DETERMINE for the global rows/columns to let PETSc determine the sum.
function Matrix [3/3]
Wraps the given UniqueEntity in a Matrix .
explicit ae108::cpppetsc::Matrix::Matrix (
UniqueEntity< Mat > mat
)
function addAlphaX
Adds alpha * X to the matrix.
void ae108::cpppetsc::Matrix::addAlphaX (
const value_type alpha,
const Matrix & X
)
Precondition:
The matrix X has the same nonzero pattern.
function assemblyView
Returns a view that permits assembly of the matrix (e.g. adding entries to the current entries). If the number of nonzero entries per row is known in advance, then it is more efficient to use the method preallocatedAssemblyView(...).
AssemblyView ae108::cpppetsc::Matrix::assemblyView ()
Remark:
Final assembly is automatically performed when the view and its inserters have been destroyed.
function blockSize
Returns the block size of the matrix.
size_type ae108::cpppetsc::Matrix::blockSize () const
function data
Returns the internal matrix.
Mat ae108::cpppetsc::Matrix::data () const
Remark:
For internal use only.
function finalize
Perform final assembly of the matrix.
void ae108::cpppetsc::Matrix::finalize ()
function localColumnRange
Returns the local column indices in the format: (first column index, last column index + 1).
std::pair< size_type, size_type > ae108::cpppetsc::Matrix::localColumnRange () const
function localRowRange
Returns the local row indices in the format: (first row index, last row index + 1).
std::pair< size_type, size_type > ae108::cpppetsc::Matrix::localRowRange () const
function localSize
Returns a pair of (number of rows, number of columns).
std::pair< size_type, size_type > ae108::cpppetsc::Matrix::localSize () const
function norm
Computes the Frobenius norm of the matrix.
real_type ae108::cpppetsc::Matrix::norm () const
function operator()
Returns the entry at the provided coordinates.
value_type ae108::cpppetsc::Matrix::operator() (
size_type row,
size_type column
) const
Remark:
Does not check bounds.
function preallocatedAssemblyView
Returns a view that permits assembly of the matrix (e.g. adding entries to the current entries).
AssemblyView ae108::cpppetsc::Matrix::preallocatedAssemblyView (
const size_type nonZeroesPerRow
)
Parameters:
nonZeroesPerRow
The maximum number of nonzero values per row to expect. If the matrix is of type MPIAIJ, then this parameter defines the number of nonzero values in the diagonal block. The same number of values is permitted in the off-diagonal block.
Exception:
- InvalidParametersException if matrix is not of type SEQAIJ or MPIAIJ.
Remark:
Final assembly is automatically performed when the view and its inserters have been destroyed.
function print
Print the matrix to world stdout.
void ae108::cpppetsc::Matrix::print () const
function replaceRowsByEye
Replaces row indices by the corresponding rows of the identity matrix. Does not change the nonzero pattern.
void ae108::cpppetsc::Matrix::replaceRowsByEye (
const std::vector< size_type > & rows
)
Parameters:
rows
Contains the global row indices to replace.
Remark:
Must be called on every process in the parallel setting.
function scale
Scale the matrix by factor.
void ae108::cpppetsc::Matrix::scale (
const value_type factor
)
function setZero
Replaces all entries by zero, but keeps the structure for sparse matrices.
void ae108::cpppetsc::Matrix::setZero ()
function size
Returns a pair of (number of rows, number of columns).
std::pair< size_type, size_type > ae108::cpppetsc::Matrix::size () const
Public Static Functions Documentation
function clone
Returns a deep copy of the matrix.
static Matrix ae108::cpppetsc::Matrix::clone (
const Matrix & matrix
)
function fromAtB
creates a matrix A^t * B.
static Matrix ae108::cpppetsc::Matrix::fromAtB (
const Matrix & A,
const Matrix & B
)
function fromInvertedBlockDiagonalOf
Creates a matrix from the inverted block diagonal.
static Matrix ae108::cpppetsc::Matrix::fromInvertedBlockDiagonalOf (
const Matrix & matrix
)
function fromLayoutOf
Creates a matrix with the same layout and nonzero pattern as the provided matrix.
static Matrix ae108::cpppetsc::Matrix::fromLayoutOf (
const Matrix & matrix
)
function fromList
Creates a Matrix from a list of rows.
static Matrix ae108::cpppetsc::Matrix::fromList (
const std::initializer_list< std::initializer_list< value_type > > list
)
Exception:
- InvalidParametersException if the number of columns is inconsistent.
function fromMesh
Creates a matrix with a structure defined by the given mesh.
static Matrix ae108::cpppetsc::Matrix::fromMesh (
const Mesh < Policy > & mesh
)
function fromPAPt
Creates a matrix P * A * P^t.
static Matrix ae108::cpppetsc::Matrix::fromPAPt (
const Matrix & P,
const Matrix & A
)
function fromProduct [1/2]
Creates a matrix A * B * C.
static Matrix ae108::cpppetsc::Matrix::fromProduct (
const Matrix & A,
const Matrix & B,
const Matrix & C
)
function fromProduct [2/2]
Creates a matrix A * B.
static Matrix ae108::cpppetsc::Matrix::fromProduct (
const Matrix & A,
const Matrix & B
)
function fromPtAP
Creates a matrix P^t * A * P.
static Matrix ae108::cpppetsc::Matrix::fromPtAP (
const Matrix & P,
const Matrix & A
)
The documentation for this class was generated from the following file cpppetsc/src/include/ae108/cpppetsc/Matrix.h