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 | Protected Attributes | Private Member Functions | Private Attributes | List of all members
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 >:
[legend]

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)
 
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)
 

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
 

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 1208 of file quadrature_generator.h.

Constructor & Destructor Documentation

◆ QGenerator()

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 1186 of file quadrature_generator.cc.

Member Function Documentation

◆ generate()

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 1202 of file quadrature_generator.cc.

◆ set_1D_quadrature()

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 1263 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 1234 of file quadrature_generator.cc.

◆ clear_quadratures()

void NonMatching::internal::QuadratureGeneratorImplementation::QGeneratorBase< dim, spacedim >::clear_quadratures
inherited

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

Definition at line 968 of file quadrature_generator.cc.

◆ get_quadratures()

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

Return the created quadratures.

Definition at line 974 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 1249 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 1254 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 1260 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 1266 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 1267 of file quadrature_generator.h.

◆ additional_data

Stores options/settings for the algorithm.

Definition at line 980 of file quadrature_generator.h.

◆ q_index

unsigned int NonMatching::internal::QuadratureGeneratorImplementation::QGeneratorBase< dim, spacedim >::q_index
protectedinherited

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

Definition at line 986 of file quadrature_generator.h.

◆ q_collection1D

const SmartPointer<const hp::QCollection<1> > NonMatching::internal::QuadratureGeneratorImplementation::QGeneratorBase< dim, spacedim >::q_collection1D
protectedinherited

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

Definition at line 992 of file quadrature_generator.h.

◆ q_partitioning

QPartitioning<dim> NonMatching::internal::QuadratureGeneratorImplementation::QGeneratorBase< dim, spacedim >::q_partitioning
protectedinherited

Quadratures that the derived classes create.

Definition at line 997 of file quadrature_generator.h.


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