|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_dof_handler_h
17 #define dealii_dof_handler_h
41 #include <boost/serialization/split_member.hpp>
52 template <
int dim,
int spacedim>
54 template <
int dim,
int spacedim>
59 namespace DoFHandlerImplementation
61 struct Implementation;
65 template <
int dim,
int spacedim>
67 struct Implementation;
71 namespace DoFAccessorImplementation
73 struct Implementation;
76 namespace DoFCellAccessorImplementation
78 struct Implementation;
204 template <
int dim,
int spacedim = dim>
208 Iterators<DoFHandler<dim, spacedim>,
false>;
210 Iterators<DoFHandler<dim, spacedim>,
true>;
622 renumber_dofs(
const std::vector<types::global_dof_index> &new_numbers);
630 const std::vector<types::global_dof_index> &new_numbers);
725 end(
const unsigned int level)
const;
948 template <
typename number>
952 &boundary_ids)
const;
959 n_boundary_dofs(
const std::set<types::boundary_id> &boundary_ids)
const;
1074 get_fe(
const unsigned int index = 0)
const;
1106 template <
class Archive>
1108 save(Archive &ar,
const unsigned int version)
const;
1114 template <
class Archive>
1116 load(Archive &ar,
const unsigned int version);
1123 template <
class Archive>
1125 serialize(Archive &archive,
const unsigned int version);
1129 BOOST_SERIALIZATION_SPLIT_MEMBER()
1148 <<
"The given list of new dof indices is not consecutive: "
1149 <<
"the index " << arg1 <<
" does not exist.");
1156 <<
"The given level " << arg1
1157 <<
" is not in the valid range!");
1169 <<
"You tried to do something on level " << arg1
1170 <<
", but this level is empty.");
1196 std::unique_ptr<::internal::DoFHandlerImplementation::Policy::
1197 PolicyBase<dim, spacedim>>
1212 std::vector<::internal::DoFHandlerImplementation::NumberCache>
1237 const unsigned int dofs_per_vertex);
1257 const unsigned int dof_number,
1258 const unsigned int dofs_per_vertex)
const;
1266 const unsigned int dof_number,
1267 const unsigned int dofs_per_vertex,
1290 std::unique_ptr<types::global_dof_index[]>
indices;
1305 template <
int structdim>
1308 const unsigned int obj_index,
1309 const unsigned int fe_index,
1310 const unsigned int local_index)
const;
1312 template <
int structdim>
1315 const unsigned int obj_index,
1316 const unsigned int fe_index,
1317 const unsigned int local_index,
1336 std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel<dim>>>
1340 std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel<dim>>>
1348 std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces<dim>>
1351 std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces<dim>>
1355 template <
int,
class,
bool>
1357 template <
class,
bool>
1359 friend struct ::internal::DoFAccessorImplementation::Implementation;
1360 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1362 friend struct ::internal::DoFHandlerImplementation::Implementation;
1363 friend struct ::internal::DoFHandlerImplementation::Policy::
1368 #ifndef DEAL_II_MSVC
1370 "The dimension <dim> of a DoFHandler must be less than or "
1371 "equal to the space dimension <spacedim> in which it lives.");
1383 template <
int dim,
int spacedim>
1387 return mg_number_cache.size() > 0;
1392 template <
int dim,
int spacedim>
1396 return number_cache.n_global_dofs > 0;
1401 template <
int dim,
int spacedim>
1405 return number_cache.n_global_dofs;
1410 template <
int dim,
int spacedim>
1416 "n_dofs(level) can only be called after distribute_mg_dofs()"));
1418 return mg_number_cache[
level].n_global_dofs;
1423 template <
int dim,
int spacedim>
1427 return number_cache.n_locally_owned_dofs;
1432 template <
int dim,
int spacedim>
1436 return number_cache.locally_owned_dofs;
1441 template <
int dim,
int spacedim>
1445 Assert(level < this->get_triangulation().n_global_levels(),
1446 ExcMessage(
"The given level index exceeds the number of levels "
1447 "present in the triangulation"));
1449 mg_number_cache.size() == this->get_triangulation().n_global_levels(),
1451 "The level dofs are not set up properly! Did you call distribute_mg_dofs()?"));
1452 return mg_number_cache[
level].locally_owned_dofs;
1457 template <
int dim,
int spacedim>
1458 const std::vector<types::global_dof_index> &
1461 if (number_cache.n_locally_owned_dofs_per_processor.empty() &&
1462 number_cache.n_global_dofs > 0)
1468 &this->get_triangulation()));
1472 comm = MPI_COMM_SELF;
1476 .n_locally_owned_dofs_per_processor =
1479 return number_cache.n_locally_owned_dofs_per_processor;
1484 template <
int dim,
int spacedim>
1485 const std::vector<IndexSet> &
1488 if (number_cache.locally_owned_dofs_per_processor.empty() &&
1489 number_cache.n_global_dofs > 0)
1495 &this->get_triangulation()));
1499 comm = MPI_COMM_SELF;
1503 .locally_owned_dofs_per_processor =
1506 return number_cache.locally_owned_dofs_per_processor;
1511 template <
int dim,
int spacedim>
1512 const std::vector<IndexSet> &
1514 const unsigned int level)
const
1516 Assert(level < this->get_triangulation().n_global_levels(),
1517 ExcMessage(
"The given level index exceeds the number of levels "
1518 "present in the triangulation"));
1520 mg_number_cache.size() == this->get_triangulation().n_global_levels(),
1522 "The level dofs are not set up properly! Did you call distribute_mg_dofs()?"));
1523 if (mg_number_cache[
level].locally_owned_dofs_per_processor.empty() &&
1524 mg_number_cache[
level].n_global_dofs > 0)
1530 &this->get_triangulation()));
1534 comm = MPI_COMM_SELF;
1537 mg_number_cache[
level])
1538 .locally_owned_dofs_per_processor =
1539 mg_number_cache[
level].get_locally_owned_dofs_per_processor(comm);
1546 template <
int dim,
int spacedim>
1553 "There is only one FiniteElement stored. The index must be zero!"));
1554 return get_fe_collection()[0];
1559 template <
int dim,
int spacedim>
1564 fe_collection.size() > 0,
1566 "You are trying to access the DoFHandler's FECollection object before it has been initialized."));
1567 return fe_collection;
1572 template <
int dim,
int spacedim>
1577 ExcMessage(
"This DoFHandler object has not been associated "
1578 "with a triangulation."));
1584 template <
int dim,
int spacedim>
1588 return block_info_object;
1593 template <
int dim,
int spacedim>
1594 template <
typename number>
1598 &boundary_ids)
const
1602 std::set<types::boundary_id> boundary_ids_only;
1605 boundary_ids.begin();
1606 p != boundary_ids.end();
1608 boundary_ids_only.insert(p->first);
1611 return n_boundary_dofs(boundary_ids_only);
1625 template <
int dim,
int spacedim>
1628 PolicyBase<dim, spacedim> &policy);
1633 template <
int dim,
int spacedim>
1634 template <
class Archive>
1638 ar &block_info_object;
1645 unsigned int n_levels = levels.size();
1647 for (
unsigned int i = 0; i < levels.size(); ++i)
1653 bool faces_is_nullptr = (faces.get() ==
nullptr);
1654 ar & faces_is_nullptr;
1655 if (!faces_is_nullptr)
1661 unsigned int n_cells = tria->n_cells();
1662 std::string fe_name = this->get_fe(0).get_name();
1665 ar &
n_cells &fe_name &policy_name;
1670 template <
int dim,
int spacedim>
1671 template <
class Archive>
1675 ar &block_info_object;
1691 levels.resize(size);
1692 for (
unsigned int i = 0; i < levels.size(); ++i)
1694 std::unique_ptr<internal::DoFHandlerImplementation::DoFLevel<dim>>
level;
1696 levels[i] = std::move(
level);
1700 bool faces_is_nullptr =
true;
1701 ar & faces_is_nullptr;
1702 if (!faces_is_nullptr)
1708 std::string fe_name;
1709 std::string policy_name;
1711 ar &
n_cells &fe_name &policy_name;
1715 "The object being loaded into does not match the triangulation "
1716 "that has been stored previously."));
1718 fe_name == this->get_fe(0).get_name(),
1720 "The finite element associated with this DoFHandler does not match "
1721 "the one that was associated with the DoFHandler previously stored."));
1724 "The policy currently associated with this DoFHandler (" +
1726 ") does not match the one that was associated with the "
1727 "DoFHandler previously stored (" +
1728 policy_name +
")."));
1733 template <
int dim,
int spacedim>
1736 const unsigned int level,
1737 const unsigned int dof_number,
1738 const unsigned int dofs_per_vertex)
const
1747 template <
int dim,
int spacedim>
1750 const unsigned int level,
1751 const unsigned int dof_number,
1752 const unsigned int dofs_per_vertex,
1757 indices[dofs_per_vertex * (
level - coarsest_level) + dof_number] = index;
std::vector<::internal::DoFHandlerImplementation::NumberCache > mg_number_cache
BlockInfo block_info_object
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
static const unsigned int default_fe_index
std::unique_ptr< types::global_dof_index[]> indices
IteratorRange< level_cell_iterator > mg_cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > faces
typename ActiveSelector::cell_iterator cell_iterator
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
typename ActiveSelector::active_cell_iterator active_cell_iterator
typename ActiveSelector::CellAccessor cell_accessor
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > levels
static ::ExceptionBase & ExcFacesHaveNoLevel()
typename LevelSelector::FaceAccessor level_face_accessor
level_cell_iterator begin_mg(const unsigned int level=0) const
typename LevelSelector::cell_iterator level_cell_iterator
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
const hp::FECollection< dim, spacedim > & get_fe_collection() const
typename ActiveSelector::face_iterator face_iterator
active_cell_iterator end_active(const unsigned int level) const
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
unsigned int get_finest_level() const
typename ActiveSelector::FaceAccessor face_accessor
static ::ExceptionBase & ExcInvalidLevel(int arg1)
void initialize_local_block_info()
static const unsigned int dim
const std::vector< IndexSet > & locally_owned_dofs_per_processor() const
static const unsigned int spacedim
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
unsigned int get_coarsest_level() const
unsigned int max_couplings_between_boundary_dofs() const
types::global_dof_index get_index(const unsigned int level, const unsigned int dof_number, const unsigned int dofs_per_vertex) const
virtual void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
static ::ExceptionBase & ExcNewNumbersNotConsecutive(types::global_dof_index arg1)
const IndexSet & locally_owned_dofs() const
static ::ExceptionBase & ExcEmptyLevel(int arg1)
IteratorRange< active_cell_iterator > active_cell_iterators() const
cell_iterator begin(const unsigned int level=0) const
IteratorRange< level_cell_iterator > mg_cell_iterators() const
virtual ~DoFHandler() override
types::global_dof_index get_dof_index(const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index) const
void load(Archive &ar, const unsigned int version)
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
DoFHandler & operator=(const DoFHandler &)=delete
bool has_level_dofs() const
std::vector< MGVertexDoFs > mg_vertex_dofs
const Triangulation< dim, spacedim > & get_triangulation() const
static const unsigned int space_dimension
#define DEAL_II_NAMESPACE_OPEN
void save(Archive &ar, const unsigned int version) const
typename ActiveSelector::active_line_iterator active_line_iterator
typename ActiveSelector::active_face_iterator active_face_iterator
std::vector< types::global_dof_index > vertex_dofs
const BlockInfo & block_info() const
#define DEAL_II_DEPRECATED
#define DeclException1(Exception1, type1, outsequence)
static const unsigned int dimension
level_cell_iterator end_mg() const
void set_index(const unsigned int level, const unsigned int dof_number, const unsigned int dofs_per_vertex, const types::global_dof_index index)
IteratorRange< cell_iterator > cell_iterators() const
std::vector< types::global_dof_index > get_n_locally_owned_dofs_per_processor(const MPI_Comm mpi_communicator) const
cell_iterator end() const
unsigned int max_couplings_between_dofs() const
typename ActiveSelector::hex_iterator hex_iterator
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcGridsDoNotMatch()
unsigned int coarsest_level
void init(const unsigned int coarsest_level, const unsigned int finest_level, const unsigned int dofs_per_vertex)
virtual const MPI_Comm & get_communicator() const
#define Assert(cond, exc)
typename LevelSelector::face_iterator level_face_iterator
static const bool is_hp_dof_handler
virtual std::size_t memory_consumption() const
::internal::DoFHandlerImplementation::NumberCache number_cache
#define DeclException0(Exception0)
typename ActiveSelector::active_hex_iterator active_hex_iterator
std::string policy_to_string(const ::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > &policy)
typename ActiveSelector::active_quad_iterator active_quad_iterator
bool has_active_dofs() const
void set_dof_index(const unsigned int obj_level, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index) const
virtual void distribute_mg_dofs()
types::global_dof_index n_locally_owned_dofs() const
std::vector< IndexSet > get_locally_owned_dofs_per_processor(const MPI_Comm mpi_communicator) const
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
const std::vector< types::global_dof_index > & n_locally_owned_dofs_per_processor() const
const std::vector< IndexSet > & locally_owned_mg_dofs_per_processor(const unsigned int level) const
types::global_dof_index n_boundary_dofs() const
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
#define DEAL_II_NAMESPACE_CLOSE
unsigned int finest_level
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
std::vector< IndexSet > locally_owned_dofs_per_processor
typename ActiveSelector::quad_iterator quad_iterator
active_cell_iterator begin_active(const unsigned int level=0) const
typename ActiveSelector::line_iterator line_iterator
hp::FECollection< dim, spacedim > fe_collection
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
void initialize(const Triangulation< dim, spacedim > &tria, const FiniteElement< dim, spacedim > &fe)
void serialize(Archive &archive, const unsigned int version)
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
#define AssertThrow(cond, exc)
virtual void set_fe(const FiniteElement< dim, spacedim > &fe)
typename LevelSelector::CellAccessor level_cell_accessor
types::global_dof_index n_dofs() const