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 Member Functions | Public Attributes | Private Attributes | Friends | List of all members
FEInterfaceValues< dim, spacedim > Class Template Reference

#include <deal.II/fe/fe_interface_values.h>

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

Public Member Functions

 FEInterfaceValues (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
 
 FEInterfaceValues (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim - 1 > &quadrature, const UpdateFlags update_flags)
 
 FEInterfaceValues (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &quadrature, const UpdateFlags update_flags)
 
 FEInterfaceValues (const hp::MappingCollection< dim, spacedim > &mapping_collection, const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &quadrature_collection, const UpdateFlags update_flags)
 
 FEInterfaceValues (const hp::FECollection< dim, spacedim > &fe_collection, const hp::QCollection< dim - 1 > &quadrature_collection, const UpdateFlags update_flags)
 
template<class CellIteratorType , class CellNeighborIteratorType >
void reinit (const CellIteratorType &cell, const unsigned int face_no, const unsigned int sub_face_no, const CellNeighborIteratorType &cell_neighbor, const unsigned int face_no_neighbor, const unsigned int sub_face_no_neighbor, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
 
template<class CellIteratorType >
void reinit (const CellIteratorType &cell, const unsigned int face_no, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int, const unsigned int fe_index=numbers::invalid_unsigned_int)
 
const FEFaceValuesBase< dim, spacedim > & get_fe_face_values (const unsigned int cell_index) const
 
const Mapping< dim, spacedim > & get_mapping () const
 
const FiniteElement< dim, spacedim > & get_fe () const
 
const Quadrature< dim - 1 > & get_quadrature () const
 
const hp::MappingCollection< dim, spacedim > & get_mapping_collection () const
 
const hp::FECollection< dim, spacedim > & get_fe_collection () const
 
const hp::QCollection< dim - 1 > & get_quadrature_collection () const
 
bool has_hp_capabilities () const
 
UpdateFlags get_update_flags () const
 
const Triangulation< dim, spacedim >::cell_iterator get_cell (const unsigned int cell_index) const
 
unsigned int get_face_number (const unsigned int cell_index) const
 
Functions to query information on a given interface
bool at_boundary () const
 
double JxW (const unsigned int quadrature_point) const
 
const std::vector< double > & get_JxW_values () const
 
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intquadrature_point_indices () const
 
const std::vector< Point< spacedim > > & get_quadrature_points () const
 
unsigned n_current_interface_dofs () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intdof_indices () const
 
std::vector< types::global_dof_indexget_interface_dof_indices () const
 
std::array< unsigned int, 2 > interface_dof_to_dof_indices (const unsigned int interface_dof_index) const
 
Tensor< 1, spacedim > normal (const unsigned int q_point_index) const
 
Access to shape functions
double shape_value (const bool here_or_there, const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Access to jumps in shape function values and their derivatives
double jump_in_shape_values (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
double jump (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 1, spacedim > jump_in_shape_gradients (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 1, spacedim > jump_gradient (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 2, spacedim > jump_in_shape_hessians (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 2, spacedim > jump_hessian (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 3, spacedim > jump_in_shape_3rd_derivatives (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 3, spacedim > jump_3rd_derivative (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Access to the average of shape function values and their derivatives
double average_of_shape_values (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
double average (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 1, spacedim > average_of_shape_gradients (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 1, spacedim > average_gradient (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 2, spacedim > average_of_shape_hessians (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Tensor< 2, spacedim > average_hessian (const unsigned int interface_dof_index, const unsigned int q_point, const unsigned int component=0) const
 
Access to jumps in the function values and derivatives
template<class InputVector >
void get_jump_in_function_values (const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
 
template<class InputVector >
void get_jump_in_function_gradients (const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
 
template<class InputVector >
void get_jump_in_function_hessians (const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
 
template<class InputVector >
void get_jump_in_function_third_derivatives (const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type > > &third_derivatives) const
 
Access to the average of the function values and derivatives
template<class InputVector >
void get_average_of_function_values (const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
 
template<class InputVector >
void get_average_of_function_gradients (const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &gradients) const
 
template<class InputVector >
void get_average_of_function_hessians (const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &hessians) const
 
Extractors Methods to extract individual components
const FEInterfaceViews::Scalar< dim, spacedim > operator[] (const FEValuesExtractors::Scalar &scalar) const
 
const FEInterfaceViews::Vector< dim, spacedim > operator[] (const FEValuesExtractors::Vector &vector) const
 

Public Attributes

const unsigned int n_quadrature_points
 

Private Attributes

std::vector< types::global_dof_indexinterface_dof_indices
 
std::vector< std::array< unsigned int, 2 > > dofmap
 
FEFaceValuesBase< dim, spacedim > * fe_face_values
 
FEFaceValuesBase< dim, spacedim > * fe_face_values_neighbor
 
Data that supports the standard FE implementation
std::unique_ptr< FEFaceValues< dim, spacedim > > internal_fe_face_values
 
std::unique_ptr< FESubfaceValues< dim, spacedim > > internal_fe_subface_values
 
std::unique_ptr< FEFaceValues< dim, spacedim > > internal_fe_face_values_neighbor
 
std::unique_ptr< FESubfaceValues< dim, spacedim > > internal_fe_subface_values_neighbor
 

Friends

template<int , int >
class FEInterfaceViews::Scalar
 
template<int , int >
class FEInterfaceViews::Vector
 

Data that supports the hp-FE implementation

std::unique_ptr< hp::FEFaceValues< dim, spacedim > > internal_hp_fe_face_values
 
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > internal_hp_fe_subface_values
 
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > internal_hp_fe_face_values_neighbor
 
std::unique_ptr< hp::FESubfaceValues< dim, spacedim > > internal_hp_fe_subface_values_neighbor
 
static ::ExceptionBaseExcOnlyAvailableWithoutHP ()
 
static ::ExceptionBaseExcOnlyAvailableWithHP ()
 

Detailed Description

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

FEInterfaceValues is a data structure to access and assemble finite element data on interfaces between two cells of a mesh.

It provides a way to access averages, jump terms, and similar operations used in Discontinuous Galerkin methods on a face between two neighboring cells. This allows the computation of typical mesh-dependent linear or bilinear forms in a similar way as FEValues does for cells and FEFaceValues does for faces. In the literature, the faces between neighboring cells are called "inner interfaces" or "facets".

Internally, this class provides an abstraction for two FEFaceValues objects (or FESubfaceValues when using adaptive refinement). The class introduces a new "interface dof index" that walks over the union of the dof indices of the two FEFaceValues objects. Helper functions allow translating between the new "interface dof index" and the corresponding "cell index" (0 for the first cell, 1 for the second cell) and "dof index" within that cell.

The class is made to be used inside MeshWorker::mesh_loop(). It is intended to be a low level replacement for MeshWorker and LocalIntegrators and a higher level abstraction compared to assembling face terms manually.

Definition at line 1397 of file fe_interface_values.h.

Constructor & Destructor Documentation

◆ FEInterfaceValues() [1/5]

template<int dim, int spacedim = dim>
FEInterfaceValues< dim, spacedim >::FEInterfaceValues ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim - 1 > &  quadrature,
const UpdateFlags  update_flags 
)

Construct the FEInterfaceValues with a single FiniteElement (same on both sides of the facet). The FEFaceValues objects will be initialized with the given mapping, quadrature, and update_flags.

◆ FEInterfaceValues() [2/5]

template<int dim, int spacedim = dim>
FEInterfaceValues< dim, spacedim >::FEInterfaceValues ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const hp::QCollection< dim - 1 > &  quadrature,
const UpdateFlags  update_flags 
)

The same as above but taking a collection of quadrature rules so that different quadrature rules can be assigned to different faces.

◆ FEInterfaceValues() [3/5]

template<int dim, int spacedim = dim>
FEInterfaceValues< dim, spacedim >::FEInterfaceValues ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim - 1 > &  quadrature,
const UpdateFlags  update_flags 
)

Construct the FEInterfaceValues with a single FiniteElement and a Q1 Mapping.

See the constructor above.

◆ FEInterfaceValues() [4/5]

template<int dim, int spacedim = dim>
FEInterfaceValues< dim, spacedim >::FEInterfaceValues ( const hp::MappingCollection< dim, spacedim > &  mapping_collection,
const hp::FECollection< dim, spacedim > &  fe_collection,
const hp::QCollection< dim - 1 > &  quadrature_collection,
const UpdateFlags  update_flags 
)

Construct the FEInterfaceValues object with different FiniteElements assigned to different faces.

◆ FEInterfaceValues() [5/5]

template<int dim, int spacedim = dim>
FEInterfaceValues< dim, spacedim >::FEInterfaceValues ( const hp::FECollection< dim, spacedim > &  fe_collection,
const hp::QCollection< dim - 1 > &  quadrature_collection,
const UpdateFlags  update_flags 
)

Same as above, but using the default linear mapping.

Member Function Documentation

◆ reinit() [1/2]

template<int dim, int spacedim = dim>
template<class CellIteratorType , class CellNeighborIteratorType >
void FEInterfaceValues< dim, spacedim >::reinit ( const CellIteratorType &  cell,
const unsigned int  face_no,
const unsigned int  sub_face_no,
const CellNeighborIteratorType &  cell_neighbor,
const unsigned int  face_no_neighbor,
const unsigned int  sub_face_no_neighbor,
const unsigned int  q_index = numbers::invalid_unsigned_int,
const unsigned int  mapping_index = numbers::invalid_unsigned_int,
const unsigned int  fe_index = numbers::invalid_unsigned_int 
)

Re-initialize this object to be used on a new interface given by two faces of two neighboring cells. The cell and cell_neighbor cells will be referred to through cell_index zero and one after this call in all places where one needs to identify the two cells adjacent to the interface.

Use numbers::invalid_unsigned_int for sub_face_no or sub_face_no_neighbor to indicate that you want to work on the entire face, not a sub-face.

The arguments (including their order) are identical to the face_worker arguments in MeshWorker::mesh_loop().

In order to do integration on a face or sub-face, this object will have to choose what quadrature formula to use. This is simple if you initialized the current FEInterfaceValues object with just a single FiniteElement and a single (face) Quadrature object, because in that case there is only one element and quadrature that will be used for all interfaces. But it is not so simple in the hp case where there may be different elements used on different cells, and different quadrature formulas should be used on different faces; one may also want to use a different mapping for different faces. As a consequence, you would have initialized the current object with a hp::FECollection, hp::QCollection, and possible an hp::MappingCollection object. The question then is: For a given face or subface, which quadrature and mapping should be used? The following decision tree will then be used:

  1. If the q_index and mapping_index arguments to this function are explicitly specified (rather than leaving them at their default values), then these indices will be used to select which element of the hp::QCollection and hp::MappingCollection passed to the constructor should serve as the quadrature and mapping to be used.
  2. If one of these arguments is left at its default value, then the function will need to choose a quadrature and/or mapping that is appropriate for the two finite element spaces used on the two cells adjacent to the current interface. As the first choice, if the quadrature or mapping collection we are considering has only one element, then that is clearly the one that should be used.
  3. If the quadrature or mapping collection have multiple elements, then we need to dig further. For quadrature objects, we can compare whether the two quadrature objects that correspond to the active_fe_index values of the two adjacent cells are identical (i.e., have quadrature points at the same locations, and have the same weights). If this is so, then it does not matter which one of the two we take, and we choose one or the other.
  4. If this has still not helped, we try to find out which of the two finite element spaces on the two adjacent cells is "larger" (say, if you had used \(Q_2\) and \(Q_4\) elements on the two adjacent cells, then the \(Q_4\) element is the larger one); the determination of which space is "larger" is made using the hp::FECollection::find_dominated_fe() function, which is not necessarily intended for this kind of query, but yields a result that serves just fine for our purposes here. We then operate on the assumption that the quadrature object associated with the "larger" of the two spaces is the appropriate one to use for the face that separates these two spaces.
    • If this function returns that one of the two elements in question is dominated by the other, then presumably it is "larger" one and we take the quadrature formula and mapping that corresponds to this "larger" element is. For example, for the \(Q_2\) element mentioned above, one would generally use a QGauss(3) quadrature formula, whereas for the \(Q_4\) element, one would use QGauss(5). To integrate jump and average terms on the interface between cells using these two elements, QGauss(5) is appropriate. Because, typically, people will order elements in the hp::FECollection in the same order as the quadrature and mapping objects in hp::QCollection and hp::MappingCollection, this function will use the index of the "larger" element in the hp::FECollection to also index into the hp::QCollection and hp::MappingCollection to retrieve quadrature and mapping objects appropriate for the current face.
    • There are cases where neither element dominates the other. For example, if one uses \(Q_2\times Q_1\) and \(Q_1\times Q_2\) elements on neighboring cells, neither of the two spaces dominates the other – or, in the context of the current function, neither space is "larger" than the other. In that case, there is no way for the current function to determine quadrature and mapping objects associated with the two elements are the appropriate ones. If that happens, you will get an error – and the only way to avoid the error is to explicitly specify for these interfaces which quadrature and mapping objects you want to use, by providing non-default values for the q_index and mapping_index arguments to this function.
Parameters
[in]cellAn iterator to the first cell adjacent to the interface.
[in]face_noAn integer identifying which face of the first cell the interface is on.
[in]sub_face_noAn integer identifying the subface (child) of the face (identified by the previous two arguments) that the interface corresponds to. If equal to numbers::invalid_unsigned_int, then the interface is considered to be the entire face.
[in]cell_neighborAn iterator to the second cell adjacent to the interface. The type of this iterator does not have to equal that of cell, but must be convertible to it. This allows using an active cell iterator for cell, and cell->neighbor(f) for cell_neighbor, since the return type of cell->neighbor(f) is simply a cell iterator (not necessarily an active cell iterator).
[in]face_no_neighborLike face_no, just for the neighboring cell.
[in]sub_face_no_neighborLike sub_face_no, just for the neighboring cell.
[in]q_indexThe index of the quadrature object within the hp::QCollection passed to the constructor to use on the current interface. See the documentation above what happens if this argument is not explicitly provided but left at its default value.
[in]mapping_indexThe index of the mapping object within the hp::MappingCollection passed to the constructor to use on the current interface. See the documentation above what happens if this argument is not explicitly provided but left at its default value.
[in]fe_indexIf left at its default, use the finite element within the hp::FECollection passed to the constructor as given by the dominating finite element across the interface (only used if the FEInterface object is initialized with an hp::FECollection, an hp::QCollection, and possibly an hp::MappingCollection).

◆ reinit() [2/2]

template<int dim, int spacedim = dim>
template<class CellIteratorType >
void FEInterfaceValues< dim, spacedim >::reinit ( const CellIteratorType &  cell,
const unsigned int  face_no,
const unsigned int  q_index = numbers::invalid_unsigned_int,
const unsigned int  mapping_index = numbers::invalid_unsigned_int,
const unsigned int  fe_index = numbers::invalid_unsigned_int 
)

Re-initialize this object to be used on an interface given by a single face face_no of the cell cell. This is useful to use FEInterfaceValues on boundaries of the domain.

As a consequence, members like jump() will assume a value of zero for the values on the "other" side. Note that no sub_face_number is needed as a boundary face can not neighbor a finer cell.

After calling this function at_boundary() will return true.

Parameters
[in]cellAn iterator to the first cell adjacent to the interface.
[in]face_noAn integer identifying which face of the first cell the interface is on.
[in]q_indexThis argument selects which quadrature formula to use See the discussion in the documentation of the other reinit() function for what happens when this argument is left at its default value.
[in]mapping_indexThis argument selects which mapping to use See the discussion in the documentation of the other reinit() function for what happens when this argument is left at its default value.
[in]fe_indexIf left at its default, use the finite element within the hp::FECollection passed to the constructor as given by the dominating finite element across the interface (only used if the FEInterface object is initialized with an hp::FECollection, an hp::QCollection, and possibly an hp::MappingCollection).

◆ get_fe_face_values()

template<int dim, int spacedim = dim>
const FEFaceValuesBase< dim, spacedim > & FEInterfaceValues< dim, spacedim >::get_fe_face_values ( const unsigned int  cell_index) const

Return a reference to the FEFaceValues or FESubfaceValues object of the specified cell of the interface.

The cell_index is either 0 or 1 and corresponds to the cell index returned by interface_dof_to_cell_and_dof_index().

◆ get_mapping()

template<int dim, int spacedim = dim>
const Mapping< dim, spacedim > & FEInterfaceValues< dim, spacedim >::get_mapping ( ) const

Constant reference to the selected mapping object.

◆ get_fe()

template<int dim, int spacedim = dim>
const FiniteElement< dim, spacedim > & FEInterfaceValues< dim, spacedim >::get_fe ( ) const

Constant reference to the selected finite element object.

◆ get_quadrature()

template<int dim, int spacedim = dim>
const Quadrature< dim - 1 > & FEInterfaceValues< dim, spacedim >::get_quadrature ( ) const

Return a reference to the quadrature object in use.

◆ get_mapping_collection()

template<int dim, int spacedim = dim>
const hp::MappingCollection< dim, spacedim > & FEInterfaceValues< dim, spacedim >::get_mapping_collection ( ) const

Return a reference to the used mapping.

◆ get_fe_collection()

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

Return a reference to the selected finite element object.

◆ get_quadrature_collection()

template<int dim, int spacedim = dim>
const hp::QCollection< dim - 1 > & FEInterfaceValues< dim, spacedim >::get_quadrature_collection ( ) const

Return a reference to the face quadrature object in use.

◆ has_hp_capabilities()

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

Returns a boolean indicating whether or not this FEInterfaceValues object has hp-capabilities enabled.

◆ get_update_flags()

template<int dim, int spacedim = dim>
UpdateFlags FEInterfaceValues< dim, spacedim >::get_update_flags ( ) const

Return the update flags set.

◆ get_cell()

template<int dim, int spacedim = dim>
const Triangulation< dim, spacedim >::cell_iterator FEInterfaceValues< dim, spacedim >::get_cell ( const unsigned int  cell_index) const

Return a triangulation iterator to the current cell of the interface.

The cell_index is either 0 or 1 and corresponds to the cell index returned by interface_dof_to_cell_and_dof_index().

◆ get_face_number()

template<int dim, int spacedim = dim>
unsigned int FEInterfaceValues< dim, spacedim >::get_face_number ( const unsigned int  cell_index) const

Return the number of the face on the interface selected the last time the reinit() function was called.

The cell_index is either 0 or 1 and corresponds to the cell index returned by interface_dof_to_cell_and_dof_index().

◆ at_boundary()

template<int dim, int spacedim = dim>
bool FEInterfaceValues< dim, spacedim >::at_boundary ( ) const

Return if the current interface is a boundary face or an internal face with two adjacent cells.

See the corresponding reinit() functions for details.

◆ JxW()

template<int dim, int spacedim = dim>
double FEInterfaceValues< dim, spacedim >::JxW ( const unsigned int  quadrature_point) const

Mapped quadrature weight. This value equals the mapped surface element times the weight of the quadrature point.

You can think of the quantity returned by this function as the surface element \(ds\) in the integral that we implement here by quadrature.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_JxW_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_JxW_values()

template<int dim, int spacedim = dim>
const std::vector< double > & FEInterfaceValues< dim, spacedim >::get_JxW_values ( ) const

Return the vector of JxW values for each quadrature point.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_JxW_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_normal_vectors()

template<int dim, int spacedim = dim>
const std::vector< Tensor< 1, spacedim > > & FEInterfaceValues< dim, spacedim >::get_normal_vectors ( ) const

Return the normal vector of the interface in each quadrature point.

The return value is identical to get_fe_face_values(0).get_normal_vectors() and therefore, are outside normal vectors from the perspective of the first cell of this interface.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_normal_vectors flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ quadrature_point_indices()

template<int dim, int spacedim = dim>
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FEInterfaceValues< dim, spacedim >::quadrature_point_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero to n_quadrature_points. This allows to write code using range-based for loops.

See also
deal.II and Modern C++ standards

◆ get_quadrature_points()

template<int dim, int spacedim = dim>
const std::vector< Point< spacedim > > & FEInterfaceValues< dim, spacedim >::get_quadrature_points ( ) const

Return a reference to the quadrature points in real space.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_quadrature_points flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ n_current_interface_dofs()

template<int dim, int spacedim = dim>
unsigned FEInterfaceValues< dim, spacedim >::n_current_interface_dofs ( ) const

Return the number of DoFs (or shape functions) on the current interface.

Note
This number is only available after a call to reinit() and can change from one call to reinit() to the next. For example, on a boundary interface it is equal to the number of dofs of the single FEFaceValues object, while it is twice that for an interior interface for a DG element. For a continuous element, it is slightly smaller because the two cells on the interface share some of the dofs.

◆ dof_indices()

template<int dim, int spacedim = dim>
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FEInterfaceValues< dim, spacedim >::dof_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero (inclusive) to n_current_interface_dofs() (exclusive). This allows one to write code using range-based for loops of the following kind:

FullMatrix<double> cell_matrix (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
cell_matrix = 0;
for (const auto &face : cell->face_iterators())
{
fe_iv.values.reinit(cell, face, ...);
for (const auto q : fe_iv.quadrature_point_indices())
for (const auto i : fe_iv.dof_indices())
for (const auto j : fe_iv.dof_indices())
cell_matrix(i,j) += ...; // Do something for DoF indices
// (i,j) at quadrature point q
}
}
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const

Here, we are looping over all degrees of freedom on all cell interfaces, with i and j taking on all valid indices for interface degrees of freedom, as defined by the finite element passed to fe_iv.

◆ get_interface_dof_indices()

template<int dim, int spacedim = dim>
std::vector< types::global_dof_index > FEInterfaceValues< dim, spacedim >::get_interface_dof_indices ( ) const

Return the set of joint DoF indices. This includes indices from both cells. If reinit was called with an active cell iterator, the indices are based on the active indices (returned by DoFCellAccessor::get_dof_indices() ), in case of level cell (that is, if is_level_cell() return true ) the mg dof indices are returned.

Note
This function is only available after a call to reinit() and can change from one call to reinit() to the next.

◆ interface_dof_to_dof_indices()

template<int dim, int spacedim = dim>
std::array< unsigned int, 2 > FEInterfaceValues< dim, spacedim >::interface_dof_to_dof_indices ( const unsigned int  interface_dof_index) const

Convert an interface dof index into the corresponding local DoF indices of the two cells. If an interface DoF is only active on one of the cells, the other index will be numbers::invalid_unsigned_int.

For discontinuous finite elements, each interface dof is located on exactly one side of the interface and, consequently, only one of the two values returned is valid (i.e., different from numbers::invalid_unsigned_int).

Note
This function is only available after a call to reinit() and the returned values may change from one call to reinit() to the next.

◆ normal()

template<int dim, int spacedim = dim>
Tensor< 1, spacedim > FEInterfaceValues< dim, spacedim >::normal ( const unsigned int  q_point_index) const

Return the normal in a given quadrature point.

The normal points in outwards direction as seen from the first cell of this interface.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_normal_vectors flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_value()

template<int dim, int spacedim = dim>
double FEInterfaceValues< dim, spacedim >::shape_value ( const bool  here_or_there,
const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return component component of the value of the shape function with interface dof index interface_dof_index in quadrature point q_point.

The argument here_or_there selects between the value on cell 0 (here, true) and cell 1 (there, false). You can also interpret it as "upstream" (true) and "downstream" (false) as defined by the direction of the normal vector in this quadrature point. If here_or_there is true, the shape functions from the first cell of the interface is used.

In other words, this function returns the limit of the value of the shape function in the given quadrature point when approaching it from one of the two cells of the interface.

Note
This function is typically used to pick the upstream or downstream value based on a direction. This can be achieved by using (direction * normal)>0 as the first argument of this function.

◆ jump_in_shape_values()

template<int dim, int spacedim = dim>
double FEInterfaceValues< dim, spacedim >::jump_in_shape_values ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return the jump \(\jump{u}=u_{\text{cell0}} - u_{\text{cell1}}\) on the interface for the shape function interface_dof_index at the quadrature point q_point of component component.

Note that one can define the jump in different ways (the value "there" minus the value "here", or the other way around; both are used in the finite element literature). The definition here uses "value here minus value there", as seen from the first cell.

If this is a boundary face (at_boundary() returns true), then \(\jump{u}=u_{\text{cell0}}\), that is "the value here (minus zero)".

Note
The name of the function is supposed to be read as "the jump (singular) in the values (plural: one or two possible values) of the shape function (singular)".

◆ jump()

template<int dim, int spacedim = dim>
double FEInterfaceValues< dim, spacedim >::jump ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

The same as above.

Deprecated:
Use the jump_in_shape_values() function instead.

◆ jump_in_shape_gradients()

template<int dim, int spacedim = dim>
Tensor< 1, spacedim > FEInterfaceValues< dim, spacedim >::jump_in_shape_gradients ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return the jump in the gradient \(\jump{\nabla u}=\nabla u_{\text{cell0}} - \nabla u_{\text{cell1}}\) on the interface for the shape function interface_dof_index at the quadrature point q_point of component component.

If this is a boundary face (at_boundary() returns true), then \(\jump{\nabla u}=\nabla u_{\text{cell0}}\).

Note
The name of the function is supposed to be read as "the jump (singular) in the gradients (plural: one or two possible gradients) of the shape function (singular)".

◆ jump_gradient()

template<int dim, int spacedim = dim>
Tensor< 1, spacedim > FEInterfaceValues< dim, spacedim >::jump_gradient ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

The same as above.

Deprecated:
Use the jump_in_shape_gradients() function instead.

◆ jump_in_shape_hessians()

template<int dim, int spacedim = dim>
Tensor< 2, spacedim > FEInterfaceValues< dim, spacedim >::jump_in_shape_hessians ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return the jump in the Hessian \(\jump{\nabla^2 u} = \nabla^2 u_{\text{cell0}} - \nabla^2 u_{\text{cell1}}\) on the interface for the shape function interface_dof_index at the quadrature point q_point of component component.

If this is a boundary face (at_boundary() returns true), then \(\jump{\nabla^2 u} = \nabla^2 u_{\text{cell0}}\).

Note
The name of the function is supposed to be read as "the jump (singular) in the Hessians (plural: one or two possible values for the derivative) of the shape function (singular)".

◆ jump_hessian()

template<int dim, int spacedim = dim>
Tensor< 2, spacedim > FEInterfaceValues< dim, spacedim >::jump_hessian ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

The same as above.

Deprecated:
Use the jump_in_shape_hessians() function instead.

◆ jump_in_shape_3rd_derivatives()

template<int dim, int spacedim = dim>
Tensor< 3, spacedim > FEInterfaceValues< dim, spacedim >::jump_in_shape_3rd_derivatives ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return the jump in the third derivative \(\jump{\nabla^3 u} = \nabla^3 u_{\text{cell0}} - \nabla^3 u_{\text{cell1}}\) on the interface for the shape function interface_dof_index at the quadrature point q_point of component component.

If this is a boundary face (at_boundary() returns true), then \(\jump{\nabla^3 u} = \nabla^3 u_{\text{cell0}}\).

Note
The name of the function is supposed to be read as "the jump (singular) in the third derivatives (plural: one or two possible values for the derivative) of the shape function (singular)".

◆ jump_3rd_derivative()

template<int dim, int spacedim = dim>
Tensor< 3, spacedim > FEInterfaceValues< dim, spacedim >::jump_3rd_derivative ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

The same as above.

Deprecated:
Use the jump_in_shape_3rd_derivatives() function instead.

◆ average_of_shape_values()

template<int dim, int spacedim = dim>
double FEInterfaceValues< dim, spacedim >::average_of_shape_values ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return the average \(\average{u}=\frac{1}{2}u_{\text{cell0}} + \frac{1}{2}u_{\text{cell1}}\) on the interface for the shape function interface_dof_index at the quadrature point q_point of component component.

If this is a boundary face (at_boundary() returns true), then \(\average{u}=u_{\text{cell0}}\).

Note
The name of the function is supposed to be read as "the average (singular) of the values (plural: one or two possible values) of the shape function (singular)".

◆ average()

template<int dim, int spacedim = dim>
double FEInterfaceValues< dim, spacedim >::average ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

The same as above.

Deprecated:
Use the average_of_shape_values() function instead.

◆ average_of_shape_gradients()

template<int dim, int spacedim = dim>
Tensor< 1, spacedim > FEInterfaceValues< dim, spacedim >::average_of_shape_gradients ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return the average of the gradient \(\average{\nabla u} = \frac{1}{2}\nabla u_{\text{cell0}} + \frac{1}{2} \nabla u_{\text{cell1}}\) on the interface for the shape function interface_dof_index at the quadrature point q_point of component component.

If this is a boundary face (at_boundary() returns true), then \(\average{\nabla u}=\nabla u_{\text{cell0}}\).

Note
The name of the function is supposed to be read as "the average (singular) of the gradients (plural: one or two possible values for the gradient) of the shape function (singular)".

◆ average_gradient()

template<int dim, int spacedim = dim>
Tensor< 1, spacedim > FEInterfaceValues< dim, spacedim >::average_gradient ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

The same as above.

Deprecated:
Use the average_of_shape_gradients() function instead.

◆ average_of_shape_hessians()

template<int dim, int spacedim = dim>
Tensor< 2, spacedim > FEInterfaceValues< dim, spacedim >::average_of_shape_hessians ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

Return the average of the Hessian \(\average{\nabla^2 u} = \frac{1}{2}\nabla^2 u_{\text{cell0}} + \frac{1}{2} \nabla^2 u_{\text{cell1}}\) on the interface for the shape function interface_dof_index at the quadrature point q_point of component component.

If this is a boundary face (at_boundary() returns true), then \(\average{\nabla^2 u}=\nabla^2 u_{\text{cell0}}\).

Note
The name of the function is supposed to be read as "the average (singular) of the Hessians (plural: one or two possible values for the second derivatives) of the shape function (singular)".

◆ average_hessian()

template<int dim, int spacedim = dim>
Tensor< 2, spacedim > FEInterfaceValues< dim, spacedim >::average_hessian ( const unsigned int  interface_dof_index,
const unsigned int  q_point,
const unsigned int  component = 0 
) const

The same as above.

Deprecated:
Use the average_of_shape_hessians() function instead.

◆ get_jump_in_function_values()

template<int dim, int spacedim = dim>
template<class InputVector >
void FEInterfaceValues< dim, spacedim >::get_jump_in_function_values ( const InputVector &  fe_function,
std::vector< typename InputVector::value_type > &  values 
) const

Return the jump in the values of the finite element function characterized by fe_function at the quadrature points of the cell interface selected the last time the reinit function of the FEInterfaceValues object was called.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jump_in_function_gradients()

template<int dim, int spacedim = dim>
template<class InputVector >
void FEInterfaceValues< dim, spacedim >::get_jump_in_function_gradients ( const InputVector &  fe_function,
std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &  gradients 
) const

Return the jump in the gradients of the finite element function characterized by fe_function at the quadrature points of the cell interface selected the last time the reinit function of the FEInterfaceValues object was called.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jump_in_function_hessians()

template<int dim, int spacedim = dim>
template<class InputVector >
void FEInterfaceValues< dim, spacedim >::get_jump_in_function_hessians ( const InputVector &  fe_function,
std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &  hessians 
) const

Return the jump in the Hessians of the finite element function characterized by fe_function at the quadrature points of the cell interface selected the last time the reinit function of the FEInterfaceValues object was called.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jump_in_function_third_derivatives()

template<int dim, int spacedim = dim>
template<class InputVector >
void FEInterfaceValues< dim, spacedim >::get_jump_in_function_third_derivatives ( const InputVector &  fe_function,
std::vector< Tensor< 3, spacedim, typename InputVector::value_type > > &  third_derivatives 
) const

Return the jump in the third derivatives of the finite element function characterized by fe_function at the quadrature points of the cell interface selected the last time the reinit function of the FEInterfaceValues object was called.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_third_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_average_of_function_values()

template<int dim, int spacedim = dim>
template<class InputVector >
void FEInterfaceValues< dim, spacedim >::get_average_of_function_values ( const InputVector &  fe_function,
std::vector< typename InputVector::value_type > &  values 
) const

Return the average of the values of the finite element function characterized by fe_function at the quadrature points of the cell interface selected the last time the reinit function of the FEInterfaceValues object was called.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_average_of_function_gradients()

template<int dim, int spacedim = dim>
template<class InputVector >
void FEInterfaceValues< dim, spacedim >::get_average_of_function_gradients ( const InputVector &  fe_function,
std::vector< Tensor< 1, spacedim, typename InputVector::value_type > > &  gradients 
) const

Return the average of the gradients of the finite element function characterized by fe_function at the quadrature points of the cell interface selected the last time the reinit function of the FEInterfaceValues object was called.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_average_of_function_hessians()

template<int dim, int spacedim = dim>
template<class InputVector >
void FEInterfaceValues< dim, spacedim >::get_average_of_function_hessians ( const InputVector &  fe_function,
std::vector< Tensor< 2, spacedim, typename InputVector::value_type > > &  hessians 
) const

Return the average of the Hessians of the finite element function characterized by fe_function at the quadrature points of the cell interface selected the last time the reinit function of the FEInterfaceValues object was called.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ operator[]() [1/2]

template<int dim, int spacedim = dim>
const FEInterfaceViews::Scalar< dim, spacedim > FEInterfaceValues< dim, spacedim >::operator[] ( const FEValuesExtractors::Scalar scalar) const

Create a view of the current FEInterfaceValues object that represents a particular scalar component of the possibly vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews.

◆ operator[]() [2/2]

template<int dim, int spacedim = dim>
const FEInterfaceViews::Vector< dim, spacedim > FEInterfaceValues< dim, spacedim >::operator[] ( const FEValuesExtractors::Vector vector) const

Create a view of the current FEInterfaceValues object that represents a set of dim scalar components (i.e. a vector) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews.

Friends And Related Symbol Documentation

◆ FEInterfaceViews::Scalar

template<int dim, int spacedim = dim>
template<int , int >
friend class FEInterfaceViews::Scalar
friend

Definition at line 2366 of file fe_interface_values.h.

◆ FEInterfaceViews::Vector

template<int dim, int spacedim = dim>
template<int , int >
friend class FEInterfaceViews::Vector
friend

Definition at line 2368 of file fe_interface_values.h.

Member Data Documentation

◆ n_quadrature_points

template<int dim, int spacedim = dim>
const unsigned int FEInterfaceValues< dim, spacedim >::n_quadrature_points

Number of quadrature points.

Definition at line 1403 of file fe_interface_values.h.

◆ interface_dof_indices

template<int dim, int spacedim = dim>
std::vector<types::global_dof_index> FEInterfaceValues< dim, spacedim >::interface_dof_indices
private

The list of DoF indices for the current interface, filled in reinit().

Definition at line 2257 of file fe_interface_values.h.

◆ dofmap

template<int dim, int spacedim = dim>
std::vector<std::array<unsigned int, 2> > FEInterfaceValues< dim, spacedim >::dofmap
private

The mapping from interface dof to the two local dof indices of the FEFaceValues objects. If an interface DoF is only active on one of the cells, the other one will have numbers::invalid_unsigned_int.

Definition at line 2264 of file fe_interface_values.h.

◆ fe_face_values

template<int dim, int spacedim = dim>
FEFaceValuesBase<dim, spacedim>* FEInterfaceValues< dim, spacedim >::fe_face_values
private

Pointer to internal_fe_face_values or internal_fe_subface_values, respectively as determined in reinit().

Definition at line 2270 of file fe_interface_values.h.

◆ fe_face_values_neighbor

template<int dim, int spacedim = dim>
FEFaceValuesBase<dim, spacedim>* FEInterfaceValues< dim, spacedim >::fe_face_values_neighbor
private

Pointer to internal_fe_face_values_neighbor, internal_fe_subface_values_neighbor, or nullptr, respectively as determined in reinit().

Definition at line 2277 of file fe_interface_values.h.

◆ internal_fe_face_values

template<int dim, int spacedim = dim>
std::unique_ptr<FEFaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_fe_face_values
private

The FEFaceValues object for the current cell.

Definition at line 2287 of file fe_interface_values.h.

◆ internal_fe_subface_values

template<int dim, int spacedim = dim>
std::unique_ptr<FESubfaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_fe_subface_values
private

The FEFaceValues object for the current cell if the cell is refined.

Definition at line 2292 of file fe_interface_values.h.

◆ internal_fe_face_values_neighbor

template<int dim, int spacedim = dim>
std::unique_ptr<FEFaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_fe_face_values_neighbor
private

The FEFaceValues object for the neighboring cell.

Definition at line 2297 of file fe_interface_values.h.

◆ internal_fe_subface_values_neighbor

template<int dim, int spacedim = dim>
std::unique_ptr<FESubfaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_fe_subface_values_neighbor
private

The FEFaceValues object for the neighboring cell if the cell is refined.

Definition at line 2303 of file fe_interface_values.h.

◆ internal_hp_fe_face_values

template<int dim, int spacedim = dim>
std::unique_ptr<hp::FEFaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_hp_fe_face_values
private

An hp::FEValues object for the FEFaceValues on the present cell.

Definition at line 2316 of file fe_interface_values.h.

◆ internal_hp_fe_subface_values

template<int dim, int spacedim = dim>
std::unique_ptr<hp::FESubfaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_hp_fe_subface_values
private

An hp::FEValues object for the FESubfaceValues on the present cell.

Definition at line 2323 of file fe_interface_values.h.

◆ internal_hp_fe_face_values_neighbor

template<int dim, int spacedim = dim>
std::unique_ptr<hp::FEFaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_hp_fe_face_values_neighbor
private

An hp::FEValues object for the FEFaceValues on the neighbor of the present cell.

Definition at line 2330 of file fe_interface_values.h.

◆ internal_hp_fe_subface_values_neighbor

template<int dim, int spacedim = dim>
std::unique_ptr<hp::FESubfaceValues<dim, spacedim> > FEInterfaceValues< dim, spacedim >::internal_hp_fe_subface_values_neighbor
private

An hp::FEValues object for the FESubfaceValues on the neighboring cell.

Definition at line 2337 of file fe_interface_values.h.


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