petsc_dm_wrapper.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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/libmesh_common.h"
22 
23 #ifdef LIBMESH_HAVE_PETSC
24 
25 #include <vector>
26 #include <memory>
27 #include <unordered_map>
28 #include <map>
29 
30 // PETSc includes
32 #include <petsc.h>
34 
35 namespace libMesh
36 {
37  // Forward declarations
38  class System;
39  class DofObject;
40 
52 {
53 public:
54 
55  PetscDMWrapper() = default;
56 
58 
60  void clear();
61 
62  void init_and_attach_petscdm(System & system, SNES & snes);
63 
64 private:
65 
67  std::vector<std::unique_ptr<DM>> _dms;
68 
70  std::vector<std::unique_ptr<PetscSection>> _sections;
71 
73  std::vector<std::unique_ptr<PetscSF>> _star_forests;
74 
76  void init_dm_data(unsigned int n_levels);
77 
79 
82  DM & get_dm(unsigned int level)
83  { libmesh_assert(level < _dms.size());
84  return *(_dms[level].get()); }
85 
87 
90  PetscSection & get_section(unsigned int level)
91  { libmesh_assert(level < _sections.size());
92  return *(_sections[level].get()); }
93 
95 
98  PetscSF & get_star_forest(unsigned int level)
99  { libmesh_assert(level < _star_forests.size());
100  return *(_star_forests[level].get()); }
101 
103 
111  void build_section(const System & system, PetscSection & section);
112 
114 
125  void build_sf( const System & system, PetscSF & star_forest );
126 
128 
135  void set_point_range_in_section( const System & system,
136  PetscSection & section,
137  std::unordered_map<dof_id_type,dof_id_type> & node_map,
138  std::unordered_map<dof_id_type,dof_id_type> & elem_map,
139  std::map<dof_id_type,unsigned int> & scalar_map);
140 
142 
145  void add_dofs_to_section (const System & system,
146  PetscSection & section,
147  const std::unordered_map<dof_id_type,dof_id_type> & node_map,
148  const std::unordered_map<dof_id_type,dof_id_type> & elem_map,
149  const std::map<dof_id_type,unsigned int> & scalar_map);
150 
152 
156  dof_id_type check_section_n_dofs( const System & system, PetscSection & section );
157 
159  void add_dofs_helper (const System & system,
160  const DofObject & dof_object,
161  dof_id_type local_id,
162  PetscSection & section);
163 
164 };
165 
166 }
167 
168 #endif // #ifdef LIBMESH_HAVE_PETSC
169 
170 #endif // LIBMESH_PETSC_DM_WRAPPER_H
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.
dof_id_type check_section_n_dofs(const System &system, PetscSection &section)
Helper function to sanity check PetscSection construction.
DM & get_dm(unsigned int level)
Get reference to DM for the given mesh level.
std::vector< std::unique_ptr< PetscSF > > _star_forests
Vector of star forests for all grid levels.
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.
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:92
std::vector< std::unique_ptr< PetscSection > > _sections
Vector of PETScSections for all grid levels.
unsigned int n_levels(const MeshBase &mesh)
Definition: mesh_tools.C:653
PetscSection & get_section(unsigned int level)
Get reference to PetscSection for the given mesh level.
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)
Helper function to reduce code duplication when setting dofs in section.
std::vector< std::unique_ptr< DM > > _dms
Vector of DMs for all grid levels.
void init_dm_data(unsigned int n_levels)
Init all the n_mesh_level dependent data structures.
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.
void clear()
Destroys and clears all build DM-related data.
void init_and_attach_petscdm(System &system, SNES &snes)
uint8_t dof_id_type
Definition: id_types.h:64