16 #ifndef dealii_hp_dof_handler_h 17 #define dealii_hp_dof_handler_h 21 #include <deal.II/base/config.h> 23 #include <deal.II/base/exceptions.h> 24 #include <deal.II/base/function.h> 25 #include <deal.II/base/iterator_range.h> 26 #include <deal.II/base/smartpointer.h> 27 #include <deal.II/base/template_constraints.h> 29 #include <deal.II/distributed/cell_data_transfer.templates.h> 31 #include <deal.II/dofs/deprecated_function_map.h> 32 #include <deal.II/dofs/dof_accessor.h> 33 #include <deal.II/dofs/dof_iterator_selector.h> 34 #include <deal.II/dofs/number_cache.h> 36 #include <deal.II/hp/dof_faces.h> 37 #include <deal.II/hp/dof_level.h> 38 #include <deal.II/hp/fe_collection.h> 45 DEAL_II_NAMESPACE_OPEN
47 template <
int dim,
int spacedim>
52 namespace DoFHandlerImplementation
54 struct Implementation;
58 template <
int dim,
int spacedim>
60 struct Implementation;
68 namespace DoFHandlerImplementation
70 struct Implementation;
77 namespace DoFAccessorImplementation
79 struct Implementation;
82 namespace DoFCellAccessorImplementation
84 struct Implementation;
193 template <
int dim,
int spacedim = dim>
196 using ActiveSelector = ::internal::DoFHandlerImplementation::
197 Iterators<DoFHandler<dim, spacedim>,
false>;
198 using LevelSelector = ::internal::DoFHandlerImplementation::
199 Iterators<DoFHandler<dim, spacedim>,
true>;
309 using level_cell_iterator =
typename LevelSelector::cell_iterator;
334 using level_cell_accessor =
typename LevelSelector::CellAccessor;
335 using level_face_accessor =
typename LevelSelector::FaceAccessor;
337 using level_face_iterator =
typename LevelSelector::face_iterator;
545 renumber_dofs(
const std::vector<types::global_dof_index> &new_numbers);
590 begin(
const unsigned int level = 0)
const;
623 end(
const unsigned int level)
const;
769 n_dofs(
const unsigned int level)
const;
787 template <
typename number>
791 &boundary_ids)
const;
798 n_boundary_dofs(
const std::set<types::boundary_id> &boundary_ids)
const;
841 const std::vector<IndexSet> &
860 const std::vector<types::global_dof_index> &
878 const std::vector<IndexSet> &
896 get_fe(
const unsigned int index)
const;
964 template <
class Archive>
966 save(Archive &ar,
const unsigned int version)
const;
972 template <
class Archive>
974 load(Archive &ar,
const unsigned int version);
976 BOOST_SERIALIZATION_SPLIT_MEMBER()
995 << "The matrix has the wrong
dimension " << arg1);
1004 types::global_dof_index,
1005 << "The given list of new dof indices is not consecutive: "
1006 << "the index " << arg1 << " does not exist.");
1013 << "The mesh contains a cell with an active_fe_index of "
1014 << arg1 << ", but the finite element collection only has "
1015 << arg2 << " elements");
1021 << "The given level " << arg1
1022 << " is not in the valid range!");
1032 << "You tried to do something on level " << arg1
1033 << ", but this level is empty.");
1051 std::unique_ptr<::
internal::DoFHandlerImplementation::Policy::
1052 PolicyBase<dim, spacedim>>
1068 template <
int structdim>
1069 types::global_dof_index
1070 get_dof_index(const
unsigned int obj_level,
1071 const
unsigned int obj_index,
1072 const
unsigned int fe_index,
1073 const
unsigned int local_index) const;
1075 template <
int structdim>
1077 set_dof_index(const
unsigned int obj_level,
1078 const
unsigned int obj_index,
1079 const
unsigned int fe_index,
1080 const
unsigned int local_index,
1081 const
types::global_dof_index global_index) const;
1196 std::vector<::
internal::DoFHandlerImplementation::NumberCache>
1239 std::map<const cell_iterator, const unsigned int>
1252 std::map<const cell_iterator, const unsigned int>
1268 parallel::distributed::
1269 CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>
1288 template <
int,
class,
bool>
1289 friend class ::DoFAccessor;
1290 template <
class,
bool>
1291 friend class ::DoFCellAccessor;
1292 friend struct ::internal::DoFAccessorImplementation::Implementation;
1293 friend struct ::internal::DoFCellAccessorImplementation::
1301 friend class ::internal::hp::DoFIndicesOnFacesOrEdges;
1302 friend struct ::internal::hp::DoFHandlerImplementation::
1304 friend struct ::internal::DoFHandlerImplementation::Policy::
1317 template <
int dim,
int spacedim>
1318 template <
typename number>
1322 &boundary_ids)
const 1326 std::set<types::boundary_id> boundary_ids_only;
1329 p = boundary_ids.begin();
1330 p != boundary_ids.end();
1332 boundary_ids_only.insert(p->first);
1335 return n_boundary_dofs(boundary_ids_only);
1349 template <
int dim,
int spacedim>
1351 policy_to_string(const ::internal::DoFHandlerImplementation::Policy::
1352 PolicyBase<dim, spacedim> &policy);
1358 template <
int dim,
int spacedim>
1359 template <
class Archive>
1364 ar &vertex_dof_offsets;
1366 ar &mg_number_cache;
1371 const unsigned int n_levels = levels.size();
1373 for (
unsigned int i = 0; i < n_levels; ++i)
1379 bool faces_is_nullptr = (faces.get() ==
nullptr);
1380 ar & faces_is_nullptr;
1381 if (!faces_is_nullptr)
1387 const unsigned int n_cells = tria->n_cells();
1388 std::string policy_name = ::internal::policy_to_string(*policy);
1390 ar &n_cells &policy_name;
1395 template <
int dim,
int spacedim>
1396 template <
class Archive>
1401 ar &vertex_dof_offsets;
1403 ar &mg_number_cache;
1417 levels.resize(size);
1418 for (
unsigned int i = 0; i < size; ++i)
1420 std::unique_ptr<::internal::hp::DoFLevel> level;
1422 levels[i] = std::move(level);
1426 bool faces_is_nullptr =
true;
1427 ar & faces_is_nullptr;
1428 if (!faces_is_nullptr)
1433 unsigned int n_cells;
1434 std::string policy_name;
1436 ar &n_cells &policy_name;
1439 n_cells == tria->n_cells(),
1441 "The object being loaded into does not match the triangulation " 1442 "that has been stored previously."));
1443 AssertThrow(policy_name == ::internal::policy_to_string(*policy),
1445 "The policy currently associated with this DoFHandler (" +
1446 ::internal::policy_to_string(*policy) +
1447 ") does not match the one that was associated with the " 1448 "DoFHandler previously stored (" +
1449 policy_name +
")."));
1452 template <
int dim,
int spacedim>
1456 return number_cache.n_global_dofs;
1461 template <
int dim,
int spacedim>
1471 template <
int dim,
int spacedim>
1475 return number_cache.n_locally_owned_dofs;
1480 template <
int dim,
int spacedim>
1484 return number_cache.locally_owned_dofs;
1489 template <
int dim,
int spacedim>
1490 const std::vector<types::global_dof_index> &
1493 return number_cache.n_locally_owned_dofs_per_processor;
1498 template <
int dim,
int spacedim>
1499 const std::vector<IndexSet> &
1502 return number_cache.locally_owned_dofs_per_processor;
1507 template <
int dim,
int spacedim>
1510 const unsigned int level)
const 1514 Assert(level < this->get_triangulation().n_global_levels(),
1515 ExcMessage(
"invalid level in locally_owned_mg_dofs"));
1516 return mg_number_cache[0].locally_owned_dofs;
1520 template <
int dim,
int spacedim>
1521 const std::vector<IndexSet> &
1523 const unsigned int level)
const 1527 Assert(level < this->get_triangulation().n_global_levels(),
1528 ExcMessage(
"invalid level in locally_owned_mg_dofs_per_processor"));
1529 return mg_number_cache[0].locally_owned_dofs_per_processor;
1534 template <
int dim,
int spacedim>
1538 Assert(fe_collection.size() > 0,
1539 ExcMessage(
"No finite element collection is associated with " 1540 "this DoFHandler"));
1541 return fe_collection;
1546 template <
int dim,
int spacedim>
1550 Assert(fe_collection.size() > 0,
1551 ExcMessage(
"No finite element collection is associated with " 1552 "this DoFHandler"));
1553 return fe_collection[number];
1558 template <
int dim,
int spacedim>
1562 Assert(fe_collection.size() > 0,
1563 ExcMessage(
"No finite element collection is associated with " 1564 "this DoFHandler"));
1565 return fe_collection;
1570 template <
int dim,
int spacedim>
1575 ExcMessage(
"This DoFHandler object has not been associated " 1576 "with a triangulation."));
1584 DEAL_II_NAMESPACE_CLOSE
cell_iterator begin(const unsigned int level=0) const
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
void post_refinement_action()
void get_active_fe_indices(std::vector< unsigned int > &active_fe_indices) const
const std::vector< IndexSet > & locally_owned_mg_dofs_per_processor(const unsigned int level) const
types::global_dof_index n_locally_owned_dofs() const
virtual void distribute_dofs(const hp::FECollection< dim, spacedim > &fe)
static const unsigned int invalid_unsigned_int
typename ActiveSelector::face_iterator face_iterator
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcInvalidFEIndex(int arg1, int arg2)
std::vector< unsigned int > active_fe_indices
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
std::vector<::internal::DoFHandlerImplementation::NumberCache > mg_number_cache
virtual std::size_t memory_consumption() const
void set_active_fe_indices(const std::vector< unsigned int > &active_fe_indices)
::internal::DoFHandlerImplementation::NumberCache number_cache
static ::ExceptionBase & ExcGridsDoNotMatch()
std::vector< std::unique_ptr<::internal::hp::DoFLevel > > levels
const hp::FECollection< dim, spacedim > & get_fe() const
std::map< const cell_iterator, const unsigned int > refined_cells_fe_index
const std::vector< types::global_dof_index > & n_locally_owned_dofs_per_processor() const
void pre_refinement_action()
active_cell_iterator begin_active(const unsigned int level=0) const
std::unique_ptr<::internal::hp::DoFIndicesOnFaces< dim > > faces
static ::ExceptionBase & ExcFunctionNotUseful()
std::map< const cell_iterator, const unsigned int > coarsened_cells_fe_index
void setup_policy_and_listeners()
static ::ExceptionBase & ExcNewNumbersNotConsecutive(types::global_dof_index arg1)
typename ActiveSelector::CellAccessor cell_accessor
#define AssertThrow(cond, exc)
std::unique_ptr< ActiveFEIndexTransfer > active_fe_index_transfer
std::vector< boost::signals2::connection > tria_listeners
std::unique_ptr< parallel::distributed::CellDataTransfer< dim, spacedim, std::vector< unsigned int > > > cell_data_transfer
static const unsigned int dimension
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
static const unsigned int default_fe_index
void pre_active_fe_index_transfer()
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
void post_distributed_serialization_of_active_fe_indices()
typename ActiveSelector::active_cell_iterator active_cell_iterator
static ::ExceptionBase & ExcNoFESelected()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::active_quad_iterator active_quad_iterator
typename ActiveSelector::quad_iterator quad_iterator
typename ActiveSelector::cell_iterator cell_iterator
static const unsigned int space_dimension
#define DeclException1(Exception1, type1, outsequence)
typename ActiveSelector::active_line_iterator active_line_iterator
#define Assert(cond, exc)
static ::ExceptionBase & ExcFacesHaveNoLevel()
IteratorRange< active_cell_iterator > active_cell_iterators() const
void prepare_for_serialization_of_active_fe_indices()
const IndexSet & locally_owned_dofs() const
void save(Archive &ar, const unsigned int version) const
const std::vector< IndexSet > & locally_owned_dofs_per_processor() const
#define DeclException0(Exception0)
unsigned int max_couplings_between_boundary_dofs() const
virtual ~DoFHandler() override
void post_active_fe_index_transfer()
static ::ExceptionBase & ExcInvalidLevel(int arg1)
const hp::FECollection< dim, spacedim > & get_fe_collection() const
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
static const bool is_hp_dof_handler
types::global_dof_index n_boundary_dofs() const
active_cell_iterator end_active(const unsigned int level) const
unsigned int global_dof_index
void initialize(const Triangulation< dim, spacedim > &tria, const hp::FECollection< dim, spacedim > &fe)
unsigned int max_couplings_between_dofs() const
void load(Archive &ar, const unsigned int version)
typename ActiveSelector::active_face_iterator active_face_iterator
hp::FECollection< dim, spacedim > fe_collection
static ::ExceptionBase & ExcMatrixHasWrongSize(int arg1)
const Triangulation< dim, spacedim > & get_triangulation() const
void post_distributed_active_fe_index_transfer()
types::global_dof_index n_dofs() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcNotImplemented()
std::map< const cell_iterator, const unsigned int > persisting_cells_fe_index
const types::global_dof_index invalid_dof_index
void create_active_fe_table()
std::vector< types::global_dof_index > vertex_dofs
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
void pre_distributed_active_fe_index_transfer()
virtual void set_fe(const hp::FECollection< dim, spacedim > &fe)
cell_iterator end() const
IteratorRange< cell_iterator > cell_iterators() const
void deserialize_active_fe_indices()
std::vector< unsigned int > vertex_dof_offsets
static ::ExceptionBase & ExcEmptyLevel(int arg1)
typename ActiveSelector::FaceAccessor face_accessor
typename ActiveSelector::active_hex_iterator active_hex_iterator
DoFHandler & operator=(const DoFHandler &)=delete
typename ActiveSelector::line_iterator line_iterator
static const types::global_dof_index invalid_dof_index
typename ActiveSelector::hex_iterator hex_iterator