Reference documentation for deal.II version 9.6.2
 
\(\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
NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim > Class Template Reference

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

Inheritance diagram for NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >:

Public Member Functions

 QGenerator (const hp::QCollection< 1 > &quadratures1D, const AdditionalQGeneratorData &additional_data)
 
void generate (const std::vector< std::reference_wrapper< const Function< 1 > > > &level_sets, const BoundingBox< 1 > &box, const unsigned int n_box_splits)
 
void set_1D_quadrature (const unsigned int q_index)
 
 QGenerator (const hp::QCollection< 1 > &q_collection1D, const AdditionalQGeneratorData &additional_data)
 
void generate (const std::vector< std::reference_wrapper< const Function< dim > > > &level_sets, const BoundingBox< dim > &box, const unsigned int n_box_splits)
 
void set_1D_quadrature (const unsigned int q_index)
 
void clear_quadratures ()
 
const QPartitioning< dim > & get_quadratures () const
 

Protected Attributes

const AdditionalQGeneratorData additional_data
 
unsigned int q_index
 
const SmartPointer< const hp::QCollection< 1 > > q_collection1D
 
QPartitioning< dim > q_partitioning
 

Private Member Functions

void create_surface_points (const std::vector< std::reference_wrapper< const Function< 1 > > > &level_sets)
 
void create_low_dim_quadratures (const unsigned int height_function_direction, const std::vector< std::reference_wrapper< const Function< dim > > > &level_sets, const BoundingBox< dim > &box, const unsigned int n_box_splits)
 
void create_high_dim_quadratures (const unsigned int height_function_direction, const std::vector< std::reference_wrapper< const Function< dim > > > &level_sets, const BoundingBox< dim > &box)
 
void split_box_and_recurse (const std::vector< std::reference_wrapper< const Function< dim > > > &level_sets, const BoundingBox< dim > &box, const std::optional< HeightDirectionData > &direction_data, const unsigned int n_box_splits)
 
void use_midpoint_method (const std::vector< std::reference_wrapper< const Function< dim > > > &level_sets, const BoundingBox< dim > &box)
 

Private Attributes

RootFinder root_finder
 
std::vector< double > roots
 
const unsigned int direction = 0
 
const Point< 0 > zero_dim_point
 
const double unit_weight = 1
 
QGenerator< dim - 1, spacedim > low_dim_algorithm
 
UpThroughDimensionCreator< dim, spacedim > up_through_dimension_creator
 
hp::QCollection< dim > tensor_products
 

Detailed Description

template<int spacedim>
class NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >

The 1d-base case of the recursive algorithm QGenerator<dim, spacedim>.

Let $L$ and $R$ be the left and right bounds of the one-dimensional BoundingBox. This interval is partitioned into $[x_0, x_1, ..., x_n]$ where $x_0 = L$, $x_n = R$, and the remaining $x_i$ are the roots of the level set functions in the interval $[L, R]$. In each interval, $[x_i, x_{i+1}]$, quadrature points are distributed according to a 1d-quadrature rule. These points are added to one of the regions of QPartitioning determined from the signs of the level set functions on the interval (see documentation of QPartitioning).

If spacedim = 1 the points $[x_1, x_{n-1}]$ are also added as surface quadrature points to QPartitioning::surface.

Definition at line 1276 of file quadrature_generator.h.

Constructor & Destructor Documentation

◆ QGenerator() [1/2]

template<int spacedim>
NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::QGenerator ( const hp::QCollection< 1 > & quadratures1D,
const AdditionalQGeneratorData & additional_data )

Constructor. Takes the same parameters QuadratureGenerator.

Definition at line 1209 of file quadrature_generator.cc.

◆ QGenerator() [2/2]

NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::QGenerator ( const hp::QCollection< 1 > & q_collection1D,
const AdditionalQGeneratorData & additional_data )

Constructor. Takes the same parameters QuadratureGenerator.

Definition at line 1159 of file quadrature_generator.cc.

Member Function Documentation

◆ generate() [1/2]

template<int spacedim>
void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::generate ( const std::vector< std::reference_wrapper< const Function< 1 > > > & level_sets,
const BoundingBox< 1 > & box,
const unsigned int n_box_splits )

Creates quadrature points over the interval defined by the incoming box and adds these quadrature points to the internally stored QPartitioning. These quadratures can then be obtained using the get_quadratures-function.

Definition at line 1225 of file quadrature_generator.cc.

◆ set_1D_quadrature() [1/2]

template<int spacedim>
void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::set_1D_quadrature ( const unsigned int q_index)

Set which 1d-quadrature in the collection passed to the constructor should be used to create the immersed quadratures.

Definition at line 1286 of file quadrature_generator.cc.

◆ create_surface_points()

template<int spacedim>
void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::create_surface_points ( const std::vector< std::reference_wrapper< const Function< 1 > > > & level_sets)
private

Adds the point defined by coordinate to the surface quadrature of ImmersedQuadrature with unit weight.

Definition at line 1257 of file quadrature_generator.cc.

◆ generate() [2/2]

void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::generate ( const std::vector< std::reference_wrapper< const Function< dim > > > & level_sets,
const BoundingBox< dim > & box,
const unsigned int n_box_splits )

Create immersed quadrature rules over the incoming box and add these to the internal QPartitioning<dim> object in the base class. These quadratures can then be obtained using the get_quadratures-function.

This function calls itself if the incoming box need to be split. n_box_splits counts the number of times this function has called itself.

Definition at line 1173 of file quadrature_generator.cc.

◆ set_1D_quadrature() [2/2]

void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::set_1D_quadrature ( const unsigned int q_index)

Set which 1d-quadrature in the collection passed to the constructor should be used to create the immersed quadratures.

Definition at line 1183 of file quadrature_generator.cc.

◆ create_low_dim_quadratures()

void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::create_low_dim_quadratures ( const unsigned int height_function_direction,
const std::vector< std::reference_wrapper< const Function< dim > > > & level_sets,
const BoundingBox< dim > & box,
const unsigned int n_box_splits )
private

Restricts the incoming level set functions to the top and bottom of the incoming box (w.r.t height_function_direction). Then call the lower dimensional QGenerator with the cross section of the box to generate the lower dimensional immersed quadrature rules.

Definition at line 1193 of file quadrature_generator.cc.

◆ create_high_dim_quadratures()

void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::create_high_dim_quadratures ( const unsigned int height_function_direction,
const std::vector< std::reference_wrapper< const Function< dim > > > & level_sets,
const BoundingBox< dim > & box )
private

Gets the $(dim - 1)$-dimensional quadratures from the lower dimensional algorithm and creates the $dim$-dimensional quadrature rules over the box from the lower dimensional ones.

Definition at line 1206 of file quadrature_generator.cc.

◆ split_box_and_recurse()

void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::split_box_and_recurse ( const std::vector< std::reference_wrapper< const Function< dim > > > & level_sets,
const BoundingBox< dim > & box,
const std::optional< HeightDirectionData > & direction_data,
const unsigned int n_box_splits )
private

Split the incoming box and call generate() recursively with each box. The box is split in 2 or 4 parts depending on the value of AdditionalQGeneratorData::split_in_half.

Definition at line 1218 of file quadrature_generator.cc.

◆ use_midpoint_method()

void NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::use_midpoint_method ( const std::vector< std::reference_wrapper< const Function< dim > > > & level_sets,
const BoundingBox< dim > & box )
private

Uses the midpoint-method to create a quadrature over the box. That is, add a single quadrature point at the center of the box with weight corresponding to the volume of the box.

The point is added to the region defined in QPartitioning according to the signs of the level set functions at the center of the box.

Definition at line 1235 of file quadrature_generator.cc.

◆ clear_quadratures()

Clear the quadratures created by the previous call to generate().

Definition at line 1036 of file quadrature_generator.cc.

◆ get_quadratures()

const QPartitioning< dim > & NonMatching::internal::QuadratureGeneratorImplementation::QGeneratorBase< dim, spacedim >::get_quadratures ( ) const

Return the created quadratures.

Definition at line 1042 of file quadrature_generator.cc.

Member Data Documentation

◆ root_finder

template<int spacedim>
RootFinder NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::root_finder
private

Class used to find the roots of the functions passed to generate().

Definition at line 1317 of file quadrature_generator.h.

◆ roots

template<int spacedim>
std::vector<double> NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::roots
private

Roots of the functions passed to generate().

Definition at line 1322 of file quadrature_generator.h.

◆ direction

template<int spacedim>
const unsigned int NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::direction = 0
private

This would be the height-function direction in higher dimensions, but in 1d there is only one coordinate direction.

Definition at line 1328 of file quadrature_generator.h.

◆ zero_dim_point

template<int spacedim>
const Point<0> NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::zero_dim_point
private

To reuse the distribute_points_between_roots()-function we need a zero-dimensional quadrature point with unit weight.

Definition at line 1334 of file quadrature_generator.h.

◆ unit_weight

template<int spacedim>
const double NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< 1, spacedim >::unit_weight = 1
private

Definition at line 1335 of file quadrature_generator.h.

◆ low_dim_algorithm

QGenerator<dim - 1, spacedim> NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::low_dim_algorithm
private

The same algorithm as this, but creating immersed quadratures in one dimension lower.

Definition at line 1244 of file quadrature_generator.h.

◆ up_through_dimension_creator

UpThroughDimensionCreator<dim, spacedim> NonMatching::internal::QuadratureGeneratorImplementation::QGenerator< dim, spacedim >::up_through_dimension_creator
private

Object responsible for creating the $dim$-dimensional quadratures from

Definition at line 1250 of file quadrature_generator.h.

◆ tensor_products

Stores tensor products of each of the Quadrature<1>'s in q_collection1d.

Definition at line 1256 of file quadrature_generator.h.

◆ additional_data

Stores options/settings for the algorithm.

Definition at line 1048 of file quadrature_generator.h.

◆ q_index

Which 1d-quadrature in the collection we should use to generate the immersed quadrature.

Definition at line 1054 of file quadrature_generator.h.

◆ q_collection1D

Index of the quadrature in q_collection1d that should use to generate the immersed quadrature rules.

Definition at line 1060 of file quadrature_generator.h.

◆ q_partitioning

Quadratures that the derived classes create.

Definition at line 1065 of file quadrature_generator.h.


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