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/petsc_matrix_base.h" 29 #include "libmesh/petsc_macro.h" 30 #include "libmesh/petsc_solver_exception.h" 34 #ifdef LIBMESH_HAVE_CXX11_THREAD 45 template <
typename T>
class DenseMatrix;
87 bool destroy_on_exit =
false);
150 const std::vector<numeric_index_type> & n_nz,
151 const std::vector<numeric_index_type> & n_oz,
168 virtual void zero ()
override;
170 virtual std::unique_ptr<SparseMatrix<T>>
zero_clone ()
const override;
172 virtual std::unique_ptr<SparseMatrix<T>>
clone ()
const override;
174 virtual void zero_rows (std::vector<numeric_index_type> & rows, T diag_value = 0.0)
override;
176 virtual void flush ()
override;
180 const T
value)
override;
184 const T
value)
override;
187 const std::vector<numeric_index_type> & rows,
188 const std::vector<numeric_index_type> & cols)
override;
191 const std::vector<numeric_index_type> & dof_indices)
override;
194 const std::vector<numeric_index_type> & brows,
195 const std::vector<numeric_index_type> & bcols)
override;
198 const std::vector<numeric_index_type> & dof_indices)
override 229 const std::map<numeric_index_type,numeric_index_type> & row_ltog,
230 const std::map<numeric_index_type,numeric_index_type> & col_ltog,
231 const T scalar)
override;
266 std::vector<numeric_index_type> & indices,
267 std::vector<T> & values)
const override;
278 const std::vector<numeric_index_type> & rows,
279 const std::vector<numeric_index_type> & cols)
const override;
283 #if PETSC_RELEASE_GREATER_EQUALS(3, 23, 0) 318 const std::vector<numeric_index_type> & n_nz,
319 const std::vector<numeric_index_type> & n_oz,
341 const std::vector<numeric_index_type> & rows,
342 const std::vector<numeric_index_type> & cols,
343 const bool reuse_submatrix)
const override;
346 PetscViewerType viewertype,
347 PetscFileMode filemode);
352 #ifdef LIBMESH_HAVE_CXX11_THREAD 364 template <NormType N>
Real norm ()
const;
366 friend class ::PetscMatrixTest;
371 #endif // #ifdef LIBMESH_HAVE_PETSC 372 #endif // LIBMESH_PETSC_MATRIX_H std::string name(const ElemQuality q)
This function returns a string containing some name for q.
virtual numeric_index_type n() const override
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 reset_preallocation()
Reset matrix to use the original nonzero pattern provided by users.
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
virtual void read_petsc_hdf5(const std::string &filename) override
Read the contents of the matrix from a file in PETSc's HDF5 sparse matrix format. ...
virtual void print_matlab(const std::string &name="") const override
Print the contents of the matrix in Matlab's sparse matrix format.
This class provides a nice interface to the PETSc C-based data structures for parallel, sparse matrices.
virtual void zero() override
Set all entries to 0.
virtual Real frobenius_norm() const
PetscMatrix(const Parallel::Communicator &comm_in)
Constructor; initializes the matrix to be empty, without any structure, i.e.
void init_without_preallocation(numeric_index_type m, numeric_index_type n, numeric_index_type m_l, numeric_index_type n_l, numeric_index_type blocksize)
Perform matrix initialization steps sans preallocation.
virtual void print_petsc_binary(const std::string &filename) override
Write the contents of the matrix to a file in PETSc's binary sparse matrix format.
void preallocate(numeric_index_type m_l, const std::vector< numeric_index_type > &n_nz, const std::vector< numeric_index_type > &n_oz, numeric_index_type blocksize)
virtual numeric_index_type m() const override
Threads::spin_mutex _petsc_matrix_mutex
PetscMatrix & operator=(PetscMatrix &&)=delete
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
void update_preallocation_and_zero()
Update the sparsity pattern based on dof_map, and set the matrix to zero.
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.
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 read_petsc_binary(const std::string &filename) override
Read the contents of the matrix from a file in PETSc's binary sparse matrix format.
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const override
The libMesh namespace provides an interface to certain functionality in the library.
virtual void add_sparse_matrix(const SparseMatrix< T > &spm, const std::map< numeric_index_type, numeric_index_type > &row_ltog, const std::map< numeric_index_type, numeric_index_type > &col_ltog, const T scalar) override
Add scalar* spm to the rows and cols of this matrix (A): A(rows[i], cols[j]) += scalar * spm(i...
virtual std::unique_ptr< SparseMatrix< T > > clone() const override
std::mutex _petsc_matrix_mutex
dof_id_type numeric_index_type
std::unique_ptr< PetscMatrix< T > > copy_from_hash()
Creates a copy of the current hash table matrix and then performs assembly.
virtual bool supports_hash_table() const override
virtual 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 create_submatrix_nosort(SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const override
Similar to the create_submatrix function, this function creates a submatrix which is defined by the i...
void _petsc_viewer(const std::string &filename, PetscViewerType viewertype, PetscFileMode filemode)
virtual void matrix_matrix_mult(SparseMatrix< T > &X, SparseMatrix< T > &Y, bool reuse=false) override
Compute Y = A*X for matrix X.
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
Add value to the element (i,j).
PetscMatrixType _mat_type
virtual void get_diagonal(NumericVector< T > &dest) const override
Copies the diagonal part of the matrix into dest.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void print_petsc_hdf5(const std::string &filename) override
Write the contents of the matrix to a file in PETSc's HDF5 sparse matrix format.
virtual void get_transpose(SparseMatrix< T > &dest) const override
Copies the transpose of the matrix into dest, which may be *this.
virtual Real l1_norm() const override
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.
This class provides a nice interface to the PETSc C-based AIJ data structures for parallel...
virtual void scale(const T scale) override
Scales all elements of this matrix by scale.
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 restore_original_nonzero_pattern() override
Reset the memory storage of 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.
Defines a dense matrix for use in Finite Element-type computations.
void finish_initialization()
Finish up the initialization process.
virtual void flush() override
For PETSc matrix , this function is similar to close but without shrinking memory.
virtual void init(const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const numeric_index_type n_nz=30, const numeric_index_type n_oz=10, const numeric_index_type blocksize=1) override
Initialize a PETSc matrix.
virtual Real linfty_norm() const override
ParallelType
Defines an enum for parallel data structure types.