deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
DoFHandler< dim, spacedim > Class Template Reference

#include <deal.II/dofs/dof_handler.h>

Inheritance diagram for DoFHandler< dim, spacedim >:
Inheritance graph
[legend]

Classes

struct  ActiveFEIndexTransfer
 
class  MGVertexDoFs
 

Public Types

using cell_accessor = typename ActiveSelector::CellAccessor
 
using face_accessor = typename ActiveSelector::FaceAccessor
 
using line_iterator = typename ActiveSelector::line_iterator
 
using active_line_iterator = typename ActiveSelector::active_line_iterator
 
using quad_iterator = typename ActiveSelector::quad_iterator
 
using active_quad_iterator = typename ActiveSelector::active_quad_iterator
 
using hex_iterator = typename ActiveSelector::hex_iterator
 
using active_hex_iterator = typename ActiveSelector::active_hex_iterator
 
using active_cell_iterator = typename ActiveSelector::active_cell_iterator
 
using cell_iterator = typename ActiveSelector::cell_iterator
 
using face_iterator = typename ActiveSelector::face_iterator
 
using active_face_iterator = typename ActiveSelector::active_face_iterator
 
using level_cell_accessor = typename LevelSelector::CellAccessor
 
using level_face_accessor = typename LevelSelector::FaceAccessor
 
using level_cell_iterator = typename LevelSelector::cell_iterator
 
using level_face_iterator = typename LevelSelector::face_iterator
 
using offset_type = unsigned int
 

Public Member Functions

 DoFHandler ()
 
 DoFHandler (const Triangulation< dim, spacedim > &tria)
 
 DoFHandler (const DoFHandler &)=delete
 
virtual ~DoFHandler () override
 
DoFHandleroperator= (const DoFHandler &)=delete
 
void set_active_fe_indices (const std::vector< types::fe_index > &active_fe_indices)
 
void set_active_fe_indices (const std::vector< unsigned int > &active_fe_indices)
 
std::vector< types::fe_indexget_active_fe_indices () const
 
void get_active_fe_indices (std::vector< unsigned int > &active_fe_indices) const
 
void set_future_fe_indices (const std::vector< types::fe_index > &future_fe_indices)
 
std::vector< types::fe_indexget_future_fe_indices () const
 
void reinit (const Triangulation< dim, spacedim > &tria)
 
void distribute_dofs (const FiniteElement< dim, spacedim > &fe)
 
void distribute_dofs (const hp::FECollection< dim, spacedim > &fe)
 
void distribute_mg_dofs ()
 
bool has_hp_capabilities () const
 
bool has_level_dofs () const
 
bool has_active_dofs () const
 
void initialize_local_block_info ()
 
void clear ()
 
void renumber_dofs (const std::vector< types::global_dof_index > &new_numbers)
 
void renumber_dofs (const unsigned int level, const std::vector< types::global_dof_index > &new_numbers)
 
unsigned int max_couplings_between_dofs () const
 
unsigned int max_couplings_between_boundary_dofs () const
 
types::global_dof_index n_dofs () const
 
types::global_dof_index n_dofs (const unsigned int level) const
 
types::global_dof_index n_boundary_dofs () const
 
template<typename number >
types::global_dof_index n_boundary_dofs (const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_ids) const
 
types::global_dof_index n_boundary_dofs (const std::set< types::boundary_id > &boundary_ids) const
 
const BlockInfoblock_info () const
 
types::global_dof_index n_locally_owned_dofs () const
 
const IndexSetlocally_owned_dofs () const
 
const IndexSetlocally_owned_mg_dofs (const unsigned int level) const
 
const FiniteElement< dim, spacedim > & get_fe (const types::fe_index index=0) const
 
const hp::FECollection< dim, spacedim > & get_fe_collection () const
 
const Triangulation< dim, spacedim > & get_triangulation () const
 
MPI_Comm get_mpi_communicator () const
 
MPI_Comm get_communicator () const
 
void prepare_for_serialization_of_active_fe_indices ()
 
void deserialize_active_fe_indices ()
 
virtual std::size_t memory_consumption () const
 
template<class Archive >
void save (Archive &ar, const unsigned int version) const
 
template<class Archive >
void load (Archive &ar, const unsigned int version)
 
template<class Archive >
void serialize (Archive &archive, const unsigned int version)
 
Cell iterator functions
cell_iterator begin (const unsigned int level=0) const
 
active_cell_iterator begin_active (const unsigned int level=0) const
 
cell_iterator end () const
 
cell_iterator end (const unsigned int level) const
 
active_cell_iterator end_active (const unsigned int level) const
 
level_cell_iterator begin_mg (const unsigned int level=0) const
 
level_cell_iterator end_mg (const unsigned int level) const
 
level_cell_iterator end_mg () const
 
Cell iterator functions returning ranges of iterators
IteratorRange< cell_iteratorcell_iterators () const
 
IteratorRange< active_cell_iteratoractive_cell_iterators () const
 
IteratorRange< level_cell_iteratormg_cell_iterators () const
 
IteratorRange< cell_iteratorcell_iterators_on_level (const unsigned int level) const
 
IteratorRange< active_cell_iteratoractive_cell_iterators_on_level (const unsigned int level) const
 
IteratorRange< level_cell_iteratormg_cell_iterators_on_level (const unsigned int level) const
 
EnableObserverPointer functionality

Classes derived from EnableObserverPointer provide a facility to subscribe to this object. This is mostly used by the ObserverPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcNoFESelected ()
 
static ::ExceptionBaseExcInvalidBoundaryIndicator ()
 
static ::ExceptionBaseExcInvalidLevel (int arg1)
 
static ::ExceptionBaseExcNewNumbersNotConsecutive (types::global_dof_index arg1)
 
static ::ExceptionBaseExcInvalidFEIndex (int arg1, int arg2)
 
static ::ExceptionBaseExcOnlyAvailableWithHP ()
 
static ::ExceptionBaseExcNotImplementedWithHP ()
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 
static constexpr unsigned int space_dimension = spacedim
 
static const types::fe_index default_fe_index = 0
 

Private Types

using ActiveSelector = ::internal::DoFHandlerImplementation::Iterators< dim, spacedim, false >
 
using LevelSelector = ::internal::DoFHandlerImplementation::Iterators< dim, spacedim, true >
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void clear_space ()
 
void clear_mg_space ()
 
void setup_policy ()
 
void connect_to_triangulation_signals ()
 
void create_active_fe_table ()
 
void update_active_fe_table ()
 
void pre_transfer_action ()
 
void post_transfer_action ()
 
void pre_distributed_transfer_action ()
 
void post_distributed_transfer_action ()
 
void check_no_subscribers () const noexcept
 

Private Attributes

BlockInfo block_info_object
 
bool hp_capability_enabled
 
ObserverPointer< const Triangulation< dim, spacedim >, DoFHandler< dim, spacedim > > tria
 
hp::FECollection< dim, spacedim > fe_collection
 
std::unique_ptr<::internal::DoFHandlerImplementation::Policy::PolicyBase< dim, spacedim > > policy
 
::internal::DoFHandlerImplementation::NumberCache number_cache
 
std::vector<::internal::DoFHandlerImplementation::NumberCachemg_number_cache
 
std::vector< std::array< std::vector< types::global_dof_index >, dim+1 > > object_dof_indices
 
std::vector< std::array< std::vector< offset_type >, dim+1 > > object_dof_ptr
 
std::array< std::vector< types::fe_index >, dim+1 > hp_object_fe_indices
 
std::array< std::vector< offset_type >, dim+1 > hp_object_fe_ptr
 
std::vector< std::vector< types::fe_index > > hp_cell_active_fe_indices
 
std::vector< std::vector< types::fe_index > > hp_cell_future_fe_indices
 
std::vector< MGVertexDoFsmg_vertex_dofs
 
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
 
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
 
std::unique_ptr< ActiveFEIndexTransferactive_fe_index_transfer
 
std::vector< boost::signals2::connection > tria_listeners
 
std::vector< boost::signals2::connection > tria_listeners_for_transfer
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

template<int , int , int , bool >
class ::DoFAccessor
 
template<int , int , bool >
class ::DoFCellAccessor
 
struct ::internal::DoFAccessorImplementation::Implementation
 
struct ::internal::DoFCellAccessorImplementation::Implementation
 
struct ::internal::DoFHandlerImplementation::Implementation
 
struct ::internal::hp::DoFHandlerImplementation::Implementation
 
struct ::internal::DoFHandlerImplementation::Policy::Implementation
 

Detailed Description

template<int dim, int spacedim = dim>
class DoFHandler< dim, spacedim >

Given a triangulation and a description of a finite element, this class enumerates degrees of freedom on all vertices, edges, faces, and cells of the triangulation. As a result, it also provides a basis for a discrete space \(V_h\) whose elements are finite element functions defined on each cell by a FiniteElement object. This class satisfies the MeshType concept requirements.

It is first used in the step-2 tutorial program.

For each 0d, 1d, 2d, and 3d subobject, this class stores a list of the indices of degrees of freedom defined on this DoFHandler. These indices refer to the unconstrained degrees of freedom, i.e. constrained degrees of freedom are numbered in the same way as unconstrained ones, and are only later eliminated. This leads to the fact that indices in global vectors and matrices also refer to all degrees of freedom and some kind of condensation is needed to restrict the systems of equations to the unconstrained degrees of freedom only. The actual layout of storage of the indices is described in the internal::DoFHandlerImplementation::DoFLevel class documentation.

The class offers iterators to traverse all cells, in much the same way as the Triangulation class does. Using the begin() and end() functions (and companions, like begin_active()), one can obtain iterators to walk over cells, and query the degree of freedom structures as well as the triangulation data. These iterators are built on top of those of the Triangulation class, but offer the additional information on degrees of freedom functionality compared to pure triangulation iterators. The order in which dof iterators are presented by the ++ and -- operators is the same as that for the corresponding iterators traversing the triangulation on which this DoFHandler is constructed.

The spacedim parameter has to be used if one wants to solve problems on surfaces. If not specified, this parameter takes the default value =dim implying that we want to solve problems in a domain whose dimension equals the dimension of the space in which it is embedded.

Distribution of indices for degrees of freedom

The degrees of freedom (‘dofs’) are distributed on the given triangulation by the function distribute_dofs(). It gets passed a finite element object describing how many degrees of freedom are located on vertices, lines, etc. It traverses the triangulation cell by cell and numbers the dofs of that cell if not yet numbered. For non-multigrid algorithms, only active cells are considered. Active cells are defined to be those cells which have no children, i.e. they are the most refined ones.

Since the triangulation is traversed starting with the cells of the coarsest active level and going to more refined levels, the lowest numbers for dofs are given to the largest cells as well as their bounding lines and vertices, with the dofs of more refined cells getting higher numbers.

This numbering implies very large bandwidths of the resulting matrices and is thus vastly suboptimal for some solution algorithms. For this reason, the DoFRenumbering class offers several algorithms to reorder the dof numbering according. See there for a discussion of the implemented algorithms.

Interaction with distributed meshes

Upon construction, this class takes a reference to a triangulation object. In most cases, this will be a reference to an object of type Triangulation, i.e. the class that represents triangulations that entirely reside on a single processor. However, it can also be of type parallel::distributed::Triangulation (see, for example, step-32, step-40 and in particular the Parallel computing with multiple processors using distributed memory topic) in which case the DoFHandler object will proceed to only manage degrees of freedom on locally owned and ghost cells. This process is entirely transparent to the used.

User defined renumbering schemes

The DoFRenumbering class offers a number of renumbering schemes like the Cuthill-McKee scheme. Basically, the function sets up an array in which for each degree of freedom we store the new index this DoF should have after renumbering. Using this array, the renumber_dofs() function of the present class is called, which actually performs the change from old DoF indices to the ones given in the array. In some cases, however, a user may want to compute their own renumbering order; in this case, one can allocate an array with one element per degree of freedom and fill it with the number that the respective degree of freedom shall be assigned. This number may, for example, be obtained by sorting the support points of the degrees of freedom in downwind direction. Then call the renumber_dofs(vector<types::global_dof_index>) function with the array, which converts old into new degree of freedom indices.

Serializing (loading or storing) DoFHandler objects

Like many other classes in deal.II, the DoFHandler class can stream its contents to an archive using BOOST's serialization facilities. The data so stored can later be retrieved again from the archive to restore the contents of this object. This facility is frequently used to save the state of a program to disk for possible later resurrection, often in the context of checkpoint/restart strategies for long running computations or on computers that aren't very reliable (e.g. on very large clusters where individual nodes occasionally fail and then bring down an entire MPI job).

The model for doing so is similar for the DoFHandler class as it is for the Triangulation class (see the section in the general documentation of that class). In particular, the load() function does not exactly restore the same state as was stored previously using the save() function. Rather, the function assumes that you load data into a DoFHandler object that is already associated with a triangulation that has a content that matches the one that was used when the data was saved. Likewise, the load() function assumes that the current object is already associated with a finite element object that matches the one that was associated with it when data was saved; the latter can be achieved by calling DoFHandler::distribute_dofs() using the same kind of finite element before re-loading data from the serialization archive.

hp-adaptive finite element methods

Instead of only using one particular FiniteElement on all cells, this class also allows for an enumeration of degrees of freedom on different finite elements on every cells. To this end, one assigns an active_fe_index to every cell that indicates which element within a collection of finite elements (represented by an object of type hp::FECollection) is the one that lives on this cell. The class then enumerates the degree of freedom associated with these finite elements on each cell of a triangulation and, if possible, identifies degrees of freedom at the interfaces of cells if they match. If neighboring cells have degrees of freedom along the common interface that do not immediate match (for example, if you have \(Q_2\) and \(Q_3\) elements meeting at a common face), then one needs to compute constraints to ensure that the resulting finite element space on the mesh remains conforming.

The whole process of working with objects of this type is explained in step-27. Many of the algorithms this class implements are described in the hp-paper.

Active FE indices and their behavior under mesh refinement

The typical workflow for using this class is to create a mesh, assign an active FE index to every active cell, call DoFHandler::distribute_dofs(), and then assemble a linear system and solve a problem on this finite element space.

Active FE indices will be automatically transferred during mesh adaptation from the old to the new mesh. Future FE indices are meant to determine the active FE index after mesh adaptation, and are used to prepare data on the old mesh for the new one. If no future FE index is specified, the finite element prevails.

In particular, the following rules apply during adaptation:

Strategies for automatic hp-adaptation which will set future FE indices based on criteria are available in the hp::Refinement namespace.

Active FE indices and parallel meshes

When this class is used with either a parallel::shared::Triangulation or a parallel::distributed::Triangulation, you can only set active FE indices on cells that are locally owned, using a call such as cell->set_active_fe_index(...). On the other hand, setting the active FE index on ghost or artificial cells is not allowed.

Ghost cells do acquire the information what element is active on them, however: whenever you call DoFHandler::distribute_dofs(), all processors that participate in the parallel mesh exchange information in such a way that the active FE index on ghost cells equals the active FE index that was set on that processor that owned that particular ghost cell. Consequently, one can query the active_fe_index on ghost cells, just not set it by hand.

On artificial cells, no information is available about the active_fe_index used there. That's because we don't even know whether these cells exist at all, and even if they did, the current processor does not know anything specific about them. See the glossary entry on artificial cells for more information.

During refinement and coarsening, information about the active_fe_index of each cell will be automatically transferred.

However, using a parallel::distributed::Triangulation with a DoFHandler in hp-mode requires additional attention during serialization, since no information on active FE indices will be automatically transferred. This has to be done manually using the prepare_for_serialization_of_active_fe_indices() and deserialize_active_fe_indices() functions. The former has to be called before parallel::distributed::Triangulation::save() is invoked, and the latter needs to be run after parallel::distributed::Triangulation::load(). If further data will be attached to the triangulation via the parallel::distributed::CellDataTransfer, parallel::distributed::SolutionTransfer, or Particles::ParticleHandler classes, all corresponding preparation and deserialization function calls need to happen in the same order. Consult the documentation of parallel::distributed::SolutionTransfer for more information.

Note
This class, function, or variable is a template, and it can only be instantiated if the following condition is true: If your compiler supports the C++20 standard, then this constraint will be enforced by a C++20 requires clause.

Definition at line 316 of file dof_handler.h.

Member Typedef Documentation

◆ ActiveSelector

template<int dim, int spacedim = dim>
using DoFHandler< dim, spacedim >::ActiveSelector = ::internal::DoFHandlerImplementation::Iterators<dim, spacedim, false>
private

Definition at line 318 of file dof_handler.h.

◆ LevelSelector

template<int dim, int spacedim = dim>
using DoFHandler< dim, spacedim >::LevelSelector = ::internal::DoFHandlerImplementation::Iterators<dim, spacedim, true>
private

Definition at line 320 of file dof_handler.h.

◆ level_cell_accessor

template<int dim, int spacedim = dim>
using DoFHandler< dim, spacedim >::level_cell_accessor = typename LevelSelector::CellAccessor

Definition at line 501 of file dof_handler.h.

◆ level_face_accessor

template<int dim, int spacedim = dim>
using DoFHandler< dim, spacedim >::level_face_accessor = typename LevelSelector::FaceAccessor

Definition at line 502 of file dof_handler.h.

◆ level_cell_iterator

template<int dim, int spacedim = dim>
using DoFHandler< dim, spacedim >::level_cell_iterator = typename LevelSelector::cell_iterator

Definition at line 504 of file dof_handler.h.

◆ level_face_iterator

template<int dim, int spacedim = dim>
using DoFHandler< dim, spacedim >::level_face_iterator = typename LevelSelector::face_iterator

Definition at line 505 of file dof_handler.h.

◆ offset_type

template<int dim, int spacedim = dim>
using DoFHandler< dim, spacedim >::offset_type = unsigned int

The type in which we store the offsets in the CRS data structures.

Definition at line 526 of file dof_handler.h.

◆ map_value_type

using EnableObserverPointer::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 238 of file enable_observer_pointer.h.

◆ map_iterator

using EnableObserverPointer::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 243 of file enable_observer_pointer.h.

Constructor & Destructor Documentation

◆ DoFHandler() [1/3]

template<int dim, int spacedim = dim>
DoFHandler< dim, spacedim >::DoFHandler ( )

Standard constructor, not initializing any data. After constructing an object with this constructor, use reinit() to get a valid DoFHandler.

◆ DoFHandler() [2/3]

template<int dim, int spacedim = dim>
DoFHandler< dim, spacedim >::DoFHandler ( const Triangulation< dim, spacedim > &  tria)
explicit

Constructor. Take tria as the triangulation to work on.

◆ DoFHandler() [3/3]

template<int dim, int spacedim = dim>
DoFHandler< dim, spacedim >::DoFHandler ( const DoFHandler< dim, spacedim > &  )
delete

Copy constructor. DoFHandler objects are large and expensive. They should not be copied, in particular not by accident, but rather deliberately constructed. As a consequence, this constructor is explicitly removed from the interface of this class.

◆ ~DoFHandler()

template<int dim, int spacedim = dim>
virtual DoFHandler< dim, spacedim >::~DoFHandler ( )
overridevirtual

Destructor.

Member Function Documentation

◆ operator=()

template<int dim, int spacedim = dim>
DoFHandler & DoFHandler< dim, spacedim >::operator= ( const DoFHandler< dim, spacedim > &  )
delete

Copy operator. DoFHandler objects are large and expensive. They should not be copied, in particular not by accident, but rather deliberately constructed. As a consequence, this operator is explicitly removed from the interface of this class.

◆ set_active_fe_indices() [1/2]

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::set_active_fe_indices ( const std::vector< types::fe_index > &  active_fe_indices)

For each locally owned cell, set the active finite element index to the corresponding value given in active_fe_indices.

The vector active_fe_indices needs to have as many entries as there are active cells. The FE indices must be in the order in which we iterate over active cells. Vector entries corresponding to active cells that are not locally owned are ignored.

Active FE indices will only be set for locally owned cells. Ghost and artificial cells will be ignored; no active FE index will be assigned to them. To exchange active FE indices on ghost cells, call distribute_dofs() afterwards.

◆ set_active_fe_indices() [2/2]

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::set_active_fe_indices ( const std::vector< unsigned int > &  active_fe_indices)

For each locally owned cell, set the active finite element index to the corresponding value given in active_fe_indices.

The vector active_fe_indices needs to have as many entries as there are active cells. The FE indices must be in the order in which we iterate over active cells. Vector entries corresponding to active cells that are not locally owned are ignored.

Active FE indices will only be set for locally owned cells. Ghost and artificial cells will be ignored; no active FE index will be assigned to them. To exchange active FE indices on ghost cells, call distribute_dofs() afterwards.

Deprecated:
Use set_active_fe_indices() with the types::fe_index datatype.

◆ get_active_fe_indices() [1/2]

template<int dim, int spacedim = dim>
std::vector< types::fe_index > DoFHandler< dim, spacedim >::get_active_fe_indices ( ) const

For each locally relevant cell, extract the active finite element index and return them in the order in which we iterate over active cells.

As we do not know the active FE index on artificial cells, they are set to the invalid value numbers::invalid_fe_index.

For DoFHandler objects without hp-capabilities, the vector will consist of zeros, indicating that all cells use the same finite element. In hp-mode, the values may be different, though.

The returned vector has as many entries as there are active cells.

◆ get_active_fe_indices() [2/2]

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::get_active_fe_indices ( std::vector< unsigned int > &  active_fe_indices) const

For each locally relevant cell, extract the active finite element index and fill the vector active_fe_indices in the order in which we iterate over active cells. This vector is resized, if necessary.

As we do not know the active FE index on artificial cells, they are set to the invalid value numbers::invalid_fe_index.

For DoFHandler objects without hp-capabilities, the vector will consist of zeros, indicating that all cells use the same finite element. In hp-mode, the values may be different, though.

The returned vector has as many entries as there are active cells.

Deprecated:
Use get_active_fe_indices() that returns the result vector.

◆ set_future_fe_indices()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::set_future_fe_indices ( const std::vector< types::fe_index > &  future_fe_indices)

For each locally owned cell, set the future finite element index to the corresponding value given in future_fe_indices.

The vector future_fe_indices needs to have as many entries as there are active cells. The FE indices must be in the order in which we iterate over active cells. Vector entries corresponding to active cells that are not locally owned are ignored.

Future FE indices will only be set for locally owned cells. Ghost and artificial cells will be ignored; no future FE index will be assigned to them.

◆ get_future_fe_indices()

template<int dim, int spacedim = dim>
std::vector< types::fe_index > DoFHandler< dim, spacedim >::get_future_fe_indices ( ) const

For each locally owned cell, extract the future finite element index and return them in the order in which we iterate over active cells.

As we do not know the future FE index on ghost and artificial cells, they are set to the invalid value numbers::invalid_fe_index. The same applies to locally owned cells that have no future FE index assigned.

The returned vector has as many entries as there are active cells.

◆ reinit()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::reinit ( const Triangulation< dim, spacedim > &  tria)

Assign a Triangulation to the DoFHandler.

Remove all associations with the previous Triangulation object and establish connections with the new one. All information about previous degrees of freedom will be removed. Activates hp-mode.

◆ distribute_dofs() [1/2]

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::distribute_dofs ( const FiniteElement< dim, spacedim > &  fe)

Go through the triangulation and "distribute" the degrees of freedom needed for the given finite element. "Distributing" degrees of freedom involves allocating memory to store the indices on all entities on which degrees of freedom can be located (e.g., vertices, edges, faces, etc.) and to then enumerate all degrees of freedom. In other words, while the mesh and the finite element object by themselves simply define a finite element space \(V_h\), the process of distributing degrees of freedom makes sure that there is a basis for this space and that the shape functions of this basis are enumerated in an indexable, predictable way.

The exact order in which degrees of freedom on a mesh are ordered, i.e., the order in which basis functions of the finite element space are enumerated, is something that deal.II treats as an implementation detail. By and large, degrees of freedom are enumerated in the same order in which we traverse cells, but you should not rely on any specific numbering. In contrast, if you want a particular ordering, use the functions in namespace DoFRenumbering.

This function is first discussed in the introduction to the step-2 tutorial program.

Note
This function makes a copy of the finite element given as argument, and stores it as a member variable, similarly to the above function set_fe().

◆ distribute_dofs() [2/2]

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::distribute_dofs ( const hp::FECollection< dim, spacedim > &  fe)

Same as above but taking an hp::FECollection object.

◆ distribute_mg_dofs()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::distribute_mg_dofs ( )

Distribute level degrees of freedom on each level for geometric multigrid. The active DoFs need to be distributed using distribute_dofs() before calling this function.

◆ has_hp_capabilities()

template<int dim, int spacedim = dim>
bool DoFHandler< dim, spacedim >::has_hp_capabilities ( ) const

Returns whether this DoFHandler has hp-capabilities.

◆ has_level_dofs()

template<int dim, int spacedim = dim>
bool DoFHandler< dim, spacedim >::has_level_dofs ( ) const

This function returns whether this DoFHandler has DoFs distributed on each multigrid level or in other words if distribute_mg_dofs() has been called.

◆ has_active_dofs()

template<int dim, int spacedim = dim>
bool DoFHandler< dim, spacedim >::has_active_dofs ( ) const

This function returns whether this DoFHandler has active DoFs. This is equivalent to asking whether (i) distribute_dofs() has been called and (ii) the finite element for which degrees of freedom have been distributed actually has degrees of freedom (which is not the case for FE_Nothing, for example).

If this object is based on a parallel::distributed::Triangulation, then the current function returns true if any partition of the parallel DoFHandler object has any degrees of freedom. In other words, the function returns true even if the Triangulation does not own any active cells on the current MPI process, but at least one process owns cells and at least this one process has any degrees of freedom associated with it.

◆ initialize_local_block_info()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::initialize_local_block_info ( )

After distribute_dofs() with an FESystem element, the block structure of global and level vectors is stored in a BlockInfo object accessible with block_info(). This function initializes the local block structure on each cell in the same object.

◆ clear()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::clear ( )

Clear all data of this object.

◆ renumber_dofs() [1/2]

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::renumber_dofs ( const std::vector< types::global_dof_index > &  new_numbers)

Renumber degrees of freedom based on a list of new DoF indices for each of the degrees of freedom.

This function is called by the functions in DoFRenumbering function after computing a new ordering of the degree of freedom indices. However, it can of course also be called from user code.

  • new_number This array must have a size equal to the number of degrees of freedom owned by the current processor, i.e. the size must be equal to what n_locally_owned_dofs() returns. If only one processor participates in storing the current mesh, then this equals the total number of degrees of freedom, i.e. the result of n_dofs(). The contents of this array are the new global indices for each freedom listed in the IndexSet returned by locally_owned_dofs(). In the case of a sequential mesh this means that the array is a list of new indices for each of the degrees of freedom on the current mesh. In the case that we have a parallel::shared::Triangulation or parallel::distributed::Triangulation underlying this DoFHandler object, the array is a list of new indices for all the locally owned degrees of freedom, enumerated in the same order as the currently locally owned DoFs. In other words, assume that degree of freedom i is currently locally owned, then new_numbers[locally_owned_dofs().index_within_set(i)] returns the new global DoF index of i. Since the IndexSet of locally_owned_dofs() is complete in the sequential case, the latter convention for the content of the array reduces to the former in the case that only one processor participates in the mesh.
Note
While it follows from the above, it may be surprising to know that the number of locally owned degrees of freedom in a parallel computation is an invariant under renumbering, even if the indices associated with these locally owned degrees of freedom are not. At a fundamental level, this invariant exists because the decision whether a degree of freedom is locally owned or not has nothing to do with that degree of freedom's (old or new) index. Indeed, degrees of freedom are locally owned if they are on a locally owned cell and not on an interface between cells where the neighboring cell has a lower subdomain id. Since both of these conditions are independent of the index associated with the DoF, a locally owned degree of freedom will also be locally owned after renumbering. On the other hand, properties such as whether the set of indices of locally owned DoFs forms a contiguous range or not (i.e., whether the locally_owned_dofs() returns an IndexSet object for which IndexSet::is_contiguous() returns true) are of course affected by the exact renumbering performed here. For example, while the initial numbering of DoF indices done in distribute_dofs() yields a contiguous numbering, the renumbering performed by DoFRenumbering::component_wise() will, in general, not yield contiguous locally owned DoF indices.

◆ renumber_dofs() [2/2]

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::renumber_dofs ( const unsigned int  level,
const std::vector< types::global_dof_index > &  new_numbers 
)

The same function as above, but renumber the degrees of freedom of a single level of a multigrid hierarchy.

◆ max_couplings_between_dofs()

template<int dim, int spacedim = dim>
unsigned int DoFHandler< dim, spacedim >::max_couplings_between_dofs ( ) const

Return the maximum number of degrees of freedom a degree of freedom in the given triangulation with the given finite element may couple with. This is the maximum number of entries per line in the system matrix; this information can therefore be used upon construction of the SparsityPattern object.

The returned number is not really the maximum number but an estimate based on the finite element and the maximum number of cells meeting at a vertex. The number holds for the constrained matrix as well.

The determination of the number of couplings can be done by simple picture drawing. An example can be found in the implementation of this function.

Note
This function is most often used to determine the maximal row length for sparsity patterns. Unfortunately, while the estimates returned by this function are rather accurate in 1d and 2d, they are often significantly too high in 3d, leading the SparsityPattern class to allocate much too much memory in some cases. Unless someone comes around to improving the present function for 3d, there is not very much one can do about these cases. The typical way to work around this problem is to use an intermediate compressed sparsity pattern that only allocates memory on demand. Refer to the step-2 and step-11 example programs on how to do this. The problem is also discussed in the documentation of the topic on Sparsity patterns.

◆ max_couplings_between_boundary_dofs()

template<int dim, int spacedim = dim>
unsigned int DoFHandler< dim, spacedim >::max_couplings_between_boundary_dofs ( ) const

Return the number of degrees of freedom located on the boundary another dof on the boundary can couple with.

The number is the same as for max_couplings_between_dofs() in one dimension less.

Note
The same applies to this function as to max_couplings_between_dofs() as regards the performance of this function. Think about one of the dynamic sparsity pattern classes instead (see Sparsity patterns).

◆ begin()

template<int dim, int spacedim = dim>
cell_iterator DoFHandler< dim, spacedim >::begin ( const unsigned int  level = 0) const

Iterator to the first used cell on level level.

◆ begin_active()

template<int dim, int spacedim = dim>
active_cell_iterator DoFHandler< dim, spacedim >::begin_active ( const unsigned int  level = 0) const

Iterator to the first active cell on level level. If the given level does not contain any active cells (i.e., all cells on this level are further refined), then this function returns end_active(level) so that loops of the kind

for (cell=dof_handler.begin_active(level);
cell!=dof_handler.end_active(level);
++cell)
{
...
}
unsigned int level
Definition grid_out.cc:4626

have zero iterations, as may be expected if there are no active cells on this level.

◆ end() [1/2]

template<int dim, int spacedim = dim>
cell_iterator DoFHandler< dim, spacedim >::end ( ) const

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ end() [2/2]

template<int dim, int spacedim = dim>
cell_iterator DoFHandler< dim, spacedim >::end ( const unsigned int  level) const

Return an iterator which is the first iterator not on the given level. If level is the last level, then this returns end().

◆ end_active()

template<int dim, int spacedim = dim>
active_cell_iterator DoFHandler< dim, spacedim >::end_active ( const unsigned int  level) const

Return an active iterator which is the first active iterator not on the given level. If level is the last level, then this returns end().

◆ begin_mg()

template<int dim, int spacedim = dim>
level_cell_iterator DoFHandler< dim, spacedim >::begin_mg ( const unsigned int  level = 0) const

Iterator to the first used cell on level level. This returns a level_cell_iterator that returns level dofs when dof_indices() is called.

◆ end_mg() [1/2]

template<int dim, int spacedim = dim>
level_cell_iterator DoFHandler< dim, spacedim >::end_mg ( const unsigned int  level) const

Iterator past the last cell on level level. This returns a level_cell_iterator that returns level dofs when dof_indices() is called.

◆ end_mg() [2/2]

template<int dim, int spacedim = dim>
level_cell_iterator DoFHandler< dim, spacedim >::end_mg ( ) const

Iterator past the end; this iterator serves for comparisons of iterators with past-the-end or before-the-beginning states.

◆ n_dofs() [1/2]

template<int dim, int spacedim = dim>
types::global_dof_index DoFHandler< dim, spacedim >::n_dofs ( ) const

Return the global number of degrees of freedom. If the current object handles all degrees of freedom itself (even if you may intend to solve your linear system in parallel, such as in step-17 or step-18), then this number equals the number of locally owned degrees of freedom since this object doesn't know anything about what you want to do with it and believes that it owns every degree of freedom it knows about.

On the other hand, if this object operates on a parallel::distributed::Triangulation object, then this function returns the global number of degrees of freedom, accumulated over all processors.

In either case, included in the returned number are those DoFs which are constrained by hanging nodes, see Constraints on degrees of freedom.

Mathematically speaking, the number returned by this function equals the dimension of the finite element space (without taking into account constraints) that corresponds to (i) the mesh on which it is defined, and (ii) the finite element that is used by the current object. It also, of course, equals the number of shape functions that span this space.

◆ n_dofs() [2/2]

template<int dim, int spacedim = dim>
types::global_dof_index DoFHandler< dim, spacedim >::n_dofs ( const unsigned int  level) const

The (global) number of multilevel degrees of freedom on a given level.

If no level degrees of freedom have been assigned to this level, returns numbers::invalid_dof_index. Else returns the number of degrees of freedom on this level.

◆ n_boundary_dofs() [1/3]

template<int dim, int spacedim = dim>
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs ( ) const

Return the number of locally owned degrees of freedom located on the boundary.

◆ n_boundary_dofs() [2/3]

template<int dim, int spacedim = dim>
template<typename number >
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs ( const std::map< types::boundary_id, const Function< spacedim, number > * > &  boundary_ids) const

Return the number of locally owned degrees of freedom located on those parts of the boundary which have a boundary indicator listed in the given set. The reason that a map rather than a set is used is the same as described in the documentation of that variant of DoFTools::make_boundary_sparsity_pattern() that takes a map.

There is, however, another overload of this function that takes a set argument (see below).

◆ n_boundary_dofs() [3/3]

template<int dim, int spacedim = dim>
types::global_dof_index DoFHandler< dim, spacedim >::n_boundary_dofs ( const std::set< types::boundary_id > &  boundary_ids) const

Return the number of degrees of freedom located on those parts of the boundary which have a boundary indicator listed in the given set. The

◆ block_info()

template<int dim, int spacedim = dim>
const BlockInfo & DoFHandler< dim, spacedim >::block_info ( ) const

Access to an object informing of the block structure of the dof handler.

If an FESystem is used in distribute_dofs(), degrees of freedom naturally split into several blocks. For each base element as many blocks appear as its multiplicity.

At the end of distribute_dofs(), the number of degrees of freedom in each block is counted, and stored in a BlockInfo object, which can be accessed here. If you have previously called distribute_mg_dofs(), the same is done on each level of the multigrid hierarchy. Additionally, the block structure on each cell can be generated in this object by calling initialize_local_block_info().

◆ n_locally_owned_dofs()

template<int dim, int spacedim = dim>
types::global_dof_index DoFHandler< dim, spacedim >::n_locally_owned_dofs ( ) const

Return the number of degrees of freedom that belong to this process.

If this is a sequential DoFHandler, then the result equals that produced by n_dofs(). (Here, "sequential" means that either the whole program does not use MPI, or that it uses MPI but only uses a single MPI process, or that there are multiple MPI processes but the Triangulation on which this DoFHandler builds works only on one MPI process.) On the other hand, if we are operating on a parallel::distributed::Triangulation or parallel::shared::Triangulation, then it includes only the degrees of freedom that the current processor owns. Note that in this case this does not include all degrees of freedom that have been distributed on the current processor's image of the mesh: in particular, some of the degrees of freedom on the interface between the cells owned by this processor and cells owned by other processors may be theirs, and degrees of freedom on ghost cells are also not necessarily included.

◆ locally_owned_dofs()

template<int dim, int spacedim = dim>
const IndexSet & DoFHandler< dim, spacedim >::locally_owned_dofs ( ) const

Return an IndexSet describing the set of locally owned DoFs as a subset of 0..n_dofs(). The number of elements of this set equals n_locally_owned_dofs().

◆ locally_owned_mg_dofs()

template<int dim, int spacedim = dim>
const IndexSet & DoFHandler< dim, spacedim >::locally_owned_mg_dofs ( const unsigned int  level) const

Return an IndexSet describing the set of locally owned DoFs used for the given multigrid level as a subset of 0..n_dofs(level).

◆ get_fe()

template<int dim, int spacedim = dim>
const FiniteElement< dim, spacedim > & DoFHandler< dim, spacedim >::get_fe ( const types::fe_index  index = 0) const

Return a constant reference to the indexth finite element object that is used by this object.

◆ get_fe_collection()

template<int dim, int spacedim = dim>
const hp::FECollection< dim, spacedim > & DoFHandler< dim, spacedim >::get_fe_collection ( ) const

Return a constant reference to the set of finite element objects that are used by this object.

◆ get_triangulation()

template<int dim, int spacedim = dim>
const Triangulation< dim, spacedim > & DoFHandler< dim, spacedim >::get_triangulation ( ) const

Return a constant reference to the triangulation underlying this object.

◆ get_mpi_communicator()

template<int dim, int spacedim = dim>
MPI_Comm DoFHandler< dim, spacedim >::get_mpi_communicator ( ) const

Return MPI communicator used by the underlying triangulation.

◆ get_communicator()

template<int dim, int spacedim = dim>
MPI_Comm DoFHandler< dim, spacedim >::get_communicator ( ) const

Return MPI communicator used by the underlying triangulation.

Deprecated:
Use get_mpi_communicator() instead.

◆ prepare_for_serialization_of_active_fe_indices()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::prepare_for_serialization_of_active_fe_indices ( )

Whenever serialization with a parallel::distributed::Triangulation as the underlying triangulation is considered, we also need to consider storing the active FE indices on all active cells as well.

This function registers that these indices are to be stored whenever the parallel::distributed::Triangulation::save() function is called on the underlying triangulation.

Note
Currently only implemented for triangulations of type parallel::distributed::Triangulation. An assertion will be triggered if a different type is registered.
See also
The documentation of parallel::distributed::SolutionTransfer has further information on serialization.

◆ deserialize_active_fe_indices()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::deserialize_active_fe_indices ( )

Whenever serialization with a parallel::distributed::Triangulation as the underlying triangulation is considered, we also need to consider storing the active FE indices on all active cells as well.

This function deserializes and distributes the previously stored active FE indices on all active cells.

Note
Currently only implemented for triangulations of type parallel::distributed::Triangulation. An assertion will be triggered if a different type is registered.
See also
The documentation of parallel::distributed::SolutionTransfer has further information on serialization.

◆ memory_consumption()

template<int dim, int spacedim = dim>
virtual std::size_t DoFHandler< dim, spacedim >::memory_consumption ( ) const
virtual

Determine an estimate for the memory consumption (in bytes) of this object.

This function is made virtual, since a dof handler object might be accessed through a pointers to this base class, although the actual object might be a derived class.

◆ save()

template<int dim, int spacedim = dim>
template<class Archive >
void DoFHandler< dim, spacedim >::save ( Archive &  ar,
const unsigned int  version 
) const

Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.

◆ load()

template<int dim, int spacedim = dim>
template<class Archive >
void DoFHandler< dim, spacedim >::load ( Archive &  ar,
const unsigned int  version 
)

Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ serialize()

template<int dim, int spacedim = dim>
template<class Archive >
void DoFHandler< dim, spacedim >::serialize ( Archive &  archive,
const unsigned int  version 
)

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ clear_space()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::clear_space ( )
private

Free all memory used for non-multigrid data structures.

◆ clear_mg_space()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::clear_mg_space ( )
private

Free all memory used for multigrid data structures.

◆ setup_policy()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::setup_policy ( )
private

Set up DoFHandler policy.

◆ connect_to_triangulation_signals()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::connect_to_triangulation_signals ( )
private

Set up connections to signals of the underlying triangulation.

◆ create_active_fe_table()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::create_active_fe_table ( )
private

Create default tables for the active and future fe_indices.

Active indices are initialized with a zero indicator, meaning that fe[0] is going to be used by default. Future indices are initialized with an invalid indicator, meaning that no p-adaptation is scheduled by default.

This method is called upon construction and whenever the underlying triangulation gets created. This ensures that each cell has a valid active and future fe_index.

◆ update_active_fe_table()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::update_active_fe_table ( )
private

Update tables for active and future fe_indices.

Whenever the underlying triangulation changes (either by adaptation or deserialization), active and future FE index tables will be adjusted to the current structure of the triangulation. Missing values of active and future indices will be initialized with their defaults (see create_active_fe_table()).

This method is called post refinement and post deserialization. This ensures that each cell has a valid active and future fe_index.

◆ pre_transfer_action()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::pre_transfer_action ( )
private

A function that will be triggered through a triangulation signal just before the associated Triangulation or parallel::shared::Triangulation is modified.

The function that stores the active FE indices of all cells that will be refined or coarsened before the refinement happens, so that they can be set again after refinement.

◆ post_transfer_action()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::post_transfer_action ( )
private

A function that will be triggered through a triangulation signal just after the associated Triangulation or parallel::shared::Triangulation is modified.

The function that restores the active FE indices of all cells that were refined or coarsened.

◆ pre_distributed_transfer_action()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::pre_distributed_transfer_action ( )
private

A function that will be triggered through a triangulation signal just before the associated parallel::distributed::Triangulation is modified.

The function that stores all active FE indices on locally owned cells for distribution over all participating processors.

◆ post_distributed_transfer_action()

template<int dim, int spacedim = dim>
void DoFHandler< dim, spacedim >::post_distributed_transfer_action ( )
private

A function that will be triggered through a triangulation signal just after the associated parallel::distributed::Triangulation is modified.

The function that restores all active FE indices on locally owned cells that have been communicated.

◆ subscribe()

void EnableObserverPointer::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 131 of file enable_observer_pointer.cc.

◆ unsubscribe()

void EnableObserverPointer::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 151 of file enable_observer_pointer.cc.

◆ n_subscriptions()

unsigned int EnableObserverPointer::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 322 of file enable_observer_pointer.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void EnableObserverPointer::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 339 of file enable_observer_pointer.h.

◆ list_subscribers() [2/2]

void EnableObserverPointer::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 199 of file enable_observer_pointer.cc.

◆ check_no_subscribers()

void EnableObserverPointer::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 53 of file enable_observer_pointer.cc.

Friends And Related Symbol Documentation

◆ ::DoFAccessor

template<int dim, int spacedim = dim>
template<int , int , int , bool >
friend class ::DoFAccessor
friend

Definition at line 1710 of file dof_handler.h.

◆ ::DoFCellAccessor

template<int dim, int spacedim = dim>
template<int , int , bool >
friend class ::DoFCellAccessor
friend

Definition at line 1712 of file dof_handler.h.

◆ ::internal::DoFAccessorImplementation::Implementation

template<int dim, int spacedim = dim>
friend struct ::internal::DoFAccessorImplementation::Implementation
friend

Definition at line 1713 of file dof_handler.h.

◆ ::internal::DoFCellAccessorImplementation::Implementation

template<int dim, int spacedim = dim>
friend struct ::internal::DoFCellAccessorImplementation::Implementation
friend

Definition at line 1714 of file dof_handler.h.

◆ ::internal::DoFHandlerImplementation::Implementation

template<int dim, int spacedim = dim>
friend struct ::internal::DoFHandlerImplementation::Implementation
friend

Definition at line 1718 of file dof_handler.h.

◆ ::internal::hp::DoFHandlerImplementation::Implementation

template<int dim, int spacedim = dim>
friend struct ::internal::hp::DoFHandlerImplementation::Implementation
friend

Definition at line 1719 of file dof_handler.h.

◆ ::internal::DoFHandlerImplementation::Policy::Implementation

template<int dim, int spacedim = dim>
friend struct ::internal::DoFHandlerImplementation::Policy::Implementation
friend

Definition at line 1720 of file dof_handler.h.

Member Data Documentation

◆ dimension

template<int dim, int spacedim = dim>
constexpr unsigned int DoFHandler< dim, spacedim >::dimension = dim
staticconstexpr

Make the dimension available in function templates.

Definition at line 511 of file dof_handler.h.

◆ space_dimension

template<int dim, int spacedim = dim>
constexpr unsigned int DoFHandler< dim, spacedim >::space_dimension = spacedim
staticconstexpr

Make the space dimension available in function templates.

Definition at line 516 of file dof_handler.h.

◆ default_fe_index

template<int dim, int spacedim = dim>
const types::fe_index DoFHandler< dim, spacedim >::default_fe_index = 0
static

The default index of the finite element to be used on a given cell.

Definition at line 521 of file dof_handler.h.

◆ block_info_object

template<int dim, int spacedim = dim>
BlockInfo DoFHandler< dim, spacedim >::block_info_object
private

An object containing information on the block structure.

Definition at line 1483 of file dof_handler.h.

◆ hp_capability_enabled

template<int dim, int spacedim = dim>
bool DoFHandler< dim, spacedim >::hp_capability_enabled
private

Boolean indicating whether or not the current DoFHandler has hp-capabilities.

Definition at line 1489 of file dof_handler.h.

◆ tria

template<int dim, int spacedim = dim>
ObserverPointer<const Triangulation<dim, spacedim>, DoFHandler<dim, spacedim> > DoFHandler< dim, spacedim >::tria
private

Address of the triangulation to work on.

Definition at line 1495 of file dof_handler.h.

◆ fe_collection

template<int dim, int spacedim = dim>
hp::FECollection<dim, spacedim> DoFHandler< dim, spacedim >::fe_collection
private

Store a hp::FECollection object. If only a single FiniteElement is used during initialization of this object, it contains the (one) FiniteElement.

Definition at line 1502 of file dof_handler.h.

◆ policy

template<int dim, int spacedim = dim>
std::unique_ptr<::internal::DoFHandlerImplementation::Policy:: PolicyBase<dim, spacedim> > DoFHandler< dim, spacedim >::policy
private

An object that describes how degrees of freedom should be distributed and renumbered.

Definition at line 1510 of file dof_handler.h.

◆ number_cache

template<int dim, int spacedim = dim>
::internal::DoFHandlerImplementation::NumberCache DoFHandler< dim, spacedim >::number_cache
private

A structure that contains all sorts of numbers that characterize the degrees of freedom this object works on.

For most members of this structure, there is an accessor function in this class that returns its value.

Definition at line 1519 of file dof_handler.h.

◆ mg_number_cache

template<int dim, int spacedim = dim>
std::vector<::internal::DoFHandlerImplementation::NumberCache> DoFHandler< dim, spacedim >::mg_number_cache
private

Data structure like number_cache, but for each multigrid level.

Definition at line 1525 of file dof_handler.h.

◆ object_dof_indices

template<int dim, int spacedim = dim>
std::vector<std::array<std::vector<types::global_dof_index>, dim + 1> > DoFHandler< dim, spacedim >::object_dof_indices
mutableprivate

Indices of degree of freedom of each d+1 geometric object (3d: vertex, line, quad, hex) for all relevant active finite elements. Identification of the appropriate position is done via object_dof_ptr (CRS scheme).

Definition at line 1533 of file dof_handler.h.

◆ object_dof_ptr

template<int dim, int spacedim = dim>
std::vector<std::array<std::vector<offset_type>, dim + 1> > DoFHandler< dim, spacedim >::object_dof_ptr
mutableprivate

Pointer to the first cached degree of freedom of a geometric object for all relevant active finite elements.

Note
In normal mode it is possible to access this data structure directly. In hp-mode, an indirection via hp_object_fe_indices/hp_object_fe_ptr is necessary.

Definition at line 1544 of file dof_handler.h.

◆ hp_object_fe_indices

template<int dim, int spacedim = dim>
std::array<std::vector<types::fe_index>, dim + 1> DoFHandler< dim, spacedim >::hp_object_fe_indices
mutableprivate

Active FE indices of each geometric object. Identification of the appropriate position of a cell in the vectors is done via hp_object_fe_ptr (CRS scheme).

Definition at line 1552 of file dof_handler.h.

◆ hp_object_fe_ptr

template<int dim, int spacedim = dim>
std::array<std::vector<offset_type>, dim + 1> DoFHandler< dim, spacedim >::hp_object_fe_ptr
mutableprivate

Pointer to the first FE index of a geometric object.

Definition at line 1557 of file dof_handler.h.

◆ hp_cell_active_fe_indices

template<int dim, int spacedim = dim>
std::vector<std::vector<types::fe_index> > DoFHandler< dim, spacedim >::hp_cell_active_fe_indices
mutableprivate

Active FE index of an active cell (identified by level and level index). This vector is only used in hp-mode.

Definition at line 1563 of file dof_handler.h.

◆ hp_cell_future_fe_indices

template<int dim, int spacedim = dim>
std::vector<std::vector<types::fe_index> > DoFHandler< dim, spacedim >::hp_cell_future_fe_indices
mutableprivate

Future FE index of an active cell (identified by level and level index). This vector is only used in hp-mode.

Definition at line 1569 of file dof_handler.h.

◆ mg_vertex_dofs

template<int dim, int spacedim = dim>
std::vector<MGVertexDoFs> DoFHandler< dim, spacedim >::mg_vertex_dofs
private

An array to store the indices for level degrees of freedom located at vertices.

Definition at line 1575 of file dof_handler.h.

◆ mg_levels

template<int dim, int spacedim = dim>
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel<dim> > > DoFHandler< dim, spacedim >::mg_levels
private

Space to store the DoF numbers for the different multigrid levels.

Definition at line 1582 of file dof_handler.h.

◆ mg_faces

template<int dim, int spacedim = dim>
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces<dim> > DoFHandler< dim, spacedim >::mg_faces
private

Space to store DoF numbers of faces in the multigrid context.

Definition at line 1588 of file dof_handler.h.

◆ active_fe_index_transfer

template<int dim, int spacedim = dim>
std::unique_ptr<ActiveFEIndexTransfer> DoFHandler< dim, spacedim >::active_fe_index_transfer
private

We embed our data structure into a pointer to control that all transfer related data only exists during the actual transfer process.

Definition at line 1594 of file dof_handler.h.

◆ tria_listeners

template<int dim, int spacedim = dim>
std::vector<boost::signals2::connection> DoFHandler< dim, spacedim >::tria_listeners
private

A list of connections with which this object connects to the triangulation to get information about when the triangulation changes.

Definition at line 1600 of file dof_handler.h.

◆ tria_listeners_for_transfer

template<int dim, int spacedim = dim>
std::vector<boost::signals2::connection> DoFHandler< dim, spacedim >::tria_listeners_for_transfer
private

A list of connections with which this object connects to the triangulation. They get triggered specifically when data needs to be transferred due to refinement or repartitioning. Only active in hp-mode.

Definition at line 1607 of file dof_handler.h.

◆ counter

std::atomic<unsigned int> EnableObserverPointer::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 227 of file enable_observer_pointer.h.

◆ counter_map

std::map<std::string, unsigned int> EnableObserverPointer::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 233 of file enable_observer_pointer.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> EnableObserverPointer::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the ObserverPointer objects that subscribe to this class.

Definition at line 249 of file enable_observer_pointer.h.

◆ object_info

const std::type_info* EnableObserverPointer::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 257 of file enable_observer_pointer.h.

◆ mutex

std::mutex EnableObserverPointer::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 280 of file enable_observer_pointer.h.


The documentation for this class was generated from the following files: