libMesh
petsc_dm_wrapper.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 #ifndef LIBMESH_PETSC_DM_WRAPPER_H
19 #define LIBMESH_PETSC_DM_WRAPPER_H
20 
21 #include "libmesh/petsc_macro.h"
22 
23 #ifdef LIBMESH_HAVE_PETSC
24 #if !PETSC_VERSION_LESS_THAN(3,7,3)
25 #if defined(LIBMESH_ENABLE_AMR) && defined(LIBMESH_HAVE_METAPHYSICL)
26 
27 // libMesh includes
28 #include "libmesh/petsc_macro.h"
29 #include "libmesh/petsc_matrix_base.h"
30 #include "libmesh/petsc_vector.h"
31 #include "libmesh/wrapped_petsc.h"
32 
33 // PETSc includes
34 #ifdef I
35 # define LIBMESH_SAW_I
36 #endif
37 #include <petsc.h>
38 #ifndef LIBMESH_SAW_I
39 # undef I // Avoid complex.h contamination
40 #endif
41 
42 // C++ includes
43 #include <vector>
44 #include <memory>
45 #include <unordered_map>
46 #include <map>
47 
48 namespace libMesh
49 {
50  // Forward declarations
51  class System;
52  class DofObject;
53 
56  {
57  int n_dofs;
58  int mesh_dim;
59  DM * coarser_dm;
60  DM * finer_dm;
61  DM * global_dm;
66 
68  std::vector<std::vector<numeric_index_type>> dof_vec;
69 
71  std::vector<PetscInt> subfields;
72 
74  n_dofs(-12345),
75  mesh_dim(-12345),
76  coarser_dm(nullptr),
77  finer_dm(nullptr),
78  global_dm(nullptr),
79  K_interp_ptr(nullptr),
80  K_sub_interp_ptr(nullptr),
81  K_restrict_ptr(nullptr),
82  current_vec(nullptr)
83  {}
84 
85  };
86 
98  {
99  public:
100 
101  PetscDMWrapper() = default;
102 
103  ~PetscDMWrapper();
104 
106  void clear();
107 
108  void init_and_attach_petscdm(System & system, SNES snes);
109  void init_and_attach_petscdm(System & system, KSP ksp);
110 
111  private:
115  unsigned int init_petscdm(System & system);
116 
122  std::vector<WrappedPetsc<DM>> _dms;
123 
132  std::vector<WrappedPetsc<PetscSection>> _sections;
133 
141  std::vector<PetscSF> _star_forests;
142 
148  std::vector<std::unique_ptr<PetscMatrixBase<Number>>> _pmtx_vec;
149 
155  std::vector<std::unique_ptr<PetscMatrixBase<Number>>> _subpmtx_vec;
156 
162  std::vector<PetscDMContext> _ctx_vec;
163 
169  std::vector<std::unique_ptr<PetscVector<Number>>> _vec_vec;
170 
172  std::vector<unsigned int> _mesh_dof_sizes;
173 
175  std::vector<unsigned int> _mesh_dof_loc_sizes;
176 
178  void init_dm_data(unsigned int n_levels, const Parallel::Communicator & comm);
179 
184  DM & get_dm(unsigned int level)
185  { libmesh_assert_less(level, _dms.size());
186  return *_dms[level]; }
187 
192  PetscSection & get_section(unsigned int level)
193  { libmesh_assert_less(level, _sections.size());
194  return *_sections[level]; }
195 
200  PetscSF & get_star_forest(unsigned int level)
201  { libmesh_assert_less(level, _star_forests.size());
202  return _star_forests[level]; }
203 
213  void build_section(const System & system, PetscSection & section);
214 
227  void build_sf( const System & system, PetscSF & star_forest );
228 
237  void set_point_range_in_section( const System & system,
238  PetscSection & section,
239  std::unordered_map<dof_id_type,dof_id_type> & node_map,
240  std::unordered_map<dof_id_type,dof_id_type> & elem_map,
241  std::map<dof_id_type,unsigned int> & scalar_map);
242 
247  void add_dofs_to_section (const System & system,
248  PetscSection & section,
249  const std::unordered_map<dof_id_type,dof_id_type> & node_map,
250  const std::unordered_map<dof_id_type,dof_id_type> & elem_map,
251  const std::map<dof_id_type,unsigned int> & scalar_map);
252 
258  dof_id_type check_section_n_dofs( PetscSection & section );
259 
260  // Helper function to reduce code duplication when setting dofs in section
261  void add_dofs_helper (const System & system,
262  const DofObject & dof_object,
263  dof_id_type local_id,
264  PetscSection & section);
265  };
266 
267 }
268 
269 #endif // #if LIBMESH_ENABLE_AMR && LIBMESH_HAVE_METAPHYSICL
270 #endif // #if PETSC_VERSION
271 #endif // #ifdef LIBMESH_HAVE_PETSC
272 
273 #endif // LIBMESH_PETSC_DM_WRAPPER_H
void init_dm_data(unsigned int n_levels, const Parallel::Communicator &comm)
Init all the n_mesh_level dependent data structures.
dof_id_type check_section_n_dofs(PetscSection &section)
Helper function to sanity check PetscSection construction The PetscSection contains local dof informa...
std::vector< std::unique_ptr< PetscVector< Number > > > _vec_vec
Vector of solution vectors for all grid levels.
std::vector< PetscDMContext > _ctx_vec
Vector of internal PetscDM context structs for all grid levels Pointers to these C++ objects are pass...
void set_point_range_in_section(const System &system, PetscSection &section, std::unordered_map< dof_id_type, dof_id_type > &node_map, std::unordered_map< dof_id_type, dof_id_type > &elem_map, std::map< dof_id_type, unsigned int > &scalar_map)
Helper function for build_section.
DM & get_dm(unsigned int level)
Get reference to DM for the given mesh level.
std::vector< WrappedPetsc< DM > > _dms
Vector of DMs for all grid levels.
unsigned int init_petscdm(System &system)
Initialize the PETSc DM and return the number of geometric multigrid levels.
std::vector< PetscInt > subfields
Stores subfield ids for use in subprojection matrixes on coarser DMs.
void add_dofs_to_section(const System &system, PetscSection &section, const std::unordered_map< dof_id_type, dof_id_type > &node_map, const std::unordered_map< dof_id_type, dof_id_type > &elem_map, const std::map< dof_id_type, unsigned int > &scalar_map)
Helper function for build_section.
PetscMatrixBase< libMesh::Number > * K_interp_ptr
The libMesh namespace provides an interface to certain functionality in the library.
This class defines a wrapper around the PETSc DM infrastructure.
PetscMatrixBase< libMesh::Number > * K_restrict_ptr
Struct to house data regarding where in the mesh hierarchy we are located.
std::vector< std::unique_ptr< PetscMatrixBase< Number > > > _pmtx_vec
Vector of projection matrixes for all grid levels.
std::vector< WrappedPetsc< PetscSection > > _sections
Vector of PETScSections for all grid levels.
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
unsigned int n_levels(const MeshBase &mesh)
Definition: mesh_tools.C:802
std::vector< std::vector< numeric_index_type > > dof_vec
Stores local dofs for each var for use in subprojection matrixes.
PetscSection & get_section(unsigned int level)
Get reference to PetscSection for the given mesh level.
std::vector< unsigned int > _mesh_dof_loc_sizes
Stores n_local_dofs for each grid level, to be used for projection vector sizing. ...
std::vector< PetscSF > _star_forests
Vector of star forests for all grid levels.
PetscSF & get_star_forest(unsigned int level)
Get reference to PetscSF for the given mesh level.
void add_dofs_helper(const System &system, const DofObject &dof_object, dof_id_type local_id, PetscSection &section)
PetscMatrixBase< libMesh::Number > * K_sub_interp_ptr
void init_and_attach_petscdm(System &system, SNES snes)
std::vector< unsigned int > _mesh_dof_sizes
Stores n_dofs for each grid level, to be used for projection matrix sizing.
The DofObject defines an abstract base class for objects that have degrees of freedom associated with...
Definition: dof_object.h:54
void build_section(const System &system, PetscSection &section)
Takes System, empty PetscSection and populates the PetscSection.
void build_sf(const System &system, PetscSF &star_forest)
Takes System, empty PetscSF and populates the PetscSF.
PetscVector< libMesh::Number > * current_vec
std::vector< std::unique_ptr< PetscMatrixBase< Number > > > _subpmtx_vec
Vector of sub projection matrixes for all grid levels for fieldsplit.
void clear()
Destroys and clears all build DM-related data.
uint8_t dof_id_type
Definition: id_types.h:67