16#ifndef dealii_dof_handler_h
17#define dealii_dof_handler_h
40#include <boost/serialization/split_member.hpp>
51template <
int dim,
int spacedim>
53template <
int dim,
int spacedim>
58 namespace DoFHandlerImplementation
60 struct Implementation;
64 template <
int dim,
int spacedim>
66 struct Implementation;
70 namespace DoFAccessorImplementation
72 struct Implementation;
75 namespace DoFCellAccessorImplementation
77 struct Implementation;
82 namespace DoFHandlerImplementation
84 struct Implementation;
93 template <
int dim,
int spacedim,
typename VectorType>
313template <
int dim,
int spacedim = dim>
811 const std::vector<types::global_dof_index> &new_numbers);
1127 template <
typename number>
1131 &boundary_ids)
const;
1326 template <
class Archive>
1328 save(Archive &ar,
const unsigned int version)
const;
1335 template <
class Archive>
1337 load(Archive &ar,
const unsigned int version);
1345 template <
class Archive>
1351 BOOST_SERIALIZATION_SPLIT_MEMBER()
1369 <<
"The given level " << arg1
1370 <<
" is not in the valid range!");
1377 <<
"The given list of new dof indices is not consecutive: "
1378 <<
"the index " << arg1 <<
" does not exist.");
1385 <<
"The mesh contains a cell with an active FE index of "
1386 << arg1 <<
", but the finite element collection only has "
1387 << arg2 <<
" elements");
1394 "The current function doesn't make sense when used with a "
1395 "DoFHandler without hp-capabilities.");
1402 "The current function has not yet been implemented for a "
1403 "DoFHandler with hp-capabilities.");
1428 const unsigned int dofs_per_vertex);
1448 const unsigned int dof_number,
1449 const unsigned int dofs_per_vertex)
const;
1457 const unsigned int dof_number,
1458 const unsigned int dofs_per_vertex,
1481 std::unique_ptr<types::global_dof_index[]>
indices;
1523 parallel::distributed::
1524 CellDataTransfer<dim, spacedim, std::vector<unsigned int>>>
1556 std::unique_ptr<::internal::DoFHandlerImplementation::Policy::
1557 PolicyBase<dim, spacedim>>
1572 std::vector<::internal::DoFHandlerImplementation::NumberCache>
1580 mutable std::vector<std::vector<types::global_dof_index>>
1594 mutable std::vector<std::array<std::vector<types::global_dof_index>, dim + 1>>
1605 mutable std::vector<std::array<std::vector<offset_type>, dim + 1>>
1613 mutable std::array<std::vector<active_fe_index_type>, dim + 1>
1625 mutable std::vector<std::vector<active_fe_index_type>>
1632 mutable std::vector<std::vector<active_fe_index_type>>
1645 std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel<dim>>>
1651 std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces<dim>>
1688 template <
int structdim>
1691 const unsigned int obj_index,
1692 const unsigned int fe_index,
1693 const unsigned int local_index)
const;
1698 template <
int structdim>
1701 const unsigned int obj_index,
1702 const unsigned int fe_index,
1703 const unsigned int local_index,
1794 template <
int,
int,
int,
bool>
1795 friend class ::DoFAccessor;
1796 template <
int,
int,
bool>
1797 friend class ::DoFCellAccessor;
1798 friend struct ::internal::DoFAccessorImplementation::Implementation;
1799 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
1803 friend struct ::internal::DoFHandlerImplementation::Implementation;
1804 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
1805 friend struct ::internal::DoFHandlerImplementation::Policy::
1811 static_assert(dim <= spacedim,
1812 "The dimension <dim> of a DoFHandler must be less than or "
1813 "equal to the space dimension <spacedim> in which it lives.");
1821 namespace DoFHandlerImplementation
1834 template <
int dim,
int spacedim>
1862 template <
int dim,
int spacedim = dim>
1872 "No FiniteElement has been found in your FECollection that is "
1873 "dominated by all children of a cell you are trying to coarsen!");
1884template <
int dim,
int spacedim>
1893template <
int dim,
int spacedim>
1902template <
int dim,
int spacedim>
1911template <
int dim,
int spacedim>
1920template <
int dim,
int spacedim>
1926 "n_dofs(level) can only be called after distribute_mg_dofs()"));
1933template <
int dim,
int spacedim>
1942template <
int dim,
int spacedim>
1951template <
int dim,
int spacedim>
1956 ExcMessage(
"The given level index exceeds the number of levels "
1957 "present in the triangulation"));
1959 mg_number_cache.size() == this->get_triangulation().n_global_levels(),
1961 "The level dofs are not set up properly! Did you call distribute_mg_dofs()?"));
1967template <
int dim,
int spacedim>
1968const std::vector<types::global_dof_index> &
1984template <
int dim,
int spacedim>
1985const std::vector<IndexSet> &
2001template <
int dim,
int spacedim>
2002const std::vector<IndexSet> &
2004 const unsigned int level)
const
2007 ExcMessage(
"The given level index exceeds the number of levels "
2008 "present in the triangulation"));
2010 mg_number_cache.size() == this->get_triangulation().n_global_levels(),
2012 "The level dofs are not set up properly! Did you call distribute_mg_dofs()?"));
2027template <
int dim,
int spacedim>
2032 ExcMessage(
"No finite element collection is associated with "
2033 "this DoFHandler"));
2039template <
int dim,
int spacedim>
2048template <
int dim,
int spacedim>
2053 ExcMessage(
"This DoFHandler object has not been associated "
2054 "with a triangulation."));
2060template <
int dim,
int spacedim>
2065 ExcMessage(
"This DoFHandler object has not been associated "
2066 "with a triangulation."));
2072template <
int dim,
int spacedim>
2083template <
int dim,
int spacedim>
2084template <
typename number>
2088 &boundary_ids)
const
2095 std::set<types::boundary_id> boundary_ids_only;
2098 boundary_ids.begin();
2099 p != boundary_ids.end();
2101 boundary_ids_only.insert(p->first);
2118 template <
int dim,
int spacedim>
2121 PolicyBase<dim, spacedim> &policy);
2126template <
int dim,
int spacedim>
2127template <
class Archive>
2153 std::string policy_name =
2176 ar &
n_cells &fe_name &policy_name;
2182template <
int dim,
int spacedim>
2183template <
class Archive>
2208 std::string policy_name;
2215 "The object being loaded into does not match the triangulation "
2216 "that has been stored previously."));
2219 ExcMessage(
"The policy currently associated with this DoFHandler (" +
2221 ") does not match the one that was associated with the "
2222 "DoFHandler previously stored (" +
2223 policy_name +
")."));
2243 std::string fe_name;
2244 std::string policy_name;
2246 ar &
n_cells &fe_name &policy_name;
2251 "The object being loaded into does not match the triangulation "
2252 "that has been stored previously."));
2254 fe_name == this->
get_fe(0).get_name(),
2256 "The finite element associated with this DoFHandler does not match "
2257 "the one that was associated with the DoFHandler previously stored."));
2260 "The policy currently associated with this DoFHandler (" +
2262 ") does not match the one that was associated with the "
2263 "DoFHandler previously stored (" +
2264 policy_name +
")."));
2270template <
int dim,
int spacedim>
2273 const unsigned int level,
2274 const unsigned int dof_number,
2275 const unsigned int dofs_per_vertex)
const
2278 ExcInvalidLevel(
level));
2279 return indices[dofs_per_vertex * (
level - coarsest_level) + dof_number];
2284template <
int dim,
int spacedim>
2287 const unsigned int level,
2288 const unsigned int dof_number,
2289 const unsigned int dofs_per_vertex,
2293 ExcInvalidLevel(
level));
2294 indices[dofs_per_vertex * (
level - coarsest_level) + dof_number] = index;
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
unsigned int coarsest_level
void init(const unsigned int coarsest_level, const unsigned int finest_level, const unsigned int dofs_per_vertex)
unsigned int get_finest_level() const
std::unique_ptr< types::global_dof_index[]> indices
unsigned int finest_level
unsigned int get_coarsest_level() const
types::global_dof_index get_index(const unsigned int level, const unsigned int dof_number, const unsigned int dofs_per_vertex) 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)
static const unsigned int dimension
cell_iterator end() const
void pre_distributed_transfer_action()
void post_transfer_action()
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 std::size_t memory_consumption() const
std::vector< std::vector< active_fe_index_type > > hp_cell_active_fe_indices
const std::vector< IndexSet > & locally_owned_mg_dofs_per_processor(const unsigned int level) const
unsigned int max_couplings_between_dofs() const
static const unsigned int space_dimension
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
std::unique_ptr< ActiveFEIndexTransfer > active_fe_index_transfer
hp::FECollection< dim, spacedim > fe_collection
void create_active_fe_table()
SmartPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
void initialize(const Triangulation< dim, spacedim > &tria, const hp::FECollection< dim, spacedim > &fe)
const hp::FECollection< dim, spacedim > & get_fe_collection() const
void save(Archive &ar, const unsigned int version) const
void set_active_fe_indices(const std::vector< unsigned int > &active_fe_indices)
BlockInfo block_info_object
bool has_active_dofs() const
std::vector< boost::signals2::connection > tria_listeners_for_transfer
level_cell_iterator end_mg() const
void renumber_dofs(const std::vector< types::global_dof_index > &new_numbers)
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
void renumber_dofs(const unsigned int level, const std::vector< types::global_dof_index > &new_numbers)
level_cell_iterator begin_mg(const unsigned int level=0) const
void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
types::global_dof_index n_boundary_dofs() const
DoFHandler(const Triangulation< dim, spacedim > &tria)
void connect_to_triangulation_signals()
const std::vector< types::global_dof_index > & n_locally_owned_dofs_per_processor() const
std::vector< std::vector< types::global_dof_index > > cell_dof_cache_indices
std::vector< MGVertexDoFs > mg_vertex_dofs
void set_fe(const hp::FECollection< dim, spacedim > &fe)
std::vector< std::vector< active_fe_index_type > > hp_cell_future_fe_indices
const std::vector< IndexSet > & locally_owned_dofs_per_processor() const
std::vector< boost::signals2::connection > tria_listeners
void post_distributed_transfer_action()
std::vector< std::array< std::vector< types::global_dof_index >, dim+1 > > object_dof_indices
const Triangulation< dim, spacedim > & get_triangulation() const
typename LevelSelector::CellAccessor level_cell_accessor
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 pre_transfer_action()
std::array< std::vector< active_fe_index_type >, dim+1 > hp_object_fe_indices
const IndexSet & locally_owned_dofs() const
void reinit(const Triangulation< dim, spacedim > &tria)
std::vector< std::array< std::vector< offset_type >, dim+1 > > object_dof_ptr
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
void serialize(Archive &archive, const unsigned int version)
DoFHandler & operator=(const DoFHandler &)=delete
DoFHandler(const DoFHandler &)=delete
active_cell_iterator begin_active(const unsigned int level=0) const
void distribute_mg_dofs()
active_cell_iterator end_active(const unsigned int level) const
void set_fe(const FiniteElement< dim, spacedim > &fe)
bool hp_capability_enabled
bool has_hp_capabilities() const
void initialize_local_block_info()
types::global_dof_index n_dofs() const
void update_active_fe_table()
level_cell_iterator end_mg(const unsigned int level) const
std::vector< std::vector< offset_type > > cell_dof_cache_ptr
cell_iterator end(const unsigned int level) const
typename LevelSelector::cell_iterator level_cell_iterator
void prepare_for_serialization_of_active_fe_indices()
std::vector<::internal::DoFHandlerImplementation::NumberCache > mg_number_cache
void load(Archive &ar, const unsigned int version)
static const unsigned int default_fe_index
typename LevelSelector::face_iterator level_face_iterator
cell_iterator begin(const unsigned int level=0) const
std::array< std::vector< offset_type >, dim+1 > hp_object_fe_ptr
MPI_Comm get_communicator() const
static const unsigned int invalid_fe_index
void get_active_fe_indices(std::vector< unsigned int > &active_fe_indices) const
void distribute_dofs(const hp::FECollection< dim, spacedim > &fe)
types::global_dof_index n_boundary_dofs(const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_ids) const
bool has_level_dofs() const
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
typename LevelSelector::FaceAccessor level_face_accessor
types::global_dof_index n_dofs(const unsigned int level) const
unsigned int max_couplings_between_boundary_dofs() const
void initialize(const Triangulation< dim, spacedim > &tria, const FiniteElement< dim, spacedim > &fe)
types::global_dof_index n_locally_owned_dofs() const
static const active_fe_index_type invalid_active_fe_index
unsigned short int active_fe_index_type
void deserialize_active_fe_indices()
types::global_dof_index n_boundary_dofs(const std::set< types::boundary_id > &boundary_ids) const
const BlockInfo & block_info() const
::internal::DoFHandlerImplementation::NumberCache number_cache
virtual ~DoFHandler() override
virtual std::string get_name() const =0
virtual MPI_Comm get_communicator() const
unsigned int n_cells() const
unsigned int size() const
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< level_cell_iterator > mg_cell_iterators() const
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
IteratorRange< cell_iterator > cell_iterators() const
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
IteratorRange< level_cell_iterator > mg_cell_iterators_on_level(const unsigned int level) const
#define DeclException0(Exception0)
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
static ::ExceptionBase & ExcOnlyAvailableWithHP()
#define Assert(cond, exc)
static ::ExceptionBase & ExcNotImplementedWithHP()
static ::ExceptionBase & ExcNewNumbersNotConsecutive(types::global_dof_index arg1)
static ::ExceptionBase & ExcInvalidFEIndex(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcNoDominatedFiniteElementOnChildren()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcInvalidLevel(int arg1)
static ::ExceptionBase & ExcNoFESelected()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::FaceAccessor face_accessor
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::line_iterator line_iterator
typename ActiveSelector::face_iterator face_iterator
typename ActiveSelector::active_line_iterator active_line_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
typename ActiveSelector::CellAccessor cell_accessor
typename ActiveSelector::quad_iterator quad_iterator
typename ActiveSelector::active_hex_iterator active_hex_iterator
typename ActiveSelector::active_quad_iterator active_quad_iterator
typename ActiveSelector::hex_iterator hex_iterator
typename ActiveSelector::active_face_iterator active_face_iterator
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
unsigned int dominated_future_fe_on_children(const typename DoFHandler< dim, spacedim >::cell_iterator &parent)
void communicate_future_fe_indices(DoFHandler< dim, spacedim > &dof_handler)
std::string policy_to_string(const ::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > &policy)
static const unsigned int invalid_unsigned_int
std::map< const cell_iterator, const unsigned int > coarsened_cells_fe_index
std::vector< unsigned int > active_fe_indices
std::map< const cell_iterator, const unsigned int > refined_cells_fe_index
std::unique_ptr< parallel::distributed::CellDataTransfer< dim, spacedim, std::vector< unsigned int > > > cell_data_transfer
std::map< const cell_iterator, const unsigned int > persisting_cells_fe_index
std::vector< types::global_dof_index > n_locally_owned_dofs_per_processor
std::vector< IndexSet > locally_owned_dofs_per_processor
types::global_dof_index n_global_dofs
std::vector< IndexSet > get_locally_owned_dofs_per_processor(const MPI_Comm &mpi_communicator) const
IndexSet locally_owned_dofs
std::vector< types::global_dof_index > get_n_locally_owned_dofs_per_processor(const MPI_Comm &mpi_communicator) const
types::global_dof_index n_locally_owned_dofs