Go to the documentation of this file.
20 #ifndef LIBMESH_PETSC_MATRIX_H
21 #define LIBMESH_PETSC_MATRIX_H
23 #include "libmesh/libmesh_common.h"
25 #ifdef LIBMESH_HAVE_PETSC
28 #include "libmesh/sparse_matrix.h"
29 #include "libmesh/petsc_macro.h"
33 #ifdef LIBMESH_HAVE_CXX11_THREAD
42 #undef semiparallel_only
46 #define semiparallel_only() do { if (this->initialized()) { const char * mytype; \
47 MatGetType(_mat,&mytype); \
48 if (!strcmp(mytype, MATSEQAIJ)) \
49 parallel_object_only(); } } while (0)
51 #define semiparallel_only()
57 # define LIBMESH_SAW_I
61 # undef I // Avoid complex.h contamination
70 template <
typename T>
class DenseMatrix;
101 PetscMatrix (
const Parallel::Communicator & comm_in);
112 const Parallel::Communicator & comm_in);
150 const std::vector<numeric_index_type> & n_nz,
151 const std::vector<numeric_index_type> & n_oz,
154 virtual void init ()
override;
168 virtual void clear ()
override;
170 virtual void zero ()
override;
172 virtual void zero_rows (std::vector<numeric_index_type> & rows, T diag_value = 0.0)
override;
174 virtual void close ()
override;
176 virtual void flush ()
override;
202 const T
value)
override;
206 const T
value)
override;
209 const std::vector<numeric_index_type> & rows,
210 const std::vector<numeric_index_type> & cols)
override;
213 const std::vector<numeric_index_type> & dof_indices)
override;
216 const std::vector<numeric_index_type> & brows,
217 const std::vector<numeric_index_type> & bcols)
override;
220 const std::vector<numeric_index_type> & dof_indices)
override
245 virtual bool closed()
const override;
284 std::vector<numeric_index_type> & indices,
285 std::vector<T> & values)
const override;
299 const std::vector<numeric_index_type> & rows,
300 const std::vector<numeric_index_type> & cols,
301 const bool reuse_submatrix)
const override;
318 #ifdef LIBMESH_HAVE_CXX11_THREAD
327 #endif // #ifdef LIBMESH_HAVE_PETSC
328 #endif // LIBMESH_PETSC_MATRIX_H
Threads::spin_mutex _petsc_matrix_mutex
virtual void clear() override
Restores the SparseMatrix<T> to a pristine state.
virtual void _get_submatrix(SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols, const bool reuse_submatrix) const override
This function either creates or re-initializes a matrix called submatrix which is defined by the indi...
virtual void print_matlab(const std::string &name="") const override
Print the contents of the matrix in Matlab's sparse matrix format.
numeric_index_type local_m() const final
Get the number of rows owned by this process.
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
The libMesh namespace provides an interface to certain functionality in the library.
virtual numeric_index_type row_stop() const override
PetscMatrix(const Parallel::Communicator &comm_in)
Constructor; initializes the matrix to be empty, without any structure, i.e.
Defines a dense matrix for use in Finite Element-type computations.
virtual void init() override
Initialize this matrix using the sparsity structure computed by dof_map.
void update_preallocation_and_zero()
Update the sparsity pattern based on dof_map, and set the matrix to zero.
virtual Real linfty_norm() const override
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
virtual void zero() override
Set all entries to 0.
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
Add value to the element (i,j).
void destroy(triangulateio &t, IO_Type)
Frees any memory which has been dynamically allocated by Triangle.
void get_row(numeric_index_type i, std::vector< numeric_index_type > &indices, std::vector< T > &values) const override
Get a row from the matrix.
virtual void zero_rows(std::vector< numeric_index_type > &rows, T diag_value=0.0) override
Sets all row entries to 0 then puts diag_value in the diagonal entry.
std::mutex _petsc_matrix_mutex
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &brows, const std::vector< numeric_index_type > &bcols) override
Add the full matrix dm to the SparseMatrix.
numeric_index_type local_n() const
Get the number of columns owned by this process.
virtual void get_diagonal(NumericVector< T > &dest) const override
Copies the diagonal part of the matrix into dest.
virtual void print_personal(std::ostream &os=libMesh::out) const override
Print the contents of the matrix to the screen with the PETSc viewer.
virtual void set_destroy_mat_on_exit(bool destroy=true)
If set to false, we don't delete the Mat on destruction and allow instead for PETSc to manage it.
PetscMatrixType _mat_type
void swap(PetscMatrix< T > &)
Swaps the internal data pointers of two PetscMatrices, no actual values are swapped.
dof_id_type numeric_index_type
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
Set the element (i,j) to value.
virtual numeric_index_type m() const override
This class provides a nice interface to the PETSc C-based data structures for parallel,...
void reset_preallocation()
Reset matrix to use the original nonzero pattern provided by users.
virtual Real l1_norm() const override
virtual void close() override
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
bool _destroy_mat_on_exit
This boolean value should only be set to false for the constructor which takes a PETSc Mat object.
void get_local_size(numeric_index_type &m, numeric_index_type &n) const
Get the number of rows and columns owned by this process.
virtual bool closed() const override
virtual numeric_index_type n() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void flush() override
For PETSc matrix , this function is similar to close but without shrinking memory.
PetscMatrix & operator=(const PetscMatrix &)=delete
virtual void get_transpose(SparseMatrix< T > &dest) const override
Copies the transpose of the matrix into dest, which may be *this.
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
Add the full matrix dm to the SparseMatrix.
virtual numeric_index_type row_start() const override
Mat _mat
PETSc matrix datatype to store values.
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices) override
Same as add_block_matrix(), but assumes the row and column maps are the same.
void set_matrix_type(PetscMatrixType mat_type)