libMesh
petsc_matrix.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_PETSC_MATRIX_H
21 #define LIBMESH_PETSC_MATRIX_H
22 
23 #include "libmesh/libmesh_common.h"
24 
25 #ifdef LIBMESH_HAVE_PETSC
26 
27 // Local includes
28 #include "libmesh/petsc_matrix_base.h"
29 #include "libmesh/petsc_macro.h"
30 #include "libmesh/petsc_solver_exception.h"
31 
32 // C++ includes
33 #include <algorithm>
34 #ifdef LIBMESH_HAVE_CXX11_THREAD
35 #include <atomic>
36 #include <mutex>
37 #endif
38 
39 class PetscMatrixTest;
40 
41 namespace libMesh
42 {
43 
44 // Forward Declarations
45 template <typename T> class DenseMatrix;
46 
47 enum PetscMatrixType : int {
48  AIJ=0,
50 
60 template <typename T>
61 class PetscMatrix final : public PetscMatrixBase<T>
62 {
63 public:
74  explicit
75  PetscMatrix (const Parallel::Communicator & comm_in);
76 
84  explicit
85  PetscMatrix (Mat m,
86  const Parallel::Communicator & comm_in,
87  bool destroy_on_exit = false);
88 
93  explicit
94  PetscMatrix (const Parallel::Communicator & comm_in,
95  const numeric_index_type m,
96  const numeric_index_type n,
97  const numeric_index_type m_l,
98  const numeric_index_type n_l,
99  const numeric_index_type n_nz=30,
100  const numeric_index_type n_oz=10,
101  const numeric_index_type blocksize=1);
102 
108  PetscMatrix (PetscMatrix &&) = delete;
109  PetscMatrix (const PetscMatrix &) = delete;
110  PetscMatrix & operator= (PetscMatrix &&) = delete;
111  virtual ~PetscMatrix ();
112 
114  virtual SparseMatrix<T> & operator= (const SparseMatrix<T> & v) override;
115 
127  virtual void init (const numeric_index_type m,
128  const numeric_index_type n,
129  const numeric_index_type m_l,
130  const numeric_index_type n_l,
131  const numeric_index_type n_nz=30,
132  const numeric_index_type n_oz=10,
133  const numeric_index_type blocksize=1) override;
134 
146  void init (const numeric_index_type m,
147  const numeric_index_type n,
148  const numeric_index_type m_l,
149  const numeric_index_type n_l,
150  const std::vector<numeric_index_type> & n_nz,
151  const std::vector<numeric_index_type> & n_oz,
152  const numeric_index_type blocksize=1);
153 
154  virtual void init (ParallelType = PARALLEL) override;
155 
162 
166  void reset_preallocation();
167 
168  virtual void zero () override;
169 
170  virtual std::unique_ptr<SparseMatrix<T>> zero_clone () const override;
171 
172  virtual std::unique_ptr<SparseMatrix<T>> clone () const override;
173 
174  virtual void zero_rows (std::vector<numeric_index_type> & rows, T diag_value = 0.0) override;
175 
176  virtual void flush () override;
177 
178  virtual void set (const numeric_index_type i,
179  const numeric_index_type j,
180  const T value) override;
181 
182  virtual void add (const numeric_index_type i,
183  const numeric_index_type j,
184  const T value) override;
185 
186  virtual void add_matrix (const DenseMatrix<T> & dm,
187  const std::vector<numeric_index_type> & rows,
188  const std::vector<numeric_index_type> & cols) override;
189 
190  virtual void add_matrix (const DenseMatrix<T> & dm,
191  const std::vector<numeric_index_type> & dof_indices) override;
192 
193  virtual void add_block_matrix (const DenseMatrix<T> & dm,
194  const std::vector<numeric_index_type> & brows,
195  const std::vector<numeric_index_type> & bcols) override;
196 
197  virtual void add_block_matrix (const DenseMatrix<T> & dm,
198  const std::vector<numeric_index_type> & dof_indices) override
199  { this->add_block_matrix (dm, dof_indices, dof_indices); }
200 
214  virtual void add (const T a, const SparseMatrix<T> & X) override;
215 
221  virtual void matrix_matrix_mult (SparseMatrix<T> & X, SparseMatrix<T> & Y, bool reuse = false) override;
222 
227  virtual
228  void add_sparse_matrix (const SparseMatrix<T> & spm,
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;
232 
233  virtual T operator () (const numeric_index_type i,
234  const numeric_index_type j) const override;
235 
236  virtual Real l1_norm () const override;
237 
238  virtual Real frobenius_norm () const;
239 
240  virtual Real linfty_norm () const override;
241 
248  virtual void print_personal(std::ostream & os=libMesh::out) const override;
249 
250  virtual void print_matlab(const std::string & name = "") const override;
251 
252  virtual void print_petsc_binary(const std::string & filename) override;
253 
254  virtual void print_petsc_hdf5(const std::string & filename) override;
255 
256  virtual void read_petsc_binary(const std::string & filename) override;
257 
258  virtual void read_petsc_hdf5(const std::string & filename) override;
259 
260  virtual void get_diagonal (NumericVector<T> & dest) const override;
261 
262  virtual void get_transpose (SparseMatrix<T> & dest) const override;
263 
264  virtual
266  std::vector<numeric_index_type> & indices,
267  std::vector<T> & values) const override;
268 
277  virtual void create_submatrix_nosort(SparseMatrix<T> & submatrix,
278  const std::vector<numeric_index_type> & rows,
279  const std::vector<numeric_index_type> & cols) const override;
280 
281  virtual void scale(const T scale) override;
282 
283 #if PETSC_RELEASE_GREATER_EQUALS(3, 23, 0)
284 
289  std::unique_ptr<PetscMatrix<T>> copy_from_hash();
290 #endif
291 
292  virtual bool supports_hash_table() const override;
293 
294  virtual void restore_original_nonzero_pattern() override;
295 
296 protected:
307  numeric_index_type m_l,
308  numeric_index_type n_l,
309  numeric_index_type blocksize);
310 
311  /*
312  * Performs matrix preallcation
313  * \param m_l The local number of rows.
314  * \param n_nz array containing the number of nonzeros in each row of the DIAGONAL portion of the local submatrix.
315  * \param n_oz Array containing the number of nonzeros in each row of the OFF-DIAGONAL portion of the local submatrix.
316  * \param blocksize Optional value indicating dense coupled blocks for systems with multiple variables all of the same */
318  const std::vector<numeric_index_type> & n_nz,
319  const std::vector<numeric_index_type> & n_oz,
320  numeric_index_type blocksize);
321 
328  void finish_initialization();
329 
340  virtual void _get_submatrix(SparseMatrix<T> & submatrix,
341  const std::vector<numeric_index_type> & rows,
342  const std::vector<numeric_index_type> & cols,
343  const bool reuse_submatrix) const override;
344 
345  void _petsc_viewer(const std::string & filename,
346  PetscViewerType viewertype,
347  PetscFileMode filemode);
348 
350 
351 private:
352 #ifdef LIBMESH_HAVE_CXX11_THREAD
353  mutable std::mutex _petsc_matrix_mutex;
354 #else
356 #endif
357 
364  template <NormType N> Real norm () const;
365 
366  friend class ::PetscMatrixTest;
367 };
368 
369 } // namespace libMesh
370 
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.
Definition: elem_quality.C:42
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.
Definition: petsc_matrix.h:197
void reset_preallocation()
Reset matrix to use the original nonzero pattern provided by users.
Definition: petsc_matrix.C:385
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&#39;s HDF5 sparse matrix format. ...
Definition: petsc_matrix.C:723
virtual void print_matlab(const std::string &name="") const override
Print the contents of the matrix in Matlab&#39;s sparse matrix format.
Definition: petsc_matrix.C:510
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.
Definition: petsc_matrix.C:401
virtual Real frobenius_norm() const
Definition: petsc_matrix.C:497
PetscMatrix(const Parallel::Communicator &comm_in)
Constructor; initializes the matrix to be empty, without any structure, i.e.
Definition: petsc_matrix.C:86
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.
Definition: petsc_matrix.C:135
virtual void print_petsc_binary(const std::string &filename) override
Write the contents of the matrix to a file in PETSc&#39;s binary sparse matrix format.
Definition: petsc_matrix.C:693
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)
Definition: petsc_matrix.C:260
virtual numeric_index_type m() const override
Threads::spin_mutex _petsc_matrix_mutex
Definition: petsc_matrix.h:355
PetscMatrix & operator=(PetscMatrix &&)=delete
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:44
void update_preallocation_and_zero()
Update the sparsity pattern based on dof_map, and set the matrix to zero.
Definition: petsc_matrix.C:379
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.
Definition: petsc_matrix.C:759
virtual void print_personal(std::ostream &os=libMesh::out) const override
Print the contents of the matrix to the screen with the PETSc viewer.
Definition: petsc_matrix.C:567
virtual void read_petsc_binary(const std::string &filename) override
Read the contents of the matrix from a file in PETSc&#39;s binary sparse matrix format.
Definition: petsc_matrix.C:713
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const override
Definition: petsc_matrix.C:435
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
Definition: petsc_matrix.C:454
Generic sparse matrix.
Definition: vector_fe_ex5.C:46
std::mutex _petsc_matrix_mutex
Definition: petsc_matrix.h:353
dof_id_type numeric_index_type
Definition: id_types.h:99
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...
Definition: petsc_matrix.C:851
void _petsc_viewer(const std::string &filename, PetscViewerType viewertype, PetscFileMode filemode)
Definition: petsc_matrix.C:663
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).
Definition: petsc_matrix.C:987
PetscMatrixType _mat_type
Definition: petsc_matrix.h:349
virtual void get_diagonal(NumericVector< T > &dest) const override
Copies the diagonal part of the matrix into dest.
Definition: petsc_matrix.C:919
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&#39;s HDF5 sparse matrix format.
Definition: petsc_matrix.C:703
virtual void get_transpose(SparseMatrix< T > &dest) const override
Copies the transpose of the matrix into dest, which may be *this.
Definition: petsc_matrix.C:931
virtual Real l1_norm() const override
Definition: petsc_matrix.C:492
OStreamProxy out
static const bool value
Definition: xdr_io.C:54
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.
Definition: petsc_matrix.C:733
This class provides a nice interface to the PETSc C-based AIJ data structures for parallel...
Definition: petsc_matrix.h:61
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...
Definition: petsc_matrix.C:800
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.
Definition: petsc_matrix.C:416
Defines a dense matrix for use in Finite Element-type computations.
Definition: dof_map.h:75
void finish_initialization()
Finish up the initialization process.
Definition: petsc_matrix.C:329
virtual void flush() override
For PETSc matrix , this function is similar to close but without shrinking memory.
Definition: petsc_matrix.C:960
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.
Definition: petsc_matrix.C:214
virtual Real linfty_norm() const override
Definition: petsc_matrix.C:502
ParallelType
Defines an enum for parallel data structure types.