Reference documentation for deal.II version GIT relicensing-478-g3275795167 2024-04-24 07:10:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | List of all members
NonMatching::ImmersedSurfaceQuadrature< dim, spacedim > Class Template Reference

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

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

Public Member Functions

 ImmersedSurfaceQuadrature ()=default
 
 ImmersedSurfaceQuadrature (const std::vector< Point< dim > > &points, const std::vector< double > &weights, const std::vector< Tensor< 1, spacedim > > &normals)
 
void clear ()
 
void push_back (const Point< dim > &point, const double weight, const Tensor< 1, spacedim > &normal)
 
const Tensor< 1, spacedim > & normal_vector (const unsigned int i) const
 
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors () const
 

Protected Attributes

std::vector< Tensor< 1, spacedim > > normals
 

Detailed Description

template<int dim, int spacedim = dim>
class NonMatching::ImmersedSurfaceQuadrature< dim, spacedim >

This class defines a quadrature formula to integrate over the intersection between an oriented surface, \(\hat{S}\), and a cell or face. The word "immersed" in the class name reflects that the surface may intersect the cell/face in an arbitrary way.

The spacedim template parameter of this class is the dimension that the (spacedim-1)-dimensional surface is embedded in: \(\hat{S} \subset \mathbb{R}^{\text{spacedim}}\). The dim parameter describes the dimension of the "object" that the surface intersects. That is, dim = spacedim corresponds to the surface intersecting a cell and dim = spacedim - 1 corresponds to the surface intersecting a face. The quadrature formula is described by a set of quadrature points, \(\hat{x}_q \in \mathbb{R}^{\text{dim}}\), weights, \(w_q\), and normalized surface normals, \(\hat{n}_q \in \mathbb{R}^{\text{spacedim}}\).

Consider first the case dim = spacedim. We typically want to compute integrals in real space. A surface, \(S\), intersecting a cell, \(K\), in real space can be mapped onto a surface, \(\hat{S}\), intersecting the unit cell, \(\hat{K}\). Thus an integral over \(S\cap K\) in real space can be transformed to an integral over \(\hat{S} \cap \hat{K}\) according to

\[ \int_{S\cap K} f dS = \int_{S\cap K} f |d\bar{S}| = \int_{\hat{S}\cap\hat{K}} f \circ F_{K} \det(J) |\left( J^{-1} \right )^T d\hat{S}|, \]

where \(F_K\) is the mapping from reference to real space and \(J\) is its Jacobian matrix. This transformation is possible since the continuous surface elements are vectors: \(d\bar{S}, d\hat{S} \in \mathbb{R}^{spacedim}\), which are parallel to the normals of \(S\) and \(\hat{S}\). That is, the normal is needed to do the transformation. Thus, in addition to storing points and weights, this quadrature stores also the normalized normal for each quadrature point. This can be viewed as storing a discrete surface element,

\[ \Delta \hat{S}_q \dealcoloneq w_q \hat{n}_q \approx d\hat{S}(\hat{x}_q), \]

for each quadrature point. The surface integral in real space would then be approximated as

\[ \int_{S\cap K} f dS \approx \sum_{q} f \left(F_{K}(\hat{x}_{q}) \right) \det(J_q) |\left( J_q^{-1} \right)^T \hat{n}_q| w_q. \]

When dim = spacedim - 1, this class represents a (spacedim-2)-dimensional integral. That is, if spacedim = 3 we have a line integral immersed in a face. Let \(\hat{r}(t)\), \(t \in [0,T]\) be an arc-length parameterizations of \(\hat{F}\cap \hat{S}\), i.e., the part of the surface that intersects the face in reference space. This means that \(\bar{r}(t) = F_K(\hat{r}(t))\) is a parameterization of \(S\cap F\). The transformation of the line integral now reads

\[ \int_{S\cap F} f dr = \int_{0}^T f(\bar{r}(t)) \left \|\frac{d\bar{r}}{dt} \right \| dt = \int_{0}^T f(F_K(\hat{r}(t))) \left \| J \frac{d\hat{r}}{dt} \right \| dt \approx \sum_{q} f \left(F_{K}(\hat{x}_{q}) \right) \|J(\hat{x}_q) \hat{t}_q \| w_q, \]

where \(\hat{t}_q = \frac{d\hat{r}}{dt}(x_q) \) is the tangent to the curve at \(\hat{x}_q\). This tangent can also be computed as \(t_q = \hat{n}_q \times \hat{n}_F / \| \hat{n}_q \times \hat{n}_F \|\) where \(\hat{n}_F\) is the face normal. It would be possible to compute the tangent by only knowing the normal to the curve in the face plane (i.e. the dim-dimensional normal). However, when these quadratures are used, the weak form typically involves the so-called conormal, which can not be computed without knowing the surface normal in \(\mathbb{R}^{\text{spacedim}}\). The conormal is the unit vector parallel to the projection of the face normal into the surface plane. This is the same as the normalized boundary form.

Definition at line 106 of file immersed_surface_quadrature.h.

Constructor & Destructor Documentation

◆ ImmersedSurfaceQuadrature() [1/2]

template<int dim, int spacedim = dim>
NonMatching::ImmersedSurfaceQuadrature< dim, spacedim >::ImmersedSurfaceQuadrature ( )
default

Default constructor to initialize the quadrature with no quadrature points.

◆ ImmersedSurfaceQuadrature() [2/2]

template<int dim, int spacedim>
NonMatching::ImmersedSurfaceQuadrature< dim, spacedim >::ImmersedSurfaceQuadrature ( const std::vector< Point< dim > > &  points,
const std::vector< double > &  weights,
const std::vector< Tensor< 1, spacedim > > &  normals 
)

Construct a quadrature formula from vectors of points, weights and surface normals. The points, weights and normals should be with respect to reference space, and the normals should be normalized.

Definition at line 21 of file immersed_surface_quadrature.cc.

Member Function Documentation

◆ clear()

template<int dim, int spacedim>
void NonMatching::ImmersedSurfaceQuadrature< dim, spacedim >::clear ( )
inline

Clears weights, points and normals vectors.

Definition at line 42 of file immersed_surface_quadrature.cc.

◆ push_back()

template<int dim, int spacedim>
void NonMatching::ImmersedSurfaceQuadrature< dim, spacedim >::push_back ( const Point< dim > &  point,
const double  weight,
const Tensor< 1, spacedim > &  normal 
)

Extend the given formula by an additional quadrature point. The point, weight and normal should be with respect to reference space, and the normal should be normalized.

This function exists since immersed quadrature rules can be rather complicated to construct. Often the construction is done by partitioning the cell into regions and constructing points on each region separately. This can make it cumbersome to create the quadrature from the constructor since all quadrature points have to be known at time of creation of the object.

Note
This function should only be used during construction of the quadrature formula.

Definition at line 53 of file immersed_surface_quadrature.cc.

◆ normal_vector()

template<int dim, int spacedim>
const Tensor< 1, spacedim > & NonMatching::ImmersedSurfaceQuadrature< dim, spacedim >::normal_vector ( const unsigned int  i) const

Return a reference to the ith surface normal.

Definition at line 69 of file immersed_surface_quadrature.cc.

◆ get_normal_vectors()

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

Return a reference to the whole vector of normals.

Definition at line 80 of file immersed_surface_quadrature.cc.

Member Data Documentation

◆ normals

template<int dim, int spacedim = dim>
std::vector<Tensor<1, spacedim> > NonMatching::ImmersedSurfaceQuadrature< dim, spacedim >::normals
protected

Vector of surface normals at each quadrature point.

Definition at line 166 of file immersed_surface_quadrature.h.


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