Go to the documentation of this file.
   18 #ifndef LIBMESH_PETSC_DM_WRAPPER_H 
   19 #define LIBMESH_PETSC_DM_WRAPPER_H 
   21 #include "libmesh/libmesh_common.h" 
   22 #include "libmesh/petsc_macro.h" 
   24 #ifdef LIBMESH_HAVE_PETSC 
   25 #if !PETSC_VERSION_LESS_THAN(3,7,3) 
   26 #if defined(LIBMESH_ENABLE_AMR) && defined(LIBMESH_HAVE_METAPHYSICL) 
   30 #include <unordered_map> 
   34 #include "libmesh/ignore_warnings.h" 
   36 # define LIBMESH_SAW_I 
   40 # undef I // Avoid complex.h contamination 
   42 #include "libmesh/restore_warnings.h" 
   43 #include "libmesh/petsc_matrix.h" 
   44 #include "libmesh/petsc_vector.h" 
   66     std::vector<std::vector<numeric_index_type>> 
dof_vec;
 
  111     std::vector<std::unique_ptr<DM>> 
_dms;
 
  120     std::vector<std::unique_ptr<PetscMatrix<Number>>> 
_pmtx_vec;
 
  126     std::vector<std::unique_ptr<PetscDMContext>> 
_ctx_vec;
 
  129     std::vector<std::unique_ptr<PetscVector<Number>>> 
_vec_vec;
 
  145       { libmesh_assert_less(level, 
_dms.size());
 
  146         return *(
_dms[level].get()); }
 
  153       { libmesh_assert_less(level, 
_sections.size());
 
  198                                      PetscSection & section,
 
  199                                      std::unordered_map<dof_id_type,dof_id_type> & node_map,
 
  200                                      std::unordered_map<dof_id_type,dof_id_type> & elem_map,
 
  201                                      std::map<dof_id_type,unsigned int> & scalar_map);
 
  208                               PetscSection & section,
 
  209                               const std::unordered_map<dof_id_type,dof_id_type> & node_map,
 
  210                               const std::unordered_map<dof_id_type,dof_id_type> & elem_map,
 
  211                               const std::map<dof_id_type,unsigned int> & scalar_map);
 
  224                           PetscSection & section);
 
  230 #endif // #if LIBMESH_ENABLE_AMR && LIBMESH_HAVE_METAPHYSICL 
  231 #endif // #if PETSC_VERSION 
  232 #endif // #ifdef LIBMESH_HAVE_PETSC 
  234 #endif // LIBMESH_PETSC_DM_WRAPPER_H 
  
Manages consistently variables, degrees of freedom, and coefficient vectors.
 
std::vector< PetscInt > subfields
Stores subfield ids for use in subprojection matrixes on coarser DMs.
 
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.
 
PetscSection & get_section(unsigned int level)
Get reference to PetscSection for the given mesh level.
 
std::vector< std::unique_ptr< PetscDMContext > > _ctx_vec
Vector of internal PetscDM context structs for all grid levels.
 
DM & get_dm(unsigned int level)
Get reference to DM for the given mesh level.
 
void add_dofs_helper(const System &system, const DofObject &dof_object, dof_id_type local_id, PetscSection §ion)
Helper function to reduce code duplication when setting dofs in section.
 
The libMesh namespace provides an interface to certain functionality in the library.
 
void build_section(const System &system, PetscSection §ion)
Takes System, empty PetscSection and populates the PetscSection.
 
std::vector< std::unique_ptr< DM > > _dms
Vector of DMs for all grid levels.
 
PetscVector< libMesh::Number > * current_vec
 
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.
 
Struct to house data regarding where in the mesh hierarchy we are located.
 
std::vector< std::vector< numeric_index_type > > dof_vec
Stores local dofs for each var for use in subprojection matrixes.
 
std::vector< unsigned int > _mesh_dof_loc_sizes
Stores n_local_dofs for each grid level, to be used for projection vector sizing.
 
The DofObject defines an abstract base class for objects that have degrees of freedom associated with...
 
std::vector< std::unique_ptr< PetscMatrix< Number > > > _subpmtx_vec
Vector of sub projection matrixes for all grid levels for fieldsplit.
 
std::vector< std::unique_ptr< PetscSection > > _sections
Vector of PETScSections for all grid levels.
 
PetscMatrix< libMesh::Number > * K_sub_interp_ptr
 
PetscMatrix< libMesh::Number > * K_restrict_ptr
 
std::vector< std::unique_ptr< PetscSF > > _star_forests
Vector of star forests for all grid levels.
 
This class defines a wrapper around the PETSc DM infrastructure.
 
std::vector< unsigned int > _mesh_dof_sizes
Stores n_dofs for each grid level, to be used for projection matrix sizing.
 
void build_sf(const System &system, PetscSF &star_forest)
Takes System, empty PetscSF and populates the PetscSF.
 
PetscMatrix< libMesh::Number > * K_interp_ptr
 
void init_dm_data(unsigned int n_levels, const Parallel::Communicator &comm)
Init all the n_mesh_level dependent data structures.
 
void init_and_attach_petscdm(System &system, SNES &snes)
 
void clear()
Destroys and clears all build DM-related data.
 
std::vector< std::unique_ptr< PetscMatrix< Number > > > _pmtx_vec
Vector of projection matrixes for all grid levels.
 
dof_id_type check_section_n_dofs(PetscSection §ion)
Helper function to sanity check PetscSection construction.
 
PetscSF & get_star_forest(unsigned int level)
Get reference to PetscSF for the given mesh level.
 
std::vector< std::unique_ptr< PetscVector< Number > > > _vec_vec
Vector of solution vectors for all grid levels.