Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30: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
Classes | Enumerations | Functions
NonMatching::internal::QuadratureGeneratorImplementation Namespace Reference

Classes

class  ExtendableQuadrature
 
struct  FunctionBounds
 
struct  HeightDirectionData
 
class  QGenerator
 
class  QGenerator< 1, spacedim >
 
class  QGeneratorBase
 
class  QPartitioning
 
class  RootFinder
 
class  UpThroughDimensionCreator
 

Enumerations

enum class  Definiteness { negative , positive , indefinite }
 

Functions

template<int dim>
void tensor_point_with_1D_quadrature (const Point< dim - 1 > &point, const double weight, const Quadrature< 1 > &quadrature1D, const double start, const double end, const unsigned int component_in_dim, ExtendableQuadrature< dim > &quadrature)
 
template<int dim>
Definiteness pointwise_definiteness (const std::vector< std::reference_wrapper< const Function< dim > > > &functions, const Point< dim > &point)
 
template<int dim>
std::pair< double, double > find_extreme_values (const std::vector< FunctionBounds< dim > > &all_function_bounds)
 
template<int dim>
std::optional< HeightDirectionDatafind_best_height_direction (const std::vector< FunctionBounds< dim > > &all_function_bounds)
 
template<int dim>
void add_tensor_product (const Quadrature< dim - 1 > &lower, const Quadrature< 1 > &quadrature1D, const double start, const double end, const unsigned int component_in_dim, ExtendableQuadrature< dim > &quadrature)
 
template<int dim>
void take_min_max_at_vertices (const Function< dim > &function, const BoundingBox< dim > &box, std::pair< double, double > &value_bounds)
 
template<int dim>
void estimate_function_bounds (const std::vector< std::reference_wrapper< const Function< dim > > > &functions, const BoundingBox< dim > &box, std::vector< FunctionBounds< dim > > &all_function_bounds)
 
bool is_indefinite (const std::pair< double, double > &function_bounds)
 
double lower_bound_on_abs (const std::pair< double, double > &function_bounds)
 
template<int dim>
bool one_positive_one_negative_definite (const std::vector< FunctionBounds< dim > > &all_function_bounds)
 
template<int dim>
void map_quadrature_to_box (const Quadrature< dim > &unit_quadrature, const BoundingBox< dim > &box, ExtendableQuadrature< dim > &quadrature)
 
template<int dim>
void restrict_to_top_and_bottom (const std::vector< std::reference_wrapper< const Function< dim > > > &functions, const BoundingBox< dim > &box, const unsigned int direction, std::vector< Functions::CoordinateRestriction< dim - 1 > > &restrictions)
 
template<int dim>
void restrict_to_point (const std::vector< std::reference_wrapper< const Function< dim > > > &functions, const Point< dim - 1 > &point, const unsigned int open_direction, std::vector< Functions::PointRestriction< dim - 1 > > &restrictions)
 
template<int dim>
void distribute_points_between_roots (const Quadrature< 1 > &quadrature1D, const BoundingBox< 1 > &interval, const std::vector< double > &roots, const Point< dim - 1 > &point, const double weight, const unsigned int height_function_direction, const std::vector< std::reference_wrapper< const Function< 1 > > > &level_sets, const AdditionalQGeneratorData &additional_data, QPartitioning< dim > &q_partitioning)
 
template<int dim>
Point< dim > face_projection_closest_zero_contour (const Point< dim - 1 > &point, const unsigned int direction, const BoundingBox< dim > &box, const Function< dim > &level_set)
 
template<int dim>
std::optional< unsigned intdirection_of_largest_extent (const BoundingBox< dim > &box)
 
template<int dim>
unsigned int compute_split_direction (const BoundingBox< dim > &box, const std::optional< HeightDirectionData > &height_direction_data)
 
template<int dim>
BoundingBox< dim > left_half (const BoundingBox< dim > &box, const unsigned int direction)
 
template<int dim>
BoundingBox< dim > right_half (const BoundingBox< dim > &box, const unsigned int direction)
 

Enumeration Type Documentation

◆ Definiteness

Type that describes the definiteness of a function over a region.

Enumerator
negative 
positive 
indefinite 

Definition at line 777 of file quadrature_generator.h.

Function Documentation

◆ tensor_point_with_1D_quadrature()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::tensor_point_with_1D_quadrature ( const Point< dim - 1 > &  point,
const double  weight,
const Quadrature< 1 > &  quadrature1D,
const double  start,
const double  end,
const unsigned int  component_in_dim,
ExtendableQuadrature< dim > &  quadrature 
)
private

Take the tensor product between (point, weight) and quadrature1d scaled over [start, end] and add the resulting dim-dimensional quadrature points to quadrature.

component_in_dim specifies which dim-dimensional coordinate quadrature1d should be written to.

Definition at line 54 of file quadrature_generator.cc.

◆ pointwise_definiteness()

template<int dim>
Definiteness NonMatching::internal::QuadratureGeneratorImplementation::pointwise_definiteness ( const std::vector< std::reference_wrapper< const Function< dim > > > &  functions,
const Point< dim > &  point 
)
private

Checks the sign of the incoming Functions at the incoming point and returns Definiteness::positive/negative if all the functions are positive/negative at the point, otherwise returns Definiteness::indefinite.

Definition at line 106 of file quadrature_generator.cc.

◆ find_extreme_values()

template<int dim>
std::pair< double, double > NonMatching::internal::QuadratureGeneratorImplementation::find_extreme_values ( const std::vector< FunctionBounds< dim > > &  all_function_bounds)
private

Returns the max/min bounds on the value, taken over all the entries in the incoming vector of FunctionBounds. That is, given the incoming function bounds, \([L_j, U_j]\), this function returns \([L, U]\), where \(L = \min_{j} L_j\) and \(U = \max_{j} U_j\).

Definition at line 202 of file quadrature_generator.cc.

◆ find_best_height_direction()

template<int dim>
std::optional< HeightDirectionData > NonMatching::internal::QuadratureGeneratorImplementation::find_best_height_direction ( const std::vector< FunctionBounds< dim > > &  all_function_bounds)
private

Finds the best choice of height function direction, given the FunctionBounds for a number of functions \(\{\psi_j\}_{j=0}^{n-1}\). Here, "best" is meant in the sense of the implicit function theorem.

Let \(J_I\) be the index set of the indefinite functions:

\(J_I = \{0,..., n - 1\} \setminus \{ j : |\psi_j| > 0 \}\).

This function converts the incoming bounds to a lower bound, \(L_{ij}\), on the absolute value of each component of the gradient:

\(|\partial_k \psi_j| > L_{jk}\).

and then returns a coordinate direction, \(i\), and a lower bound \(L\), such that

\[ i = \arg \max_{k} \min_{j \in J_I} L_{jk}, \\ L = \max_{k} \min_{j \in J_I} L_{jk}. \]

This means \(i\) is a coordinate direction such that all functions intersected by the zero contour (i.e. those belonging to \(J_I\)) fulfill

\(|\partial_i \psi_j| > L\).

Note that the estimated lower bound, \(L\), can be zero or negative. This means that no suitable height function direction exists. If all of the incoming functions are positive or negative definite the returned std::optional is non-set.

Definition at line 276 of file quadrature_generator.cc.

◆ add_tensor_product()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::add_tensor_product ( const Quadrature< dim - 1 > &  lower,
const Quadrature< 1 > &  quadrature1D,
const double  start,
const double  end,
const unsigned int  component_in_dim,
ExtendableQuadrature< dim > &  quadrature 
)

For each (point, weight) in lower create a dim-dimensional quadrature using tensor_point_with_1D_quadrature and add the results to quadrature.

Definition at line 83 of file quadrature_generator.cc.

◆ take_min_max_at_vertices()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::take_min_max_at_vertices ( const Function< dim > &  function,
const BoundingBox< dim > &  box,
std::pair< double, double > &  value_bounds 
)

Given the incoming lower and upper bounds on the value of a function \([L, U]\), return the minimum/maximum of \([L, U]\) and the function values at the vertices. That is, this function returns

\([\min(L, L_f), \max(U, U_f)]\),

where \(L_f = \min_{v} f(x_v)\), \(U_f = \max_{v} f(x_v)|\), and \(x_v\) is a vertex.

It is assumed that the incoming function is scalar valued.

Definition at line 151 of file quadrature_generator.cc.

◆ estimate_function_bounds()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::estimate_function_bounds ( const std::vector< std::reference_wrapper< const Function< dim > > > &  functions,
const BoundingBox< dim > &  box,
std::vector< FunctionBounds< dim > > &  all_function_bounds 
)

Estimate bounds on each of the functions in the incoming vector over the incoming box.

Bounds on the functions value and the gradient components are first computed using FunctionTools::taylor_estimate_function_bounds. In addition, the function value is checked for min/max at the at the vertices of the box. The gradient is not checked at the box vertices.

Definition at line 179 of file quadrature_generator.cc.

◆ is_indefinite()

bool NonMatching::internal::QuadratureGeneratorImplementation::is_indefinite ( const std::pair< double, double > &  function_bounds)
inline

Return true if the incoming function bounds correspond to a function which is indefinite, i.e., that is not negative or positive definite.

Definition at line 223 of file quadrature_generator.cc.

◆ lower_bound_on_abs()

double NonMatching::internal::QuadratureGeneratorImplementation::lower_bound_on_abs ( const std::pair< double, double > &  function_bounds)
inline

Return a lower bound, \(L_a\), on the absolute value of a function, \(f(x)\):

\(L_a \leq |f(x)|\),

by estimating it from the incoming lower and upper bounds: \(L \leq f(x) \leq U\).

By rewriting the lower and upper bounds as \(F - C \leq f(x) \leq F + C\), where \(L = F - C\), \(U = F + C\) (or \(F = (U + L)/2\), \(C = (U - L)/2\)), we get \(|f(x) - F| \leq C\). Using the inverse triangle inequality gives \(|F| - |f(x)| \leq |f(x) - F| \leq C\). Thus, \(L_a = |F| - C\).

Note that the returned value can be negative. This is used to indicate "how far away" a function is from being definite.

Definition at line 255 of file quadrature_generator.cc.

◆ one_positive_one_negative_definite()

template<int dim>
bool NonMatching::internal::QuadratureGeneratorImplementation::one_positive_one_negative_definite ( const std::vector< FunctionBounds< dim > > &  all_function_bounds)
inline

Return true if there are exactly two incoming FunctionBounds and they corresponds to one function being positive definite and one being negative definite. Return false otherwise.

Definition at line 336 of file quadrature_generator.cc.

◆ map_quadrature_to_box()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::map_quadrature_to_box ( const Quadrature< dim > &  unit_quadrature,
const BoundingBox< dim > &  box,
ExtendableQuadrature< dim > &  quadrature 
)

Transform the points and weights of the incoming quadrature, unit_quadrature, from unit space to the incoming box and add these to quadrature.

Note that unit_quadrature should be a quadrature over [0,1]^dim.

Definition at line 365 of file quadrature_generator.cc.

◆ restrict_to_top_and_bottom()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::restrict_to_top_and_bottom ( const std::vector< std::reference_wrapper< const Function< dim > > > &  functions,
const BoundingBox< dim > &  box,
const unsigned int  direction,
std::vector< Functions::CoordinateRestriction< dim - 1 > > &  restrictions 
)

For each of the incoming dim-dimensional functions, create the restriction to the top and bottom of the incoming BoundingBox and add these two (dim-1)-dimensional functions to restrictions. Here, top and bottom is meant with respect to the incoming direction. For each function, the "bottom-restriction" will be added before the "top-restriction"

Note
restrictions will be cleared, so after this function restrictions.size() == 2 * functions.size().

Definition at line 392 of file quadrature_generator.cc.

◆ restrict_to_point()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::restrict_to_point ( const std::vector< std::reference_wrapper< const Function< dim > > > &  functions,
const Point< dim - 1 > &  point,
const unsigned int  open_direction,
std::vector< Functions::PointRestriction< dim - 1 > > &  restrictions 
)

Restrict each of the incoming functions to point, while keeping the coordinate direction open_direction open, and add the restriction to restrictions.

Note
restrictions will be cleared, so after this function restrictions.size() == functions.size().

Definition at line 428 of file quadrature_generator.cc.

◆ distribute_points_between_roots()

template<int dim>
void NonMatching::internal::QuadratureGeneratorImplementation::distribute_points_between_roots ( const Quadrature< 1 > &  quadrature1D,
const BoundingBox< 1 > &  interval,
const std::vector< double > &  roots,
const Point< dim - 1 > &  point,
const double  weight,
const unsigned int  height_function_direction,
const std::vector< std::reference_wrapper< const Function< 1 > > > &  level_sets,
const AdditionalQGeneratorData additional_data,
QPartitioning< dim > &  q_partitioning 
)

Let \(\{ y_0, ..., y_{n+1} \}\) be such that \([y_0, y_{n+1}]\) is the interval and \(\{ y_1, ..., y_n \}\) are the roots. In each subinterval, \([y_i, y_{i+1}]\), distribute point according to the 1D-quadrature rule \(\{(x_q, w_q)\}_q\) (quadrature1D). Take the tensor product with the quadrature point \((x, w)\) (point, weight) to create dim-dimensional quadrature points

\[ X_q = x_I \times (y_i + (y_{i+1} - y_i) x_q), W_q = w_I (y_{i+1} - y_i) w_q, \]

and add these points to q_partitioning.

Definition at line 463 of file quadrature_generator.cc.

◆ face_projection_closest_zero_contour()

template<int dim>
Point< dim > NonMatching::internal::QuadratureGeneratorImplementation::face_projection_closest_zero_contour ( const Point< dim - 1 > &  point,
const unsigned int  direction,
const BoundingBox< dim > &  box,
const Function< dim > &  level_set 
)

Takes a (dim-1)-dimensional point from the cross-section (orthogonal to direction) of the box. Creates the two dim-dimensional points, which are the projections from the cross section to the faces of the box and returns the point closest to the zero-contour of the incoming level set function.

Definition at line 706 of file quadrature_generator.cc.

◆ direction_of_largest_extent()

template<int dim>
std::optional< unsigned int > NonMatching::internal::QuadratureGeneratorImplementation::direction_of_largest_extent ( const BoundingBox< dim > &  box)

Return the coordinate direction of the largest side of the box. If two or more sides have the same length the returned std::optional will be non-set.

Definition at line 982 of file quadrature_generator.cc.

◆ compute_split_direction()

template<int dim>
unsigned int NonMatching::internal::QuadratureGeneratorImplementation::compute_split_direction ( const BoundingBox< dim > &  box,
const std::optional< HeightDirectionData > &  height_direction_data 
)

Return the coordinate direction that the box should be split in, assuming that the box should be split it half.

If the box is larger in one coordante direction, this direction is returned. If the box have the same extent in all directions, we choose the coordinate direction which is closest to being a height-function direction. That is, the direction \(i\) that has a least negative estimate of \(|\partial_i \psi_j|\). As a last resort, we choose the direction 0, if height_direction_data non-set.

Definition at line 1018 of file quadrature_generator.cc.

◆ left_half()

template<int dim>
BoundingBox< dim > NonMatching::internal::QuadratureGeneratorImplementation::left_half ( const BoundingBox< dim > &  box,
const unsigned int  direction 
)
inline

Split the incoming box in half with respect to the incoming coordinate direction and return the left half.

Definition at line 1045 of file quadrature_generator.cc.

◆ right_half()

template<int dim>
BoundingBox< dim > NonMatching::internal::QuadratureGeneratorImplementation::right_half ( const BoundingBox< dim > &  box,
const unsigned int  direction 
)
inline

Split the incoming box in half with respect to the incoming coordinate direction and return the right half.

Definition at line 1064 of file quadrature_generator.cc.