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 | Private Member Functions | Private Attributes | List of all members
NonMatching::FEInterfaceValues< dim > Class Template Reference

#include <deal.II/non_matching/fe_values.h>

Public Types

using AdditionalData = typename FaceQuadratureGenerator< dim >::AdditionalData
 

Public Member Functions

template<class VectorType >
 FEInterfaceValues (const hp::FECollection< dim > &fe_collection, const Quadrature< 1 > &quadrature, const RegionUpdateFlags region_update_flags, const MeshClassifier< dim > &mesh_classifier, const DoFHandler< dim > &dof_handler, const VectorType &level_set, const AdditionalData &additional_data=AdditionalData())
 
template<class VectorType >
 FEInterfaceValues (const hp::MappingCollection< dim > &mapping_collection, const hp::FECollection< dim > &fe_collection, const hp::QCollection< dim - 1 > &q_collection, const hp::QCollection< 1 > &q_collection_1d, const RegionUpdateFlags region_update_flags, const MeshClassifier< dim > &mesh_classifier, const DoFHandler< dim > &dof_handler, const VectorType &level_set, const AdditionalData &additional_data=AdditionalData())
 
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)
 
template<class CellIteratorType >
void reinit (const CellIteratorType &cell, const unsigned int face_no)
 
const std_cxx17::optional<::FEInterfaceValues< dim > > & get_inside_fe_values () const
 
const std_cxx17::optional<::FEInterfaceValues< dim > > & get_outside_fe_values () const
 

Private Member Functions

void initialize (const hp::QCollection< dim - 1 > &q_collection)
 
template<bool level_dof_access>
void do_reinit (const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > &cell, const unsigned int face_no, const std::function< void(::FEInterfaceValues< dim > &)> &call_reinit)
 

Private Attributes

const SmartPointer< const hp::MappingCollection< dim > > mapping_collection
 
const SmartPointer< const hp::FECollection< dim > > fe_collection
 
const hp::QCollection< 1 > q_collection_1D
 
LocationToLevelSet current_face_location
 
unsigned int active_fe_index
 
const RegionUpdateFlags region_update_flags
 
const SmartPointer< const MeshClassifier< dim > > mesh_classifier
 
std::deque< std_cxx17::optional<::FEInterfaceValues< dim > > > fe_values_inside_full_quadrature
 
std::deque< std_cxx17::optional<::FEInterfaceValues< dim > > > fe_values_outside_full_quadrature
 
std_cxx17::optional<::FEInterfaceValues< dim > > fe_values_inside
 
std_cxx17::optional<::FEInterfaceValues< dim > > fe_values_outside
 
DiscreteFaceQuadratureGenerator< dim > face_quadrature_generator
 

Detailed Description

template<int dim>
class NonMatching::FEInterfaceValues< dim >

This class is intended to facilitate assembling interface terms on faces in immersed (in the sense of cut) finite element methods. These types of terms occur mainly in cut discontinuous Galerkin methods. This class works analogously to NonMatching::FEValues. The domain is assumed to be described by a level set function, \(\psi : \mathbb{R}^{dim} \to \mathbb{R}\), and this class assumes that we want to integrate over two different regions of each face, \(F\):

\[ N = \{x \in F : \psi(x) < 0 \}, \\ P = \{x \in F : \psi(x) > 0 \}, \]

which we as before refer to as the "inside" and "outside" regions of the face.

When calling one of the reinit functions of this class, this class will create quadrature rules over these regions in the background and set up FEInterfaceValues objects with these quadratures. These objects can then be accessed through one of the functions: get_inside_fe_values() or get_outside_fe_values(). For the same reason as described in NonMatching::FEValues, the returned the FEInterfaceValues objects that get_inside/outside_fe_values() returns is wrapped in a std_cxx17::optional. Assembling using this class would then typically be done in the following way:

NonMatching::FEInterfaceValues<dim> fe_interface_values(...)
for (const auto &cell : dof_handler.active_cell_iterators())
{
for (const unsigned int face_index : cell->face_indices())
{
if (cell->at_boundary(face_index))
fe_interface_values.reinit(cell, face_index);
else
fe_interface_values.reinit(cell,
face_index,
subface_index,
cell->neighbor(face_index),
cell->neighbor_of_neighbor(face_index),
neighbor_subface_index);
const std_cxx17::optional<::FEInterfaceValues<dim>>
&inside_fe_values = fe_interface_values.get_inside_fe_values();
// Check if the returned optional contains a value
if (inside_fe_values)
{
// Assemble locally
}
}
}
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)

To reduce the amount of work, the reinit() function of this class uses the MeshClassifier passed to the constructor to check how the incoming cell relates to the level set function. The immersed quadrature rules are only generated if the cell is intersected. If the cell is completely inside or outside, it returns a cached FEInterfaceValues object created with a quadrature over the reference cell: \([0, 1]^{dim-1}\).

Definition at line 437 of file fe_values.h.

Member Typedef Documentation

◆ AdditionalData

template<int dim>
using NonMatching::FEInterfaceValues< dim >::AdditionalData = typename FaceQuadratureGenerator<dim>::AdditionalData

Definition at line 440 of file fe_values.h.

Constructor & Destructor Documentation

◆ FEInterfaceValues() [1/2]

template<int dim>
template<class VectorType >
FEInterfaceValues< dim >::FEInterfaceValues ( const hp::FECollection< dim > &  fe_collection,
const Quadrature< 1 > &  quadrature,
const RegionUpdateFlags  region_update_flags,
const MeshClassifier< dim > &  mesh_classifier,
const DoFHandler< dim > &  dof_handler,
const VectorType &  level_set,
const AdditionalData additional_data = AdditionalData() 
)

Constructor.

Parameters
fe_collectionCollection of FiniteElements to be used.
quadrature1-dimensional quadrature rule used to generate the immersed quadrature rules. See the QuadratureGenerator class. On the non intersected faces a tensor product of this quadrature will be used.
mesh_classifierObject used to determine when the immersed quadrature rules need to be generated.
region_update_flagsStruct storing UpdateFlags for the inside/outside region of the cell.
dof_handlerThe DoFHandler associated with the discrete level set function.
level_setThe degrees of freedom of the discrete level set function.
additional_dataAdditional data passed to the QuadratureGenerator class.
Note
Pointers to fe_collection, mesh_classifier, dof_handler, and level_set are stored internally, so these need to have a longer life span than the instance of this class.

Definition at line 274 of file fe_values.cc.

◆ FEInterfaceValues() [2/2]

template<int dim>
template<class VectorType >
FEInterfaceValues< dim >::FEInterfaceValues ( const hp::MappingCollection< dim > &  mapping_collection,
const hp::FECollection< dim > &  fe_collection,
const hp::QCollection< dim - 1 > &  q_collection,
const hp::QCollection< 1 > &  q_collection_1d,
const RegionUpdateFlags  region_update_flags,
const MeshClassifier< dim > &  mesh_classifier,
const DoFHandler< dim > &  dof_handler,
const VectorType &  level_set,
const AdditionalData additional_data = AdditionalData() 
)

Constructor.

Parameters
mapping_collectionCollection of Mappings to be used.
fe_collectionCollection of FiniteElements to be used.
q_collectionCollection of Quadrature rules over \([0, 1]^{dim-1}\) that should be used when a face is not intersected and we do not need to generate immersed quadrature rules.
q_collection_1dCollection of 1-dimensional quadrature rules used to generate the immersed quadrature rules. See the QuadratureGenerator class.
mesh_classifierObject used to determine when the immersed quadrature rules need to be generated.
region_update_flagsStruct storing UpdateFlags for the inside/outside region of the cell.
dof_handlerThe DoFHandler associated with the discrete level set function.
level_setThe degrees of freedom of the discrete level set function.
additional_dataAdditional data passed to the QuadratureGenerator class.
Note
Pointers to mapping_collection, fe_collection, mesh_classifier, dof_handler, and level_set are stored internally, so these need to have a longer life span than the instance of this class.

Definition at line 305 of file fe_values.cc.

Member Function Documentation

◆ reinit() [1/2]

template<int dim>
template<class CellIteratorType , class CellNeighborIteratorType >
void NonMatching::FEInterfaceValues< dim >::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 
)

Reinitialize on the shared face between two neighboring cells. After calling this function, an FEInterfaceValues object can be retrieved for each region using the functions get_inside_fe_values() or get_outside_fe_values().

Note
Currently cell and cell_neighbor need to have the same finite element associated with them.
Specifying sub_face_no/sub_face_no_neighbor is currently not implemented, the passed values of these must equal numbers::invalid_unsigned_int.

◆ reinit() [2/2]

template<int dim>
template<class CellIteratorType >
void NonMatching::FEInterfaceValues< dim >::reinit ( const CellIteratorType &  cell,
const unsigned int  face_no 
)

Reinitialize on a face on the boundary of the domain, when there is no neighbor and a single face face_no of the cell. After calling this function, an FEInterfaceValues object can be retrieved for each region using the functions get_inside_fe_values() or get_outside_fe_values().

◆ get_inside_fe_values()

template<int dim>
const std_cxx17::optional<::FEInterfaceValues< dim > > & FEInterfaceValues< dim >::get_inside_fe_values

Return an FEInterfaceValues object reinitialized with a quadrature for the inside region of the cell: \(\{x \in K : \psi(x) < 0 \}\).

Note
If the quadrature rule over the region is empty, e.g. because the cell is completely located in the outside domain, the returned optional will not contain a value.

Definition at line 461 of file fe_values.cc.

◆ get_outside_fe_values()

template<int dim>
const std_cxx17::optional<::FEInterfaceValues< dim > > & FEInterfaceValues< dim >::get_outside_fe_values

Return an FEInterfaceValues object reinitialized with a quadrature for the outside region of the cell: \(\{x \in K : \psi(x) > 0 \}\).

Note
If the quadrature rule over the region is empty, e.g. because the cell is completely located in the inside domain, the returned optional will not contain a value.

Definition at line 473 of file fe_values.cc.

◆ initialize()

template<int dim>
void FEInterfaceValues< dim >::initialize ( const hp::QCollection< dim - 1 > &  q_collection)
private

Do work common to the constructors. The incoming QCollection should be quadratures integrating over \([0, 1]^{dim-1}\). These will be used on the non-intersected cells.

Definition at line 332 of file fe_values.cc.

◆ do_reinit()

template<int dim>
template<bool level_dof_access>
void FEInterfaceValues< dim >::do_reinit ( const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > &  cell,
const unsigned int  face_no,
const std::function< void(::FEInterfaceValues< dim > &)> &  call_reinit 
)
private

Do work that is common two the two reinit functions of the class. call_reinit is a std::function that describes how we should call reinit on a single FEInterfaceValues object, which is what differs between the two reinit functions.

Definition at line 384 of file fe_values.cc.

Member Data Documentation

◆ mapping_collection

template<int dim>
const SmartPointer<const hp::MappingCollection<dim> > NonMatching::FEInterfaceValues< dim >::mapping_collection
private

A pointer to the collection of mappings to be used.

Definition at line 593 of file fe_values.h.

◆ fe_collection

template<int dim>
const SmartPointer<const hp::FECollection<dim> > NonMatching::FEInterfaceValues< dim >::fe_collection
private

A pointer to the collection of finite elements to be used.

Definition at line 598 of file fe_values.h.

◆ q_collection_1D

template<int dim>
const hp::QCollection<1> NonMatching::FEInterfaceValues< dim >::q_collection_1D
private

Collection of 1-dimensional quadrature rules that are used by QuadratureGenerator as base for generating the immersed quadrature rules.

Definition at line 604 of file fe_values.h.

◆ current_face_location

template<int dim>
LocationToLevelSet NonMatching::FEInterfaceValues< dim >::current_face_location
private

Location of the last face that reinit was called with.

Definition at line 609 of file fe_values.h.

◆ active_fe_index

template<int dim>
unsigned int NonMatching::FEInterfaceValues< dim >::active_fe_index
private

Active fe index of the last cell that reinit was called with.

Definition at line 614 of file fe_values.h.

◆ region_update_flags

template<int dim>
const RegionUpdateFlags NonMatching::FEInterfaceValues< dim >::region_update_flags
private

The update flags passed to the constructor.

Definition at line 619 of file fe_values.h.

◆ mesh_classifier

template<int dim>
const SmartPointer<const MeshClassifier<dim> > NonMatching::FEInterfaceValues< dim >::mesh_classifier
private

Pointer to the MeshClassifier passed to the constructor.

Definition at line 624 of file fe_values.h.

◆ fe_values_inside_full_quadrature

template<int dim>
std::deque<std_cxx17::optional<::FEInterfaceValues<dim> > > NonMatching::FEInterfaceValues< dim >::fe_values_inside_full_quadrature
private

For each element in the FECollection passed to the constructor, this object contains an FEInterfaceValues object created with a quadrature rule over the full reference cell: \([0, 1]^{dim-1}\) and UpdateFlags for the inside region. Thus, these optionals should always contain a value.

When LocationToLevelSet of the cell is INSIDE (and we do not need to generate an immersed quadrature), we return the FEInterfaceValues object in this container corresponding to the cell's active_fe_index.

This container is a std::deque, which is compatible with the FEInterfaceValues class that does not have a copy-constructor.

Definition at line 642 of file fe_values.h.

◆ fe_values_outside_full_quadrature

template<int dim>
std::deque<std_cxx17::optional<::FEInterfaceValues<dim> > > NonMatching::FEInterfaceValues< dim >::fe_values_outside_full_quadrature
private

For each element in the FECollection passed to the constructor, this object contains an FEInterfaceValues object created with a quadrature rule over the full reference cell: \([0, 1]^{dim-1}\) and UpdateFlags for the outside region. Thus, these optionals should always contain a value.

When LocationToLevelSet of the cell is OUTSIDE (and we do not need to generate an immersed quadrature), we return the FEValues object in this container corresponding to the cell's active_fe_index.

This container is a std::deque, which is compatible with the FEInterfaceValues class that does not have a copy-constructor.

Definition at line 659 of file fe_values.h.

◆ fe_values_inside

template<int dim>
std_cxx17::optional<::FEInterfaceValues<dim> > NonMatching::FEInterfaceValues< dim >::fe_values_inside
private

FEInterfaceValues object created with a quadrature rule integrating over the inside region, \(\{x \in B : \psi(x) < 0 \}\), that was generated in the last call to reinit(..). If the cell in the last call was not intersected or if 0 quadrature points were generated, this optional will not contain a value.

Definition at line 668 of file fe_values.h.

◆ fe_values_outside

template<int dim>
std_cxx17::optional<::FEInterfaceValues<dim> > NonMatching::FEInterfaceValues< dim >::fe_values_outside
private

FEInterfaceValues object created with a quadrature rule integrating over the outside region, \(\{x \in B : \psi(x) > 0 \}\), that was generated in the last call to reinit(..). If the cell in the last call was not intersected or if 0 quadrature points were generated, this optional will not contain a value.

Definition at line 677 of file fe_values.h.

◆ face_quadrature_generator

template<int dim>
DiscreteFaceQuadratureGenerator<dim> NonMatching::FEInterfaceValues< dim >::face_quadrature_generator
private

Object that generates the immersed quadrature rules.

Definition at line 682 of file fe_values.h.


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