Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/petsc_full_matrix.h>
Public Types | |
using | size_type = types::global_dof_index |
Public Types inherited from PETScWrappers::MatrixBase | |
using | const_iterator = MatrixIterators::const_iterator |
using | size_type = types::global_dof_index |
using | value_type = PetscScalar |
Public Member Functions | |
FullMatrix () | |
FullMatrix (const size_type m, const size_type n) | |
void | reinit (const size_type m, const size_type n) |
virtual const MPI_Comm & | get_mpi_communicator () const override |
Public Member Functions inherited from PETScWrappers::MatrixBase | |
MatrixBase () | |
MatrixBase (const MatrixBase &)=delete | |
MatrixBase & | operator= (const MatrixBase &)=delete |
virtual | ~MatrixBase () override |
MatrixBase & | operator= (const value_type d) |
void | clear () |
void | set (const size_type i, const size_type j, const PetscScalar value) |
void | set (const std::vector< size_type > &indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=false) |
void | set (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=false) |
void | set (const size_type row, const std::vector< size_type > &col_indices, const std::vector< PetscScalar > &values, const bool elide_zero_values=false) |
void | set (const size_type row, const size_type n_cols, const size_type *col_indices, const PetscScalar *values, const bool elide_zero_values=false) |
void | add (const size_type i, const size_type j, const PetscScalar value) |
void | add (const std::vector< size_type > &indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=true) |
void | add (const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< PetscScalar > &full_matrix, const bool elide_zero_values=true) |
void | add (const size_type row, const std::vector< size_type > &col_indices, const std::vector< PetscScalar > &values, const bool elide_zero_values=true) |
void | add (const size_type row, const size_type n_cols, const size_type *col_indices, const PetscScalar *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false) |
void | clear_row (const size_type row, const PetscScalar new_diag_value=0) |
void | clear_rows (const std::vector< size_type > &rows, const PetscScalar new_diag_value=0) |
void | compress (const VectorOperation::values operation) |
PetscScalar | operator() (const size_type i, const size_type j) const |
PetscScalar | el (const size_type i, const size_type j) const |
PetscScalar | diag_element (const size_type i) const |
size_type | m () const |
size_type | n () const |
size_type | local_size () const |
std::pair< size_type, size_type > | local_range () const |
bool | in_local_range (const size_type index) const |
size_type | n_nonzero_elements () const |
size_type | row_length (const size_type row) const |
PetscReal | l1_norm () const |
PetscReal | linfty_norm () const |
PetscReal | frobenius_norm () const |
PetscScalar | matrix_norm_square (const VectorBase &v) const |
PetscScalar | matrix_scalar_product (const VectorBase &u, const VectorBase &v) const |
PetscScalar | trace () const |
MatrixBase & | operator*= (const PetscScalar factor) |
MatrixBase & | operator/= (const PetscScalar factor) |
MatrixBase & | add (const PetscScalar factor, const MatrixBase &other) |
MatrixBase & | add (const MatrixBase &other, const PetscScalar factor) |
void | vmult (VectorBase &dst, const VectorBase &src) const |
void | Tvmult (VectorBase &dst, const VectorBase &src) const |
void | vmult_add (VectorBase &dst, const VectorBase &src) const |
void | Tvmult_add (VectorBase &dst, const VectorBase &src) const |
PetscScalar | residual (VectorBase &dst, const VectorBase &x, const VectorBase &b) const |
const_iterator | begin () const |
const_iterator | end () const |
const_iterator | begin (const size_type r) const |
const_iterator | end (const size_type r) const |
operator Mat () const | |
Mat & | petsc_matrix () |
void | transpose () |
PetscBool | is_symmetric (const double tolerance=1.e-12) |
PetscBool | is_hermitian (const double tolerance=1.e-12) |
void | write_ascii (const PetscViewerFormat format=PETSC_VIEWER_DEFAULT) |
void | print (std::ostream &out, const bool alternative_output=false) const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
void | do_reinit (const size_type m, const size_type n) |
Additional Inherited Members | |
Static Public Member Functions inherited from PETScWrappers::MatrixBase | |
static ::ExceptionBase & | ExcSourceEqualsDestination () |
static ::ExceptionBase & | ExcWrongMode (int arg1, int arg2) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions inherited from PETScWrappers::MatrixBase | |
void | prepare_action (const VectorOperation::values new_action) |
void | assert_is_compressed () |
void | prepare_add () |
void | prepare_set () |
void | mmult (MatrixBase &C, const MatrixBase &B, const VectorBase &V) const |
void | Tmmult (MatrixBase &C, const MatrixBase &B, const VectorBase &V) const |
Protected Attributes inherited from PETScWrappers::MatrixBase | |
Mat | matrix |
VectorOperation::values | last_action |
Implementation of a sequential dense matrix class based on PETSc. All the functionality is actually in the base class, except for the calls to generate a sequential dense matrix. This is possible since PETSc only works on an abstract matrix type and internally distributes to functions that do the actual work depending on the actual matrix type (much like using virtual functions). Only the functions creating a matrix of specific type differ, and are implemented in this particular class.
Definition at line 49 of file petsc_full_matrix.h.
Declare type for container size.
Definition at line 55 of file petsc_full_matrix.h.
FullMatrix< number >::FullMatrix | ( | ) |
Default constructor. Create an empty matrix.
Definition at line 27 of file petsc_full_matrix.cc.
FullMatrix< number >::FullMatrix | ( | const size_type | m, |
const size_type | n | ||
) |
Create a full matrix of dimensions m
times n
.
Definition at line 33 of file petsc_full_matrix.cc.
void FullMatrix< number >::reinit | ( | const size_type | m, |
const size_type | n | ||
) |
Throw away the present matrix and generate one that has the same properties as if it were created by the constructor of this class with the same argument list as the present function.
Definition at line 39 of file petsc_full_matrix.cc.
|
overridevirtual |
Return a reference to the MPI communicator object in use with this matrix. Since this is a sequential matrix, it returns the MPI_COMM_SELF communicator.
Implements PETScWrappers::MatrixBase.
Definition at line 61 of file petsc_full_matrix.cc.
|
private |
Do the actual work for the respective reinit() function and the matching constructor, i.e. create a matrix. Getting rid of the previous matrix is left to the caller.
Definition at line 50 of file petsc_full_matrix.cc.