18 #ifndef LIBMESH_PETSC_DM_WRAPPER_H 19 #define LIBMESH_PETSC_DM_WRAPPER_H 21 #include "libmesh/petsc_macro.h" 23 #ifdef LIBMESH_HAVE_PETSC 24 #if !PETSC_VERSION_LESS_THAN(3,7,3) 25 #if defined(LIBMESH_ENABLE_AMR) && defined(LIBMESH_HAVE_METAPHYSICL) 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" 35 # define LIBMESH_SAW_I 39 # undef I // Avoid complex.h contamination 45 #include <unordered_map> 68 std::vector<std::vector<numeric_index_type>>
dof_vec;
122 std::vector<WrappedPetsc<DM>>
_dms;
148 std::vector<std::unique_ptr<PetscMatrixBase<Number>>>
_pmtx_vec;
169 std::vector<std::unique_ptr<PetscVector<Number>>>
_vec_vec;
185 { libmesh_assert_less(level,
_dms.size());
186 return *
_dms[level]; }
193 { libmesh_assert_less(level,
_sections.size());
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);
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);
264 PetscSection & section);
269 #endif // #if LIBMESH_ENABLE_AMR && LIBMESH_HAVE_METAPHYSICL 270 #endif // #if PETSC_VERSION 271 #endif // #ifdef LIBMESH_HAVE_PETSC 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 §ion)
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 §ion, 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 §ion, 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.
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 §ion)
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...
void build_section(const System &system, PetscSection §ion)
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.