Reference documentation for deal.II version 9.1.0-pre
Public Member Functions | Protected Attributes | List of all members
NonMatching::ImmersedSurfaceQuadrature< dim > Class Template Reference

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

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

Public Member Functions

 ImmersedSurfaceQuadrature ()=default
 
 ImmersedSurfaceQuadrature (const std::vector< Point< dim >> &points, const std::vector< double > &weights, const std::vector< Tensor< 1, dim >> &normals)
 
void push_back (const Point< dim > &point, const double weight, const Tensor< 1, dim > &normal)
 
const Tensor< 1, dim > & normal_vector (const unsigned int i) const
 
const std::vector< Tensor< 1, dim > > & get_normal_vectors () const
 
- Public Member Functions inherited from Quadrature< dim >
 Quadrature (const unsigned int n_quadrature_points=0)
 
 Quadrature (const SubQuadrature &, const Quadrature< 1 > &)
 
 Quadrature (const Quadrature< dim!=1?1:0 > &quadrature_1d)
 
 Quadrature (const Quadrature< dim > &q)
 
 Quadrature (Quadrature< dim > &&) noexcept=default
 
 Quadrature (const std::vector< Point< dim >> &points, const std::vector< double > &weights)
 
 Quadrature (const std::vector< Point< dim >> &points)
 
 Quadrature (const Point< dim > &point)
 
virtual ~Quadrature () override=default
 
Quadratureoperator= (const Quadrature< dim > &)
 
Quadratureoperator= (Quadrature< dim > &&)=default
 
bool operator== (const Quadrature< dim > &p) const
 
void initialize (const std::vector< Point< dim >> &points, const std::vector< double > &weights)
 
unsigned int size () const
 
const Point< dim > & point (const unsigned int i) const
 
const std::vector< Point< dim > > & get_points () const
 
double weight (const unsigned int i) const
 
const std::vector< double > & get_weights () const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
bool is_tensor_product () const
 
std::conditional< dim==1, std::array< Quadrature< 1 >, dim >, const std::array< Quadrature< 1 >, dim > & >::type get_tensor_basis () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Protected Attributes

std::vector< Tensor< 1, dim > > normals
 
- Protected Attributes inherited from Quadrature< dim >
std::vector< Point< dim > > quadrature_points
 
std::vector< double > weights
 
bool is_tensor_product_flag
 
std::unique_ptr< std::array< Quadrature< 1 >, dim > > tensor_basis
 

Additional Inherited Members

- Public Types inherited from Quadrature< dim >
using SubQuadrature = Quadrature< dim-1 >
 
- Static Public Member Functions inherited from Subscriptor
static::ExceptionBase & ExcInUse (int arg1, std::string arg2, std::string arg3)
 
static::ExceptionBase & ExcNoSubscriber (std::string arg1, std::string arg2)
 

Detailed Description

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

Defines a quadrature formula for integration over an oriented surface, \(\hat{S}\), immersed in the unit cell. By immersed it is meant that the surface may intersect the unit cell in an arbitrary way. The quadrature formula is described by a set of quadrature points, \(\hat{x}_q\), weights, \(w_q\), and normalized surface normals, \(\hat{n}_q\).

We typically want to compute surface 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 a surface integral over \(S\cap K\) in real space can be transformed to a surface integral over \(\hat{S} \cap \hat{K}\) according to

\[ \int_{S\cap K} f(x) dS = \int_{S\cap K} f(x) |d\bar{S}| = \int_{\hat{S}\cap\hat{K}} f(F_{K}(\hat{x})) \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. This transformation is possible since the continuous surface elements are vectors: \(d\bar{S}, d\hat{S} \in \mathbb{R}^{dim}\) which are parallel to the normals of \(S\) and \(\hat{S}\). So in order to compute the integral in real space one needs information about the normal 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 := 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(x) 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. \]

immersed_surface_quadrature.svg
Author
Simon Sticko, 2017

Definition at line 75 of file immersed_surface_quadrature.h.

Constructor & Destructor Documentation

Default constructor to initialize the quadrature with no quadrature points.

template<int dim>
NonMatching::ImmersedSurfaceQuadrature< dim >::ImmersedSurfaceQuadrature ( const std::vector< Point< dim >> &  points,
const std::vector< double > &  weights,
const std::vector< Tensor< 1, dim >> &  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 22 of file immersed_surface_quadrature.cc.

Member Function Documentation

template<int dim>
void NonMatching::ImmersedSurfaceQuadrature< dim >::push_back ( const Point< dim > &  point,
const double  weight,
const Tensor< 1, dim > &  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 43 of file immersed_surface_quadrature.cc.

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

Return a reference to the ith surface normal.

Definition at line 58 of file immersed_surface_quadrature.cc.

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

Return a reference to the whole vector of normals.

Definition at line 68 of file immersed_surface_quadrature.cc.

Member Data Documentation

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

Vector of surface normals at each quadrature point.

Definition at line 129 of file immersed_surface_quadrature.h.


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