Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > Class Template Reference

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

Inheritance diagram for DoFCellAccessor< dimension_, space_dimension_, level_dof_access >:
[legend]

Public Types

using AccessorData = DoFHandler< dimension_, space_dimension_ >
 
using BaseClass = DoFAccessor< dimension_, dimension_, space_dimension_, level_dof_access >
 
using Container = DoFHandler< dimension_, space_dimension_ >
 
using face_iterator = TriaIterator< DoFAccessor< dimension_ - 1, dimension_, space_dimension_, level_dof_access > >
 

Public Member Functions

TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent () const
 
void set_dof_indices (const std::vector< types::global_dof_index > &dof_indices)
 
void set_mg_dof_indices (const std::vector< types::global_dof_index > &dof_indices)
 
const DoFHandler< dim, spacedim > & get_dof_handler () const
 
void copy_from (const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
 
void copy_from (const TriaAccessorBase< structdim, dim, spacedim > &da)
 
bool operator== (const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
 
bool operator!= (const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
 
Constructors and initialization
 DoFCellAccessor (const Triangulation< dimension_, space_dimension_ > *tria, const int level, const int index, const AccessorData *local_data)
 
template<int structdim2, int dim2, int spacedim2>
 DoFCellAccessor (const InvalidAccessor< structdim2, dim2, spacedim2 > &)
 
template<int structdim2, int dim2, int spacedim2, bool level_dof_access2>
 DoFCellAccessor (const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
 
 DoFCellAccessor (const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &)=default
 
 DoFCellAccessor (DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
 
 ~DoFCellAccessor ()=default
 
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator= (const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &da)=delete
 
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator= (DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
 
Accessing sub-objects and neighbors
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor (const unsigned int i) const
 
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor (const unsigned int i) const
 
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_or_periodic_neighbor (const unsigned int i) const
 
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child (const unsigned int i) const
 
boost::container::small_vector< TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > >, GeometryInfo< dimension_ >::max_children_per_cell > child_iterators () const
 
face_iterator face (const unsigned int i) const
 
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators () const
 
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_child_on_subface (const unsigned int face_no, const unsigned int subface_no) const
 
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor_child_on_subface (const unsigned int face_no, const unsigned int subface_no) const
 
Extracting values from global vectors
template<class InputVector , typename number >
void get_dof_values (const InputVector &values, Vector< number > &local_values) const
 
template<class InputVector , typename ForwardIterator >
void get_dof_values (const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
 
template<class InputVector , typename ForwardIterator >
void get_dof_values (const AffineConstraints< typename InputVector::value_type > &constraints, const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
 
template<class OutputVector , typename number >
void set_dof_values (const Vector< number > &local_values, OutputVector &values) const
 
template<class InputVector , typename number >
void get_interpolated_dof_values (const InputVector &values, Vector< number > &interpolated_values, const types::fe_index fe_index=numbers::invalid_fe_index) const
 
template<class OutputVector , typename number >
void set_dof_values_by_interpolation (const Vector< number > &local_values, OutputVector &values, const types::fe_index fe_index=numbers::invalid_fe_index, const bool perform_check=false) const
 
template<class OutputVector , typename number >
void distribute_local_to_global_by_interpolation (const Vector< number > &local_values, OutputVector &values, const types::fe_index fe_index=numbers::invalid_fe_index) const
 
template<typename number , typename OutputVector >
void distribute_local_to_global (const Vector< number > &local_source, OutputVector &global_destination) const
 
template<typename ForwardIterator , typename OutputVector >
void distribute_local_to_global (ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
 
template<typename ForwardIterator , typename OutputVector >
void distribute_local_to_global (const AffineConstraints< typename OutputVector::value_type > &constraints, ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
 
template<typename number , typename OutputMatrix >
void distribute_local_to_global (const FullMatrix< number > &local_source, OutputMatrix &global_destination) const
 
template<typename number , typename OutputMatrix , typename OutputVector >
void distribute_local_to_global (const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, OutputMatrix &global_matrix, OutputVector &global_vector) const
 
Accessing the DoF indices of this object
void get_active_or_mg_dof_indices (std::vector< types::global_dof_index > &dof_indices) const
 
void get_dof_indices (std::vector< types::global_dof_index > &dof_indices) const
 
void get_mg_dof_indices (std::vector< types::global_dof_index > &dof_indices) const
 
Accessing the finite element associated with this object
const FiniteElement< dimension_, space_dimension_ > & get_fe () const
 
types::fe_index active_fe_index () const
 
void set_active_fe_index (const types::fe_index i) const
 
Dealing with refinement indicators
const FiniteElement< dimension_, space_dimension_ > & get_future_fe () const
 
types::fe_index future_fe_index () const
 
void set_future_fe_index (const types::fe_index i) const
 
bool future_fe_index_set () const
 
void clear_future_fe_index () const
 
Accessing sub-objects
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator line (const unsigned int i) const
 
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, 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 types::fe_index fe_index=numbers::invalid_fe_index) const
 
void get_mg_dof_indices (const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
 
void set_mg_dof_indices (const int level, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index)
 
types::global_dof_index vertex_dof_index (const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
 
types::global_dof_index mg_vertex_dof_index (const int level, const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
 
types::global_dof_index dof_index (const unsigned int i, const types::fe_index fe_index=numbers::invalid_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
 
types::fe_index nth_active_fe_index (const unsigned int n) const
 
std::set< types::fe_indexget_active_fe_indices () const
 
bool fe_index_is_active (const types::fe_index fe_index) const
 
const FiniteElement< dim, spacedim > & get_fe (const types::fe_index fe_index) 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 dim = dimension_
 
static const unsigned int spacedim = space_dimension_
 
static constexpr unsigned int dimension
 
static constexpr unsigned int space_dimension
 

Protected Member Functions

void set_dof_handler (DoFHandler< dim, spacedim > *dh)
 
void set_dof_index (const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
 
void set_mg_dof_index (const int level, const unsigned int i, const types::global_dof_index index) const
 
void set_mg_vertex_dof_index (const int level, const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
 

Protected Attributes

DoFHandler< dim, spacedim > * dof_handler
 

Friends

struct ::internal::DoFCellAccessorImplementation::Implementation
 

Detailed Description

template<int dimension_, int space_dimension_, bool level_dof_access>
class DoFCellAccessor< dimension_, space_dimension_, level_dof_access >

Grant access to the degrees of freedom on a cell.

Note that since for the class we derive from, i.e. DoFAccessor<dim>, the two template parameters are equal, the base class is actually derived from CellAccessor, which makes the functions of this class available to the DoFCellAccessor class as well.

Definition at line 1321 of file dof_accessor.h.

Member Typedef Documentation

◆ AccessorData

template<int dimension_, int space_dimension_, bool level_dof_access>
using DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::AccessorData = DoFHandler<dimension_, space_dimension_>

Data type passed by the iterator class.

Definition at line 1341 of file dof_accessor.h.

◆ BaseClass

template<int dimension_, int space_dimension_, bool level_dof_access>
using DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::BaseClass = DoFAccessor<dimension_, dimension_, space_dimension_, level_dof_access>

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

Definition at line 1347 of file dof_accessor.h.

◆ Container

template<int dimension_, int space_dimension_, bool level_dof_access>
using DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::Container = DoFHandler<dimension_, space_dimension_>

Define the type of the container this is part of.

Definition at line 1353 of file dof_accessor.h.

◆ face_iterator

template<int dimension_, int space_dimension_, bool level_dof_access>
using DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::face_iterator = TriaIterator<DoFAccessor<dimension_ - 1, dimension_, space_dimension_, level_dof_access> >

A type for an iterator over the faces of a cell. This is what the face() function returns.

Definition at line 1359 of file dof_accessor.h.

Constructor & Destructor Documentation

◆ DoFCellAccessor() [1/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::DoFCellAccessor ( const Triangulation< dimension_, space_dimension_ > *  tria,
const int  level,
const int  index,
const AccessorData local_data 
)

Constructor

◆ DoFCellAccessor() [2/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<int structdim2, int dim2, int spacedim2>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::DoFCellAccessor ( 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).

◆ DoFCellAccessor() [3/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<int structdim2, int dim2, int spacedim2, bool level_dof_access2>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::DoFCellAccessor ( const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &  )
explicit

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

◆ DoFCellAccessor() [4/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::DoFCellAccessor ( const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &  )
default

Copy constructor.

◆ DoFCellAccessor() [5/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::DoFCellAccessor ( DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&  )
default

Move constructor.

◆ ~DoFCellAccessor()

template<int dimension_, int space_dimension_, bool level_dof_access>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::~DoFCellAccessor ( )
default

Destructor

Member Function Documentation

◆ operator=() [1/2]

template<int dimension_, int space_dimension_, bool level_dof_access>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::operator= ( const DoFCellAccessor< dimension_, space_dimension_, 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.

◆ operator=() [2/2]

template<int dimension_, int space_dimension_, bool level_dof_access>
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::operator= ( DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&  )
default

Move assignment operator.

◆ parent()

template<int dimension_, int space_dimension_, bool level_dof_access>
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::parent ( ) const

Return the parent of this cell as a DoF cell iterator. If the parent does not exist (i.e., if the object is at the coarsest level of the mesh hierarchy), an exception is generated.

This function is needed since the parent function of the base class CellAccessor returns a triangulation cell accessor without access to the DoF data.

◆ neighbor()

template<int dimension_, int space_dimension_, bool level_dof_access>
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::neighbor ( const unsigned int  i) const

Return the ith neighbor as a DoF cell iterator. This function is needed since the neighbor function of the base class returns a cell accessor without access to the DoF data.

◆ periodic_neighbor()

template<int dim, int spacedim, bool lda>
TriaIterator< DoFCellAccessor< dim, spacedim, lda > > DoFCellAccessor< dim, spacedim, lda >::periodic_neighbor ( const unsigned int  i) const
inline

Return the ith periodic neighbor as a DoF cell iterator. This function is needed since the neighbor function of the base class returns a cell accessor without access to the DoF data.

Definition at line 154 of file dof_accessor.cc.

◆ neighbor_or_periodic_neighbor()

template<int dim, int spacedim, bool lda>
TriaIterator< DoFCellAccessor< dim, spacedim, lda > > DoFCellAccessor< dim, spacedim, lda >::neighbor_or_periodic_neighbor ( const unsigned int  i) const
inline

Return the ith neighbor or periodic neighbor as a DoF cell iterator. This function is needed since the neighbor function of the base class returns a cell accessor without access to the DoF data.

Definition at line 167 of file dof_accessor.cc.

◆ child()

template<int dimension_, int space_dimension_, bool level_dof_access>
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::child ( const unsigned int  i) const

Return the ith child as a DoF cell iterator. This function is needed since the child function of the base class returns a cell accessor without access to the DoF data.

◆ child_iterators()

template<int dimension_, int space_dimension_, bool level_dof_access>
boost::container::small_vector< TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > >, GeometryInfo< dimension_ >::max_children_per_cell > DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::child_iterators ( ) const

Return an array of iterators to all children of this cell.

◆ face()

template<int dimension_, int space_dimension_, bool level_dof_access>
face_iterator DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::face ( const unsigned int  i) const

Return an iterator to the ith face of this cell.

This function returns a DoFAccessor with structdim == 0 in 1d, a DoFAccessor::line in 2d, and a DoFAccessor::quad in 3d.

◆ face_iterators()

template<int dimension_, int space_dimension_, bool level_dof_access>
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::face_iterators ( ) const

Return an array of iterators to all faces of this cell.

◆ neighbor_child_on_subface()

template<int dim, int spacedim, bool lda>
TriaIterator< DoFCellAccessor< dim, spacedim, lda > > DoFCellAccessor< dim, spacedim, lda >::neighbor_child_on_subface ( const unsigned int  face_no,
const unsigned int  subface_no 
) const

Return the result of the neighbor_child_on_subface function of the base class, but convert it so that one can also access the DoF data (the function in the base class only returns an iterator with access to the triangulation data).

Definition at line 125 of file dof_accessor.cc.

◆ periodic_neighbor_child_on_subface()

template<int dim, int spacedim, bool lda>
TriaIterator< DoFCellAccessor< dim, spacedim, lda > > DoFCellAccessor< dim, spacedim, lda >::periodic_neighbor_child_on_subface ( const unsigned int  face_no,
const unsigned int  subface_no 
) const

Return the result of the periodic_neighbor_child_on_subface function of the base class, but convert it so that one can also access the DoF data (the function in the base class only returns an iterator with access to the triangulation data).

Definition at line 139 of file dof_accessor.cc.

◆ get_dof_values() [1/3]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<class InputVector , typename number >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_dof_values ( const InputVector &  values,
Vector< number > &  local_values 
) const

Collect the values of the given vector restricted to the dofs of this cell in the standard ordering: dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc. In other words, this function implements a gather operation.

The vector has to have the right size before being passed to this function. This function is only callable for active cells.

The input vector may be either a Vector<float>, Vector<double>, or a BlockVector<double>, or a PETSc or Trilinos vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy.

◆ get_dof_values() [2/3]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<class InputVector , typename ForwardIterator >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_dof_values ( const InputVector &  values,
ForwardIterator  local_values_begin,
ForwardIterator  local_values_end 
) const

Collect the values of the given vector restricted to the dofs of this cell in the standard ordering: dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc. In other words, this function implements a gather operation.

The vector has to have the right size before being passed to this function. This function is only callable for active cells.

The input vector may be either a Vector<float>, Vector<double>, or a BlockVector<double>, or a PETSc or Trilinos vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy.

◆ get_dof_values() [3/3]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<class InputVector , typename ForwardIterator >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_dof_values ( const AffineConstraints< typename InputVector::value_type > &  constraints,
const InputVector &  values,
ForwardIterator  local_values_begin,
ForwardIterator  local_values_end 
) const

Collect the values of the given vector restricted to the dofs of this cell in the standard ordering: dofs on vertex 0, dofs on vertex 1, etc, dofs on line 0, dofs on line 1, etc, dofs on quad 0, etc. In other words, this function implements a gather operation.

The vector has to have the right size before being passed to this function. This function is only callable for active cells.

The input vector may be either a Vector<float>, Vector<double>, or a BlockVector<double>, or a PETSc or Trilinos vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy. The AffineConstraints object passed as an argument to this function makes sure that constraints are correctly distributed when the dof values are calculated.

◆ set_dof_values()

template<int dimension_, int space_dimension_, bool level_dof_access>
template<class OutputVector , typename number >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::set_dof_values ( const Vector< number > &  local_values,
OutputVector &  values 
) const

This function is the counterpart to get_dof_values(): it takes a vector of values for the degrees of freedom of the cell pointed to by this iterator and writes these values into the global data vector values. In other words, this function implements a scatter operation. This function is only callable for active cells.

Note that for continuous finite elements, calling this function affects the dof values on neighboring cells as well. It may also violate continuity requirements for hanging nodes, if neighboring cells are less refined than the present one. These requirements are not taken care of and must be enforced by the user afterwards.

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

The output vector may be either a Vector<float>, Vector<double>, or a BlockVector<double>, or a PETSc vector if deal.II is compiled to support these libraries. It is in the responsibility of the caller to assure that the types of the numbers stored in input and output vectors are compatible and with similar accuracy.

◆ get_interpolated_dof_values()

template<int dim, int spacedim, bool lda>
template<class InputVector , typename number >
void DoFCellAccessor< dim, spacedim, lda >::get_interpolated_dof_values ( const InputVector &  values,
Vector< number > &  interpolated_values,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const

Return the interpolation of the given finite element function to the present cell. In the simplest case, the cell is a terminal one, i.e., it has no children; then, the returned value is the vector of nodal values on that cell. You could as well get the desired values through the get_dof_values function. In the other case, when the cell has children, we use the restriction matrices provided by the finite element class to compute the interpolation from the children to the present cell.

If the cell is part of a DoFHandler with hp-capabilities, cells only have an associated finite element space if they are active. However, this function is supposed to also provide information on inactive cells with children. Consequently, it carries a third argument that can be used in the hp-context that denotes the finite element space we are supposed to interpolate onto. If the cell is active, this function then obtains the finite element function from the values vector on this cell and interpolates it onto the space described by the fe_indexth element of the hp::FECollection associated with the DoFHandler of which this cell is a part of. If the cell is not active, then we first perform this interpolation on all of its terminal children and then interpolate this function down to the cell requested keeping the function space the same.

It is assumed that both input vectors already have the right size beforehand.

Note
Unlike the get_dof_values() function, this function is only available on cells, rather than on lines, quads, and hexes, since interpolation is presently only provided for cells by the finite element classes.

Definition at line 46 of file dof_accessor_get.cc.

◆ set_dof_values_by_interpolation()

template<int dim, int spacedim, bool lda>
template<class OutputVector , typename number >
void DoFCellAccessor< dim, spacedim, lda >::set_dof_values_by_interpolation ( const Vector< number > &  local_values,
OutputVector &  values,
const types::fe_index  fe_index = numbers::invalid_fe_index,
const bool  perform_check = false 
) const

This function is the counterpart to get_interpolated_dof_values(): you specify the dof values on a cell and these are interpolated to the children of the present cell and set on the terminal cells.

In principle, it works as follows: if the cell pointed to by this object is terminal (i.e., has no children), then the dof values are set in the global data vector by calling the set_dof_values() function; otherwise, the values are prolonged to each of the children and this function is called for each of them.

Using the get_interpolated_dof_values() and this function, you can compute the interpolation of a finite element function to a coarser grid by first getting the interpolated solution on a cell of the coarse grid and afterwards redistributing it using this function.

Note that for continuous finite elements, calling this function affects the dof values on neighboring cells as well. It may also violate continuity requirements for hanging nodes, if neighboring cells are less refined than the present one, or if their children are less refined than the children of this cell. These requirements are not taken care of and must be enforced by the user afterward.

If the cell is part of a DoFHandler with hp-capabilities, cells only have an associated finite element space if they are active. However, this function is supposed to also work on inactive cells with children. Consequently, it carries a third argument that can be used in the hp- context that denotes the finite element space we are supposed to interpret the input vector of this function in. If the cell is active, this function then interpolates the input vector interpreted as an element of the space described by the fe_indexth element of the hp::FECollection associated with the DoFHandler of which this cell is a part of, and interpolates it into the space that is associated with this cell. On the other hand, if the cell is not active, then we first perform this interpolation from this cell to its children using the given fe_index until we end up on an active cell, at which point we follow the procedure outlined at the beginning of the paragraph.

It is assumed that both vectors already have the right size beforehand. This function relies on the existence of a natural interpolation property of finite element spaces of a cell to its children, denoted by the prolongation matrices of finite element classes. For some elements, the spaces on coarse and fine grids are not nested, in which case the interpolation to a child is not the identity; refer to the documentation of the respective finite element class for a description of what the prolongation matrices represent in this case.

Note
Cells set the values of DoFs independently and might overwrite previously set values in the global vector, for example when calling the same function earlier from a different cell. By setting perform_check, you can enable a check that the previous value and the one to be set here are at least roughly the same. In practice, they might be slightly different because they are computed in a way that theoretically ensures that they are the same, but in practice they are only equal up to round-off.
Unlike the get_dof_values() function, this function is only available on cells, rather than on lines, quads, and hexes, since interpolation is presently only provided for cells by the finite element classes.

Definition at line 272 of file dof_accessor_set.cc.

◆ distribute_local_to_global_by_interpolation()

template<int dim, int spacedim, bool lda>
template<class OutputVector , typename number >
void DoFCellAccessor< dim, spacedim, lda >::distribute_local_to_global_by_interpolation ( const Vector< number > &  local_values,
OutputVector &  values,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const

Similar to set_dof_values_by_interpolation() with the difference that values are added into the vector.

Note
In parallel::distributed::SolutionTransfer, this function is used to accumulate the contributions of all cells to a DoF; with a subsequent multiplication with the inverse of the valence, finally, the average value is obtained.

Definition at line 294 of file dof_accessor_set.cc.

◆ distribute_local_to_global() [1/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<typename number , typename OutputVector >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::distribute_local_to_global ( const Vector< number > &  local_source,
OutputVector &  global_destination 
) const

Distribute a local (cell based) vector to a global one by mapping the local numbering of the degrees of freedom to the global one and entering the local values into the global vector. In other words, this function implements a scatter operation.

The elements are added to the existing elements in the global vector, rather than just set, since this is usually what one wants. You may also want to take a look at the AffineConstraints::distribute_local_to_global() function if you need to deal with constraints.

◆ distribute_local_to_global() [2/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<typename ForwardIterator , typename OutputVector >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::distribute_local_to_global ( ForwardIterator  local_source_begin,
ForwardIterator  local_source_end,
OutputVector &  global_destination 
) const

Distribute a local (cell based) vector in iterator format to a global one by mapping the local numbering of the degrees of freedom to the global one and entering the local values into the global vector. In other words, this function implements a scatter operation.

The elements are added to the existing elements in the global vector, rather than just set, since this is usually what one wants. You may also want to take a look at the AffineConstraints::distribute_local_to_global() function if you need to deal with constraints.

◆ distribute_local_to_global() [3/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<typename ForwardIterator , typename OutputVector >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::distribute_local_to_global ( const AffineConstraints< typename OutputVector::value_type > &  constraints,
ForwardIterator  local_source_begin,
ForwardIterator  local_source_end,
OutputVector &  global_destination 
) const

Distribute a local (cell based) vector in iterator format to a global one by mapping the local numbering of the degrees of freedom to the global one and entering the local values into the global vector. In other words, this function implements a scatter operation.

The elements are added up to the elements in the global vector, rather than just set, since this is usually what one wants. Moreover, the AffineConstraints object passed to this function makes sure that also constraints are eliminated in this process.

◆ distribute_local_to_global() [4/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<typename number , typename OutputMatrix >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::distribute_local_to_global ( const FullMatrix< number > &  local_source,
OutputMatrix &  global_destination 
) const

This function does much the same as the distribute_local_to_global(Vector,Vector) function, but operates on matrices instead of vectors. If the matrix type is a sparse matrix then it is supposed to have non-zero entry slots where required.

◆ distribute_local_to_global() [5/5]

template<int dimension_, int space_dimension_, bool level_dof_access>
template<typename number , typename OutputMatrix , typename OutputVector >
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::distribute_local_to_global ( const FullMatrix< number > &  local_matrix,
const Vector< number > &  local_vector,
OutputMatrix &  global_matrix,
OutputVector &  global_vector 
) const

This function does what the two distribute_local_to_global functions with vector and matrix argument do, but all at once.

◆ get_active_or_mg_dof_indices()

template<int dimension_, int space_dimension_, bool level_dof_access>
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_active_or_mg_dof_indices ( std::vector< types::global_dof_index > &  dof_indices) const

Obtain the global indices of the local degrees of freedom on this cell.

If this object accesses a level cell (indicated by the third template argument or is_level_cell), then return the result of get_mg_dof_indices(), else return get_dof_indices().

You will get a level_cell_iterator when calling begin_mg() and a normal one otherwise.

Examples for this use are in the implementation of DoFRenumbering.

◆ get_dof_indices() [1/2]

template<int dimension_, int space_dimension_, bool level_dof_access>
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_dof_indices ( std::vector< types::global_dof_index > &  dof_indices) const
inline

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).

Parameters
[out]dof_indicesThe vector into which the indices will be written. It has to have the right size (namely, fe.n_dofs_per_cell(), fe.dofs_per_face, or fe.dofs_per_line, depending on which kind of object this function is called) before being passed to this function.

This function reimplements the same function in the base class. In contrast to the function in the base class, we do not need the fe_index here because there is always a unique finite element index on cells.

This is a function which requires that the cell is active.

Also see get_active_or_mg_dof_indices().

Note
In many places in the tutorial and elsewhere in the library, the argument to this function is called local_dof_indices by convention. The name is not meant to indicate the local numbers of degrees of freedom (which are always between zero and fe.n_dofs_per_cell()) but instead that the returned values are the global indices of those degrees of freedom that are located locally on the current cell.

Definition at line 282 of file dof_accessor.cc.

◆ get_mg_dof_indices() [1/2]

template<int dimension_, int space_dimension_, bool level_dof_access>
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_mg_dof_indices ( std::vector< types::global_dof_index > &  dof_indices) const
inline

Retrieve the global indices of the degrees of freedom on this cell in the level vector associated to the level of the cell.

Definition at line 299 of file dof_accessor.cc.

◆ get_fe() [1/2]

template<int dimension_, int space_dimension_, bool level_dof_access>
const FiniteElement< dimension_, space_dimension_ > & DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_fe ( ) const

Return the finite element that is used on the cell pointed to by this iterator. For DoFHandler objects without hp-capabilities, this is of course always the same element, independent of the cell we are presently on, but for hp-DoFHandler objects this may change from cell to cell.

Note
Since degrees of freedom only exist on active cells for DoFHandler objects with hp-capabilities (i.e., there is currently no implementation of multilevel such objects), it does not make sense to query the finite element on non-active cells since they do not have finite element spaces associated with them without having any degrees of freedom. Consequently, this function will produce an exception when called on non-active cells.

◆ active_fe_index()

template<int dimension_, int space_dimension_, bool level_dof_access>
types::fe_index DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::active_fe_index ( ) const

Return the index inside the hp::FECollection of the FiniteElement used for this cell. This function is only useful if the DoFHandler object associated with the current cell has hp-capabilities enabled.

Note
Since degrees of freedom only exist on active cells for DoFHandler objects with hp-capabilities (i.e., there is currently no implementation of multilevel such objects), it does not make sense to query the finite element on non-active cells since they do not have finite element spaces associated with them without having any degrees of freedom. Consequently, this function will produce an exception when called on non-active cells.
When using parallel meshes, either through the parallel::shared::Triangulation or parallel::distributed::Triangulation classes, it is only allowed to call this function on locally owned or ghost cells. No information is available on artificial cells. Furthermore, active_fe_index information is only exchanged from locally owned cells on one processor to other processors where they may be ghost cells, during the call to DoFHandler::set_fe() and DoFHandler::distribute_dofs(). Be aware that if you call set_active_fe_index() on a cell after calling one of these functions, then this information will not be propagated to other processors who may have this cell as a ghost cell. See the documentation of DoFHandler for more information.

◆ set_active_fe_index()

template<int dimension_, int space_dimension_, bool level_dof_access>
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::set_active_fe_index ( const types::fe_index  i) const

Set the index of the FiniteElement used for this cell. This determines which element in an hp::FECollection to use. This function is only useful if the DoF handler object associated with the current cell has hp- capabilities enabled.

Note
Since degrees of freedom only exist on active cells for DoFHandler objects with hp-capabilities (i.e., there is currently no implementation of multilevel such objects), it does not make sense to query the finite element on non-active cells since they do not have finite element spaces associated with them without having any degrees of freedom. Consequently, this function will produce an exception when called on non-active cells.
When using parallel meshes, either through the parallel::shared::Triangulation or parallel::distributed::Triangulation classes, it is only allowed to call this function on locally owned or ghost cells. No information is available on artificial cells. Furthermore, active_fe_index information is only exchanged from locally owned cells on one processor to other processors where they may be ghost cells, during the call to DoFHandler::set_fe() and DoFHandler::distribute_dofs(). Be aware that if you call set_active_fe_index() on a cell after calling one of these functions, then this information will not be propagated to other processors who may have this cell as a ghost cell. See the documentation of DoFHandler for more information.

◆ set_dof_indices()

template<int dim, int spacedim, bool lda>
void DoFCellAccessor< dim, spacedim, lda >::set_dof_indices ( const std::vector< types::global_dof_index > &  dof_indices)

Set the DoF indices of this cell to the given values. This function bypasses the DoF cache, if one exists for the given DoF handler class.

Definition at line 106 of file dof_accessor.cc.

◆ set_mg_dof_indices() [1/2]

template<int dimension_, int space_dimension_, bool level_dof_access>
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::set_mg_dof_indices ( const std::vector< types::global_dof_index > &  dof_indices)
inline

Set the Level DoF indices of this cell to the given values.

Definition at line 313 of file dof_accessor.cc.

◆ get_future_fe()

template<int dimension_, int space_dimension_, bool level_dof_access>
const FiniteElement< dimension_, space_dimension_ > & DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::get_future_fe ( ) const

Return the finite element that will be assigned to this cell next time the triangulation gets refined and coarsened. If no future finite element has been specified for this cell via the set_future_fe_index() function, the active one will remain unchanged, in which case the active finite element will be returned.

For DoFHandlers without hp-capabilities enabled, this is of course always the same element, independent of the cell we are presently on, but for hp- DoFHandler objects this may change from cell to cell.

Note
Since degrees of freedom only exist on active cells for DoFHandler objects with hp-capabilities (i.e., there is currently no implementation of multilevel such objects), it does not make sense to query the finite element on non-active cells since they do not have finite element spaces associated with them without having any degrees of freedom. Consequently, this function will produce an exception when called on non-active cells.

◆ future_fe_index()

template<int dimension_, int space_dimension_, bool level_dof_access>
types::fe_index DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::future_fe_index ( ) const

Return the fe_index of the finite element that will be assigned to this cell next time the triangulation gets refined and coarsened. If no future finite element has been specified for this cell via the set_future_fe_index() function, the active one will remain unchanged, in which case the fe_index of the active finite element will be returned.

Note
Since degrees of freedom only exist on active cells for DoFHandler objects with hp-capabilities (i.e., there is currently no implementation of multilevel such objects), it does not make sense to query the finite element on non-active cells since they do not have finite element spaces associated with them without having any degrees of freedom. Consequently, this function will produce an exception when called on non-active cells.
When using parallel meshes, either through the parallel::shared::Triangulation or parallel::distributed::Triangulation classes, it is only allowed to call this function on locally owned cells.

◆ set_future_fe_index()

template<int dimension_, int space_dimension_, bool level_dof_access>
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::set_future_fe_index ( const types::fe_index  i) const

Set the fe_index of the finite element that will be assigned to this cell next time the triangulation gets refined and coarsened. A previously assigned future finite element will be overwritten.

See notes of future_fe_index() for information about restrictions on this functionality.

◆ future_fe_index_set()

template<int dimension_, int space_dimension_, bool level_dof_access>
bool DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::future_fe_index_set ( ) const

Return whether a future finite element has been set.

See notes of future_fe_index() for information about restrictions on this functionality.

◆ clear_future_fe_index()

template<int dimension_, int space_dimension_, bool level_dof_access>
void DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::clear_future_fe_index ( ) const

Revoke the future finite element assigned. Thus, the active finite element will remain unchanged next time the triangulation gets refined and coarsened.

See notes on future_fe_index() for information about restrictions on this functionality.

◆ get_dof_handler()

const DoFHandler< dim, spacedim > & DoFAccessor< structdim, dim, spacedim, level_dof_access >::get_dof_handler ( ) const
inherited

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

◆ copy_from() [1/2]

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::copy_from ( const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &  a)
inherited

Implement the copy operator needed for the iterator classes.

◆ copy_from() [2/2]

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::copy_from ( const TriaAccessorBase< structdim, dim, spacedim > &  da)
inherited

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

◆ is_level_cell()

bool DoFAccessor< structdim, dim, spacedim, level_dof_access >::is_level_cell
inlinestaticinherited

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

Definition at line 366 of file dof_accessor.h.

◆ line()

typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator DoFAccessor< structdim, dim, spacedim, level_dof_access >::line ( const unsigned int  i) const
inherited

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()

typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator DoFAccessor< structdim, dim, spacedim, level_dof_access >::quad ( const unsigned int  i) const
inherited

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() [2/2]

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::get_dof_indices ( std::vector< types::global_dof_index > &  dof_indices,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const
inlineinherited

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, when the relevant DoFHandler object has hp-capabilities enabled, 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.

Definition at line 445 of file dof_accessor.cc.

◆ get_mg_dof_indices() [2/2]

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::get_mg_dof_indices ( const int  level,
std::vector< types::global_dof_index > &  dof_indices,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const
inlineinherited

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.

Definition at line 456 of file dof_accessor.cc.

◆ set_mg_dof_indices() [2/2]

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::set_mg_dof_indices ( const int  level,
const std::vector< types::global_dof_index > &  dof_indices,
const types::fe_index  fe_index = numbers::invalid_fe_index 
)
inlineinherited

Set the level DoF indices that are returned by get_mg_dof_indices.

Definition at line 465 of file dof_accessor.cc.

◆ vertex_dof_index()

types::global_dof_index DoFAccessor< structdim, dim, spacedim, level_dof_access >::vertex_dof_index ( const unsigned int  vertex,
const unsigned int  i,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const
inherited

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, when hp-capabilities are enabled, 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 cell->active_fe_index(). Alternatively, if fe_index is left to its default value when this function is called on a cell, then this is interpreted as equal to cell->active_fe_index().

◆ mg_vertex_dof_index()

types::global_dof_index DoFAccessor< structdim, dim, spacedim, level_dof_access >::mg_vertex_dof_index ( const int  level,
const unsigned int  vertex,
const unsigned int  i,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const
inherited

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()

types::global_dof_index DoFAccessor< structdim, dim, spacedim, level_dof_access >::dof_index ( const unsigned int  i,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const
inherited

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, when hp-capabilities are enabled, 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()

types::global_dof_index DoFAccessor< structdim, dim, spacedim, level_dof_access >::mg_dof_index ( const int  level,
const unsigned int  i 
) const
inherited

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

◆ n_active_fe_indices()

unsigned int DoFAccessor< structdim, dim, spacedim, level_dof_access >::n_active_fe_indices ( ) const
inherited

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

When hp-capabilities are disabled the answer is, of course, always one. However, when hp-capabilities are enabled, 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()

types::fe_index DoFAccessor< structdim, dim, spacedim, level_dof_access >::nth_active_fe_index ( const unsigned int  n) const
inherited

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()

std::set< types::fe_index > DoFAccessor< structdim, dim, spacedim, level_dof_access >::get_active_fe_indices ( ) const
inherited

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()

bool DoFAccessor< structdim, dim, spacedim, level_dof_access >::fe_index_is_active ( const types::fe_index  fe_index) const
inherited

Return true if the finite element with given index is active on the present object. When the current DoFHandler does not have hp- capabilities, 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() [2/2]

const FiniteElement< dim, spacedim > & DoFAccessor< structdim, dim, spacedim, level_dof_access >::get_fe ( const types::fe_index  fe_index) const
inherited

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.

◆ ExcInvalidObject()

static ::ExceptionBase & DoFAccessor< structdim, dim, spacedim, level_dof_access >::ExcInvalidObject ( )
staticinherited

Exceptions for child classes

Note
The message that will be printed by this exception reads:
"This accessor object has not been " "associated with any DoFHandler object."

◆ ExcVectorNotEmpty()

static ::ExceptionBase & DoFAccessor< structdim, dim, spacedim, level_dof_access >::ExcVectorNotEmpty ( )
staticinherited

Exception

◆ ExcVectorDoesNotMatch()

static ::ExceptionBase & DoFAccessor< structdim, dim, spacedim, level_dof_access >::ExcVectorDoesNotMatch ( )
staticinherited

Exception

◆ ExcMatrixDoesNotMatch()

static ::ExceptionBase & DoFAccessor< structdim, dim, spacedim, level_dof_access >::ExcMatrixDoesNotMatch ( )
staticinherited

Exception

◆ ExcNotActive()

static ::ExceptionBase & DoFAccessor< structdim, dim, spacedim, level_dof_access >::ExcNotActive ( )
staticinherited

A function has been called for a cell which should be active, but is refined.

◆ ExcCantCompareIterators()

static ::ExceptionBase & DoFAccessor< structdim, dim, spacedim, level_dof_access >::ExcCantCompareIterators ( )
staticinherited

Exception

◆ operator==()

bool DoFAccessor< structdim, dim, spacedim, level_dof_access >::operator== ( const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &  ) const
inherited

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!=()

bool DoFAccessor< structdim, dim, spacedim, level_dof_access >::operator!= ( const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &  ) const
inherited

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

◆ set_dof_handler()

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::set_dof_handler ( DoFHandler< dim, spacedim > *  dh)
protectedinherited

Reset the DoF handler pointer.

◆ set_dof_index()

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::set_dof_index ( const unsigned int  i,
const types::global_dof_index  index,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const
protectedinherited

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, when the relevant DoFHandler has hp-capabilities, 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_mg_dof_index()

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::set_mg_dof_index ( const int  level,
const unsigned int  i,
const types::global_dof_index  index 
) const
protectedinherited

◆ set_mg_vertex_dof_index()

void DoFAccessor< structdim, dim, spacedim, level_dof_access >::set_mg_vertex_dof_index ( const int  level,
const unsigned int  vertex,
const unsigned int  i,
const types::global_dof_index  index,
const types::fe_index  fe_index = numbers::invalid_fe_index 
) const
protectedinherited

Friends And Related Symbol Documentation

◆ ::internal::DoFCellAccessorImplementation::Implementation

template<int dimension_, int space_dimension_, bool level_dof_access>
friend struct ::internal::DoFCellAccessorImplementation::Implementation
friend

Definition at line 2113 of file dof_accessor.h.

Member Data Documentation

◆ dim

template<int dimension_, int space_dimension_, bool level_dof_access>
const unsigned int DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::dim = dimension_
static

Extract dimension from DoFHandler.

Definition at line 1330 of file dof_accessor.h.

◆ spacedim

template<int dimension_, int space_dimension_, bool level_dof_access>
const unsigned int DoFCellAccessor< dimension_, space_dimension_, level_dof_access >::spacedim = space_dimension_
static

Extract space dimension from DoFHandler.

Definition at line 1335 of file dof_accessor.h.

◆ dimension

const unsigned int DoFAccessor< structdim, dim, spacedim, level_dof_access >::dimension
staticconstexprinherited

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

Definition at line 217 of file dof_accessor.h.

◆ space_dimension

const unsigned int DoFAccessor< structdim, dim, spacedim, level_dof_access >::space_dimension
staticconstexprinherited

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

Definition at line 223 of file dof_accessor.h.

◆ dof_handler

DoFHandler<dim, spacedim>* DoFAccessor< structdim, dim, spacedim, level_dof_access >::dof_handler
protectedinherited

Store the address of the DoFHandler object to be accessed.

Definition at line 658 of file dof_accessor.h.


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