libMesh
petsc_matrix_base.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_BASE_H
21 #define LIBMESH_PETSC_MATRIX_BASE_H
22 
23 #include "libmesh/libmesh_common.h"
24 
25 #ifdef LIBMESH_HAVE_PETSC
26 
27 // Local includes
28 #include "libmesh/sparse_matrix.h"
29 #include "libmesh/petsc_macro.h"
30 #include "libmesh/petsc_solver_exception.h"
31 #include "libmesh/wrapped_petsc.h"
32 
33 // Petsc include files.
34 #ifdef I
35 # define LIBMESH_SAW_I
36 #endif
37 #include <petscmat.h>
38 #ifndef LIBMESH_SAW_I
39 # undef I // Avoid complex.h contamination
40 #endif
41 
42 // Macro to identify and debug functions which should be called in
43 // parallel on parallel matrices but which may be called in serial on
44 // serial matrices. This macro will only be valid inside non-static
45 // PetscMatrix methods
46 #undef semiparallel_only
47 #undef exceptionless_semiparallel_only
48 #ifndef NDEBUG
49 #include <cstring>
50 
51 #define semiparallel_only() do { if (this->initialized()) { const char * mytype; \
52  LibmeshPetscCall(MatGetType(this->_mat,&mytype)); \
53  if (!strcmp(mytype, MATSEQAIJ)) \
54  parallel_object_only(); } } while (0)
55 #define exceptionless_semiparallel_only() do { if (this->initialized()) { const char * mytype; \
56  auto semiparallel_only_ierr = MatGetType(this->_mat,&mytype); \
57  libmesh_ignore(semiparallel_only_ierr); \
58  if (!strcmp(mytype, MATSEQAIJ)) \
59  exceptionless_parallel_object_only(); } } while (0)
60 #else
61 #define semiparallel_only()
62 #define exceptionless_semiparallel_only()
63 #endif
64 
65 
66 namespace libMesh
67 {
68 
78 template <typename T>
79 class PetscMatrixBase : public SparseMatrix<T>
80 {
81 public:
82  explicit
83  PetscMatrixBase (const Parallel::Communicator & comm_in);
84 
92  explicit
93  PetscMatrixBase (Mat m,
94  const Parallel::Communicator & comm_in,
95  bool destroy_on_exit = false);
96 
102  PetscMatrixBase (PetscMatrixBase &&) = delete;
103  PetscMatrixBase (const PetscMatrixBase &) = delete;
105  virtual ~PetscMatrixBase ();
106 
107  virtual SolverPackage solver_package() override
108  {
109  return PETSC_SOLVERS;
110  }
111 
120  Mat mat () { libmesh_assert (_mat); return _mat; }
121  Mat mat() const { libmesh_assert(_mat); return _mat; }
122 
126  virtual void clear() noexcept override;
127 
132  void set_destroy_mat_on_exit(bool destroy = true);
133 
138  void swap (PetscMatrixBase<T> &);
139 
140  using SparseMatrix<T>::operator=;
141 
145  void set_context();
146 
150  static PetscMatrixBase<T> * get_context(Mat mat, const TIMPI::Communicator & comm);
151 
152  virtual numeric_index_type m () const override;
153 
154  virtual numeric_index_type local_m () const final;
155 
156  virtual numeric_index_type n () const override;
157 
161  virtual numeric_index_type local_n () const final;
162 
163  virtual numeric_index_type row_start () const override;
164 
165  virtual numeric_index_type row_stop () const override;
166 
167  virtual numeric_index_type col_start () const override;
168 
169  virtual numeric_index_type col_stop () const override;
170 
171  virtual void close () override;
172 
173  virtual bool closed() const override;
174 
175 protected:
176 
180  Mat _mat;
181 
187 };
188 
189 } // namespace libMesh
190 
191 #endif // #ifdef LIBMESH_HAVE_PETSC
192 #endif // LIBMESH_PETSC_MATRIX_BASE_H
virtual numeric_index_type n() const override
virtual numeric_index_type local_n() const final
Get the number of columns owned by this process.
bool _destroy_mat_on_exit
This boolean value should only be set to false for the constructor which takes a PETSc Mat object...
virtual void clear() noexcept override
clear() is called from the destructor, so it should not throw.
virtual numeric_index_type row_stop() const override
virtual numeric_index_type col_stop() const override
virtual numeric_index_type local_m() const final
Get the number of rows owned by this process.
This class provides a nice interface to the PETSc C-based data structures for parallel, sparse matrices.
void swap(PetscMatrixBase< T > &)
Swaps the internal data pointers of two PetscMatrices, no actual values are swapped.
virtual numeric_index_type m() const override
const Parallel::Communicator & comm() const
static PetscMatrixBase< T > * get_context(Mat mat, const TIMPI::Communicator &comm)
The libMesh namespace provides an interface to certain functionality in the library.
void set_context()
Set the context (ourself) for _mat.
Generic sparse matrix.
Definition: vector_fe_ex5.C:46
dof_id_type numeric_index_type
Definition: id_types.h:99
PetscMatrixBase & operator=(PetscMatrixBase &&)=delete
void set_destroy_mat_on_exit(bool destroy=true)
If set to false, we don&#39;t delete the Mat on destruction and allow instead for PETSc to manage it...
libmesh_assert(ctx)
virtual void close() override
Calls the SparseMatrix&#39;s internal assembly routines, ensuring that the values are consistent across p...
virtual numeric_index_type col_start() const override
Mat _mat
PETSc matrix datatype to store values.
PetscMatrixBase(const Parallel::Communicator &comm_in)
virtual bool closed() const override
virtual SolverPackage solver_package() override
virtual numeric_index_type row_start() const override
SolverPackage
Defines an enum for various linear solver packages.
void destroy(triangulateio &t, IO_Type)
Frees any memory which has been dynamically allocated by Triangle.