Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
DoFAccessor< structdim, DoFHandlerType, level_dof_access > Class Template Reference

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

Inheritance diagram for DoFAccessor< structdim, DoFHandlerType, level_dof_access >:
[legend]

Public Types

using BaseClass = typename ::internal::DoFAccessorImplementation::Inheritance< structdim, dimension, space_dimension >::BaseClass
 
using AccessorData = DoFHandlerType
 
- Public Types inherited from TriaAccessor< structdim, dim, spacedim >
using AccessorData = typename TriaAccessorBase< structdim, dim, spacedim >::AccessorData
 
- Public Types inherited from TriaAccessorBase< structdim, dim, spacedim >
using LocalData = void *
 

Public Member Functions

const DoFHandlerType & get_dof_handler () const
 
template<bool level_dof_access2>
void copy_from (const DoFAccessor< structdim, DoFHandlerType, level_dof_access2 > &a)
 
void copy_from (const TriaAccessorBase< structdim, DoFHandlerType::dimension, DoFHandlerType::space_dimension > &da)
 
template<int dim2, class DoFHandlerType2 , bool level_dof_access2>
bool operator== (const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) const
 
template<int dim2, class DoFHandlerType2 , bool level_dof_access2>
bool operator!= (const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &) const
 
Constructors
 DoFAccessor ()
 
 DoFAccessor (const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > *tria, const int level, const int index, const DoFHandlerType *dof_handler)
 
template<int structdim2, int dim2, int spacedim2>
 DoFAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &)
 
template<int dim2, class DoFHandlerType2 , bool level_dof_access2>
 DoFAccessor (const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &)
 
template<bool level_dof_access2>
 DoFAccessor (const DoFAccessor< structdim, DoFHandlerType, level_dof_access2 > &)
 
DoFAccessor< structdim, DoFHandlerType, level_dof_access > & operator= (const DoFAccessor< structdim, DoFHandlerType, level_dof_access > &da)=delete
 
Accessing sub-objects
TriaIterator< DoFAccessor< structdim, DoFHandlerType, level_dof_access > > child (const unsigned int c) const
 
typename ::internal::DoFHandlerImplementation::Iterators< DoFHandlerType, level_dof_access >::line_iterator line (const unsigned int i) const
 
typename ::internal::DoFHandlerImplementation::Iterators< DoFHandlerType, level_dof_access >::quad_iterator quad (const unsigned int i) const
 
Accessing the DoF indices of this object
void get_dof_indices (std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
 
void get_mg_dof_indices (const int level, std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
 
void set_mg_dof_indices (const int level, const std::vector< types::global_dof_index > &dof_indices, const unsigned int fe_index=DoFHandlerType::default_fe_index)
 
types::global_dof_index vertex_dof_index (const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
 
types::global_dof_index mg_vertex_dof_index (const int level, const unsigned int vertex, const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
 
types::global_dof_index dof_index (const unsigned int i, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
 
types::global_dof_index mg_dof_index (const int level, const unsigned int i) const
 
Accessing the finite element associated with this object
unsigned int n_active_fe_indices () const
 
unsigned int nth_active_fe_index (const unsigned int n) const
 
std::set< unsigned int > get_active_fe_indices () const
 
bool fe_index_is_active (const unsigned int fe_index) const
 
const FiniteElement< DoFHandlerType::dimension, DoFHandlerType::space_dimension > & get_fe (const unsigned int fe_index) const
 
- Public Member Functions inherited from TriaAccessor< structdim, dim, spacedim >
 TriaAccessor (const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
 
template<int structdim2, int dim2, int spacedim2>
 TriaAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &)
 
template<int structdim2, int dim2, int spacedim2>
 TriaAccessor (const TriaAccessor< structdim2, dim2, spacedim2 > &)
 
void operator= (const TriaAccessor &)=delete
 
bool used () const
 
template<>
double extent_in_direction (const unsigned int axis) const
 
template<>
double extent_in_direction (const unsigned int axis) const
 
template<>
double extent_in_direction (const unsigned int axis) const
 
TriaIterator< TriaAccessor< 0, dim, spacedim > > vertex_iterator (const unsigned int i) const
 
unsigned int vertex_index (const unsigned int i) const
 
Point< spacedim > & vertex (const unsigned int i) const
 
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::line_iterator line (const unsigned int i) const
 
unsigned int line_index (const unsigned int i) const
 
typename ::internal::TriangulationImplementation::Iterators< dim, spacedim >::quad_iterator quad (const unsigned int i) const
 
unsigned int quad_index (const unsigned int i) const
 
bool face_orientation (const unsigned int face) const
 
bool face_flip (const unsigned int face) const
 
bool face_rotation (const unsigned int face) const
 
bool line_orientation (const unsigned int line) const
 
bool has_children () const
 
unsigned int n_children () const
 
unsigned int number_of_children () const
 
unsigned int max_refinement_depth () const
 
TriaIterator< TriaAccessor< structdim, dim, spacedim > > child (const unsigned int i) const
 
TriaIterator< TriaAccessor< structdim, dim, spacedim > > isotropic_child (const unsigned int i) const
 
RefinementCase< structdim > refinement_case () const
 
int child_index (const unsigned int i) const
 
int isotropic_child_index (const unsigned int i) const
 
types::boundary_id boundary_id () const
 
void set_boundary_id (const types::boundary_id) const
 
void set_all_boundary_ids (const types::boundary_id) const
 
bool at_boundary () const
 
const Manifold< dim, spacedim > & get_manifold () const
 
types::manifold_id manifold_id () const
 
void set_manifold_id (const types::manifold_id) const
 
void set_all_manifold_ids (const types::manifold_id) const
 
bool user_flag_set () const
 
void set_user_flag () const
 
void clear_user_flag () const
 
void recursively_set_user_flag () const
 
void recursively_clear_user_flag () const
 
void clear_user_data () const
 
void set_user_pointer (void *p) const
 
void clear_user_pointer () const
 
void * user_pointer () const
 
void recursively_set_user_pointer (void *p) const
 
void recursively_clear_user_pointer () const
 
void set_user_index (const unsigned int p) const
 
void clear_user_index () const
 
unsigned int user_index () const
 
void recursively_set_user_index (const unsigned int p) const
 
void recursively_clear_user_index () const
 
double diameter () const
 
std::pair< Point< spacedim >, double > enclosing_ball () const
 
BoundingBox< spacedim > bounding_box () const
 
double extent_in_direction (const unsigned int axis) const
 
double minimum_vertex_distance () const
 
Point< spacedim > intermediate_point (const Point< structdim > &coordinates) const
 
Point< structdim > real_to_unit_cell_affine_approximation (const Point< spacedim > &point) const
 
Point< spacedim > center (const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
 
Point< spacedim > barycenter () const
 
double measure () const
 
bool is_translation_of (const TriaIterator< TriaAccessor< structdim, dim, spacedim >> &o) const
 
- Public Member Functions inherited from TriaAccessorBase< structdim, dim, spacedim >
void operator= (const TriaAccessorBase *)=delete
 
int level () const
 
int index () const
 
IteratorState::IteratorStates state () const
 
const Triangulation< dim, spacedim > & get_triangulation () const
 

Static Public Member Functions

static bool is_level_cell ()
 
static ::ExceptionBaseExcInvalidObject ()
 
static ::ExceptionBaseExcVectorNotEmpty ()
 
static ::ExceptionBaseExcVectorDoesNotMatch ()
 
static ::ExceptionBaseExcMatrixDoesNotMatch ()
 
static ::ExceptionBaseExcNotActive ()
 
static ::ExceptionBaseExcCantCompareIterators ()
 

Static Public Attributes

static const unsigned int dimension = DoFHandlerType::dimension
 
static const unsigned int space_dimension = DoFHandlerType::space_dimension
 
- Static Public Attributes inherited from TriaAccessorBase< structdim, dim, spacedim >
static const unsigned int space_dimension = spacedim
 
static const unsigned int dimension = dim
 
static const unsigned int structure_dimension = structdim
 

Protected Member Functions

void set_dof_handler (DoFHandlerType *dh)
 
void set_dof_index (const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
 
void set_vertex_dof_index (const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const unsigned int fe_index=DoFHandlerType::default_fe_index) const
 
- Protected Member Functions inherited from TriaAccessorBase< structdim, dim, spacedim >
 TriaAccessorBase (const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *=nullptr)
 
 TriaAccessorBase (const TriaAccessorBase &)
 
void copy_from (const TriaAccessorBase &)
 
TriaAccessorBaseoperator= (const TriaAccessorBase &)
 
bool operator< (const TriaAccessorBase &other) const
 
bool operator== (const TriaAccessorBase &) const
 
bool operator!= (const TriaAccessorBase &) const
 
::internal::TriangulationImplementation::TriaObjects< ::internal::TriangulationImplementation::TriaObject< structdim > > & objects () const
 
void operator++ ()
 
void operator-- ()
 

Protected Attributes

DoFHandlerType * dof_handler
 
- Protected Attributes inherited from TriaAccessorBase< structdim, dim, spacedim >
typename ::internal::TriaAccessorImplementation::PresentLevelType< structdim, dim >::type present_level
 
int present_index
 
const Triangulation< dim, spacedim > * tria
 

Friends

template<typename >
class TriaRawIterator
 
template<int , class , bool >
class DoFAccessor
 
template<int dim, int spacedim>
class DoFHandler
 
template<int dim, int spacedim>
class hp::DoFHandler
 

Additional Inherited Members

- Protected Types inherited from TriaAccessorBase< structdim, dim, spacedim >
using AccessorData = void
 

Detailed Description

template<int structdim, typename DoFHandlerType, bool level_dof_access>
class DoFAccessor< structdim, DoFHandlerType, level_dof_access >

A class that gives access to the degrees of freedom stored in a DoFHandler or hp::DoFHandler object. Accessors are used to access the data that pertains to edges, faces, and cells of a triangulation. The concept is explained in more detail in connection to Iterators on mesh-like containers.

This class follows mainly the route laid out by the accessor library declared in the triangulation library (TriaAccessor). It enables the user to access the degrees of freedom on lines, quads, or hexes. The first template argument of this class determines the dimensionality of the object under consideration: 1 for lines, 2 for quads, and 3 for hexes. The second argument denotes the type of DoF handler we should work on. It can either be DoFHandler or hp::DoFHandler. From the second template argument we also deduce the dimension of the Triangulation this object refers to as well as the dimension of the space into which it is embedded. Finally, the template argument level_dof_access governs the behavior of the function get_active_or_mg_dof_indices(). See the section on Generic loops below.

Alias

Usage is best to happen through the alias to the various kinds of iterators provided by the DoFHandler and hp::DoFHandler classes, since they are more secure to changes in the class naming and template interface as well as providing easier typing (much less complicated names!).

Generic loops and the third template argument

Many loops look very similar, whether they operate on the active dofs of the active cells of the Triangulation or on the level dofs of a single level or the whole grid hierarchy. In order to use polymorphism in such loops, they access degrees of freedom through the function get_active_or_mg_dof_indices(), which changes behavior according to the third template argument. If the argument is false, then the active dofs of active cells are accessed. If it is true, the level dofs are used. DoFHandler has functions, for instance begin() and begin_mg(), which return either type or the other. Additionally, they can be cast into each other, in case this is needed, since they access the same data.

It is recommended to use the function get_active_or_mg_dof_indices() in generic loops in lieu of get_dof_indices() or get_mg_dof_indices().

Inheritance

If the structural dimension given by the first template argument equals the dimension of the DoFHandler (given as the second template argument), then we are obviously dealing with cells, rather than lower-dimensional objects. In that case, inheritance is from CellAccessor, to provide access to all the cell specific information afforded by that class. Otherwise, i.e. for lower-dimensional objects, inheritance is from TriaAccessor.

There is a DoFCellAccessor class that provides the equivalent to the CellAccessor class.

Template Parameters
structdimThe dimensionality of the objects the accessor represents. For example, points have structdim equal to zero, edges have structdim equal to one, etc.
DoFHandlerTypeThe type of the DoF handler into which accessor of this type point. This is either the DoFHandler or hp::DoFHandler class. See also the DoFHandlerType concept.
level_dof_accessIf false, then the accessor simply represents a cell, face, or edge in a DoFHandler for which degrees of freedom only exist on the finest level. Some operations are not allowed in this case, such as asking for DoF indices on non-active cells. On the other hand, if this template argument is true, then the accessor represents an object in a multilevel hierarchy of degrees of freedom. In this case, accessing DoF indices of any cell is possible, and will return the level indices (which, for active cells, may be different from the global indices).
Author
Wolfgang Bangerth, 1998, 2006, 2008, Timo Heister, Guido Kanschat, 2012, 2013

Definition at line 209 of file dof_accessor.h.

Member Typedef Documentation

◆ BaseClass

template<int structdim, typename DoFHandlerType, bool level_dof_access>
using DoFAccessor< structdim, DoFHandlerType, level_dof_access >::BaseClass = typename ::internal::DoFAccessorImplementation:: Inheritance<structdim, dimension, space_dimension>::BaseClass

Declare an alias to the base class to make accessing some of the exception classes simpler.

Definition at line 233 of file dof_accessor.h.

◆ AccessorData

template<int structdim, typename DoFHandlerType, bool level_dof_access>
using DoFAccessor< structdim, DoFHandlerType, level_dof_access >::AccessorData = DoFHandlerType

Data type passed by the iterator class.

Definition at line 238 of file dof_accessor.h.

Constructor & Destructor Documentation

◆ DoFAccessor() [1/5]

template<int structdim, typename DoFHandlerType, bool level_dof_access>
DoFAccessor< structdim, DoFHandlerType, level_dof_access >::DoFAccessor ( )

Default constructor. Provides an accessor that can't be used.

◆ DoFAccessor() [2/5]

template<int structdim, typename DoFHandlerType, bool level_dof_access>
DoFAccessor< structdim, DoFHandlerType, level_dof_access >::DoFAccessor ( const Triangulation< DoFHandlerType::dimension, DoFHandlerType::space_dimension > *  tria,
const int  level,
const int  index,
const DoFHandlerType *  dof_handler 
)

Constructor that generates an access that points to a particular cell or face or edge in a DoFHandler or hp::DoFHandler.

Parameters
triaThe triangulation into which this accessor points.
levelThe level within the mesh hierarchy of the object pointed to. For example, coarse mesh cells will have level zero, their children level one, and so on. This argument is ignored for faces and edges which do not have a level.
indexThe index of the object pointed to within the specified refinement level.
dof_handlerA pointer to the DoFHandler or hp::DoFHandler object to which the accessor shall refer. This DoFHandler object must of course be built on the same triangulation as the one specified in the first argument.

◆ DoFAccessor() [3/5]

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<int structdim2, int dim2, int spacedim2>
DoFAccessor< structdim, DoFHandlerType, level_dof_access >::DoFAccessor ( const InvalidAccessor< structdim2, dim2, spacedim2 > &  )

Conversion constructor. This constructor exists to make certain constructs simpler to write in dimension independent code. For example, it allows assigning a face iterator to a line iterator, an operation that is useful in 2d but doesn't make any sense in 3d. The constructor here exists for the purpose of making the code conform to C++ but it will unconditionally abort; in other words, assigning a face iterator to a line iterator is better put into an if-statement that checks that the dimension is two, and assign to a quad iterator in 3d (an operator that, without this constructor would be illegal if we happen to compile for 2d).

◆ DoFAccessor() [4/5]

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<int dim2, class DoFHandlerType2 , bool level_dof_access2>
DoFAccessor< structdim, DoFHandlerType, level_dof_access >::DoFAccessor ( const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &  )

Another conversion operator between objects that don't make sense, just like the previous one.

◆ DoFAccessor() [5/5]

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<bool level_dof_access2>
DoFAccessor< structdim, DoFHandlerType, level_dof_access >::DoFAccessor ( const DoFAccessor< structdim, DoFHandlerType, level_dof_access2 > &  )

Copy constructor allowing to switch level access and active access.

Member Function Documentation

◆ operator=()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
DoFAccessor<structdim, DoFHandlerType, level_dof_access>& DoFAccessor< structdim, DoFHandlerType, level_dof_access >::operator= ( const DoFAccessor< structdim, DoFHandlerType, level_dof_access > &  da)
delete

Copy operator. These operators are usually used in a context like iterator a,b; *a=*b;. Presumably, the intent here is to copy the object pointed to by b to the object pointed to by a. However, the result of dereferencing an iterator is not an object but an accessor; consequently, this operation is not useful for iterators on DoF handler objects. Consequently, this operator is declared as deleted and can not be used.

◆ get_dof_handler()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
const DoFHandlerType& DoFAccessor< structdim, DoFHandlerType, level_dof_access >::get_dof_handler ( ) const

Return a handle on the DoFHandler object which we are using.

◆ copy_from() [1/2]

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<bool level_dof_access2>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::copy_from ( const DoFAccessor< structdim, DoFHandlerType, level_dof_access2 > &  a)

Implement the copy operator needed for the iterator classes.

◆ copy_from() [2/2]

template<int structdim, typename DoFHandlerType, bool level_dof_access>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::copy_from ( const TriaAccessorBase< structdim, DoFHandlerType::dimension, DoFHandlerType::space_dimension > &  da)

Copy operator used by the iterator class. Keeps the previously set dof handler, but sets the object coordinates of the TriaAccessor.

◆ is_level_cell()

template<int sd, typename DoFHandlerType , bool level_dof_access>
bool DoFAccessor< sd, DoFHandlerType, level_dof_access >::is_level_cell ( )
inlinestatic

Tell the caller whether get_active_or_mg_dof_indices() accesses active or level dofs.

Definition at line 2023 of file dof_accessor.h.

◆ child()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
TriaIterator<DoFAccessor<structdim, DoFHandlerType, level_dof_access> > DoFAccessor< structdim, DoFHandlerType, level_dof_access >::child ( const unsigned int  c) const

Return an iterator pointing to the c-th child.

◆ line()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
typename ::internal::DoFHandlerImplementation:: Iterators<DoFHandlerType, level_dof_access>::line_iterator DoFAccessor< structdim, DoFHandlerType, level_dof_access >::line ( const unsigned int  i) const

Pointer to the ith line bounding this object. If the current object is a line itself, then the only valid index is i equals to zero, and the function returns an iterator to itself.

◆ quad()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
typename ::internal::DoFHandlerImplementation:: Iterators<DoFHandlerType, level_dof_access>::quad_iterator DoFAccessor< structdim, DoFHandlerType, level_dof_access >::quad ( const unsigned int  i) const

Pointer to the ith quad bounding this object. If the current object is a quad itself, then the only valid index is i equals to zero, and the function returns an iterator to itself.

◆ get_dof_indices()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::get_dof_indices ( std::vector< types::global_dof_index > &  dof_indices,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
) const

Return the global indices of the degrees of freedom located on this object in the standard ordering defined by the finite element (i.e., dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc.) This function is only available on active objects (see this glossary entry).

The cells needs to be an active cell (and not artificial in a parallel distributed computation).

The vector has to have the right size before being passed to this function.

The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.

However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().

For cells, there is only a single possible finite element index (namely the one for that cell, returned by cell->active_fe_index. Consequently, the derived DoFCellAccessor class has an overloaded version of this function that calls the present function with cell->active_fe_index as last argument.

◆ get_mg_dof_indices()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::get_mg_dof_indices ( const int  level,
std::vector< types::global_dof_index > &  dof_indices,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
) const

Return the global multilevel indices of the degrees of freedom that live on the current object with respect to the given level within the multigrid hierarchy. The indices refer to the local numbering for the level this line lives on.

◆ set_mg_dof_indices()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::set_mg_dof_indices ( const int  level,
const std::vector< types::global_dof_index > &  dof_indices,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
)

Set the level DoF indices that are returned by get_mg_dof_indices.

◆ vertex_dof_index()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
types::global_dof_index DoFAccessor< structdim, DoFHandlerType, level_dof_access >::vertex_dof_index ( const unsigned int  vertex,
const unsigned int  i,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
) const

Global DoF index of the i degree associated with the vertexth vertex of the present cell.

The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.

However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().

◆ mg_vertex_dof_index()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
types::global_dof_index DoFAccessor< structdim, DoFHandlerType, level_dof_access >::mg_vertex_dof_index ( const int  level,
const unsigned int  vertex,
const unsigned int  i,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
) const

Return the global DoF index of the ith degree of freedom associated with the vertexth vertex on level level. Also see vertex_dof_index().

◆ dof_index()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
types::global_dof_index DoFAccessor< structdim, DoFHandlerType, level_dof_access >::dof_index ( const unsigned int  i,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
) const

Index of the ith degree of freedom of this object.

The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.

However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().

Note
While the get_dof_indices() function returns an array that contains the indices of all degrees of freedom that somehow live on this object (i.e. on the vertices, edges or interior of this object), the current dof_index() function only considers the DoFs that really belong to this particular object's interior. In other words, as an example, if the current object refers to a quad (a cell in 2d, a face in 3d) and the finite element associated with it is a bilinear one, then the get_dof_indices() will return an array of size 4 while dof_index() will produce an exception because no degrees are defined in the interior of the face.

◆ mg_dof_index()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
types::global_dof_index DoFAccessor< structdim, DoFHandlerType, level_dof_access >::mg_dof_index ( const int  level,
const unsigned int  i 
) const

Return the dof_index on the given level. Also see dof_index.

◆ n_active_fe_indices()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
unsigned int DoFAccessor< structdim, DoFHandlerType, level_dof_access >::n_active_fe_indices ( ) const

Return the number of finite elements that are active on a given object.

For non-hp DoFHandler objects, the answer is of course always one. However, for hp::DoFHandler objects, this isn't the case: If this is a cell, the answer is of course one. If it is a face, the answer may be one or two, depending on whether the two adjacent cells use the same finite element or not. If it is an edge in 3d, the possible return value may be one or any other value larger than that.

◆ nth_active_fe_index()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
unsigned int DoFAccessor< structdim, DoFHandlerType, level_dof_access >::nth_active_fe_index ( const unsigned int  n) const

Return the n-th active fe index on this object. For cells and all non- hp objects, there is only a single active fe index, so the argument must be equal to zero. For lower-dimensional hp objects, there are n_active_fe_indices() active finite elements, and this function can be queried for their indices.

◆ get_active_fe_indices()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
std::set<unsigned int> DoFAccessor< structdim, DoFHandlerType, level_dof_access >::get_active_fe_indices ( ) const

Returns all active fe indices on this object.

The size of the returned set equals the number of finite elements that are active on this object.

◆ fe_index_is_active()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
bool DoFAccessor< structdim, DoFHandlerType, level_dof_access >::fe_index_is_active ( const unsigned int  fe_index) const

Return true if the finite element with given index is active on the present object. For non-hp DoF accessors, this is of course the case only if fe_index equals zero. For cells, it is the case if fe_index equals active_fe_index() of this cell. For faces and other lower- dimensional objects, there may be more than one fe_index that are active on any given object (see n_active_fe_indices()).

◆ get_fe()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
const FiniteElement<DoFHandlerType::dimension, DoFHandlerType::space_dimension>& DoFAccessor< structdim, DoFHandlerType, level_dof_access >::get_fe ( const unsigned int  fe_index) const

Return a reference to the finite element used on this object with the given fe_index. fe_index must be used on this object, i.e. fe_index_is_active(fe_index) must return true.

◆ operator==()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<int dim2, class DoFHandlerType2 , bool level_dof_access2>
bool DoFAccessor< structdim, DoFHandlerType, level_dof_access >::operator== ( const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &  ) const

Compare for equality. Return true if the two accessors refer to the same object.

The template parameters of this function allow for a comparison of very different objects. Therefore, some of them are disabled. Namely, if the dimension, or the dof handler of the two objects differ, an exception is generated. It can be expected that this is an unwanted comparison.

The template parameter level_dof_access2 is ignored, such that an iterator with level access can be equal to one with access to the active degrees of freedom.

◆ operator!=()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<int dim2, class DoFHandlerType2 , bool level_dof_access2>
bool DoFAccessor< structdim, DoFHandlerType, level_dof_access >::operator!= ( const DoFAccessor< dim2, DoFHandlerType2, level_dof_access2 > &  ) const

Compare for inequality. The boolean not of operator==().

◆ set_dof_handler()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::set_dof_handler ( DoFHandlerType *  dh)
protected

Reset the DoF handler pointer.

◆ set_dof_index()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::set_dof_index ( const unsigned int  i,
const types::global_dof_index  index,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
) const
protected

Set the index of the ith degree of freedom of this object to index.

The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.

However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().

◆ set_vertex_dof_index()

template<int structdim, typename DoFHandlerType, bool level_dof_access>
void DoFAccessor< structdim, DoFHandlerType, level_dof_access >::set_vertex_dof_index ( const unsigned int  vertex,
const unsigned int  i,
const types::global_dof_index  index,
const unsigned int  fe_index = DoFHandlerType::default_fe_index 
) const
protected

Set the global index of the i degree on the vertex-th vertex of the present cell to index.

The last argument denotes the finite element index. For the standard DoFHandler class, this value must be equal to its default value since that class only supports the same finite element on all cells anyway.

However, for hp objects (i.e. the hp::DoFHandler class), different finite element objects may be used on different cells. On faces between two cells, as well as vertices, there may therefore be two sets of degrees of freedom, one for each of the finite elements used on the adjacent cells. In order to specify which set of degrees of freedom to work on, the last argument is used to disambiguate. Finally, if this function is called for a cell object, there can only be a single set of degrees of freedom, and fe_index has to match the result of active_fe_index().

Friends And Related Function Documentation

◆ TriaRawIterator

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<typename >
friend class TriaRawIterator
friend

Iterator classes need to be friends because they need to access operator== and operator!=.

Definition at line 737 of file dof_accessor.h.

◆ DoFHandler

template<int structdim, typename DoFHandlerType, bool level_dof_access>
template<int dim, int spacedim>
friend class DoFHandler
friend

Make the DoFHandler class a friend so that it can call the set_xxx() functions.

Definition at line 747 of file dof_accessor.h.

Member Data Documentation

◆ dimension

template<int structdim, typename DoFHandlerType, bool level_dof_access>
const unsigned int DoFAccessor< structdim, DoFHandlerType, level_dof_access >::dimension = DoFHandlerType::dimension
static

A static variable that allows users of this class to discover the value of the second template argument.

Definition at line 220 of file dof_accessor.h.

◆ space_dimension

template<int structdim, typename DoFHandlerType, bool level_dof_access>
const unsigned int DoFAccessor< structdim, DoFHandlerType, level_dof_access >::space_dimension = DoFHandlerType::space_dimension
static

A static variable that allows users of this class to discover the value of the third template argument.

Definition at line 226 of file dof_accessor.h.

◆ dof_handler

template<int structdim, typename DoFHandlerType, bool level_dof_access>
DoFHandlerType* DoFAccessor< structdim, DoFHandlerType, level_dof_access >::dof_handler
protected

Store the address of the DoFHandler object to be accessed.

Definition at line 636 of file dof_accessor.h.


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