KASKADE 7 development version
|
Miscellaneous support routines and classes. More...
Modules | |
Combinatorics | |
Functions for counting and enumerating combinatorial structures. | |
Multithreading | |
Support routines and data structures for multithreaded execution. | |
Exceptions | |
A hierarchy of exception classes providing more or less detailed information. | |
Classes | |
struct | Kaskade::Creator< Product > |
Abstract base class of creators that can be registered at the pluggable factory. More... | |
class | Kaskade::Factory< K, P > |
A pluggable factory. More... | |
class | Kaskade::Deprecated |
Convenience class for marking deprecated functions. A static object of this class can be created inside deprecated functions. It will output a warning message to stderr if the deprecated function is called, but only once. More... | |
struct | Kaskade::ConvertTo< T, Real > |
Reports the converted type. More... | |
class | Kaskade::Timings |
Supports gathering and reporting execution times information for nested program parts. More... | |
class | Kaskade::ScopedTimingSection |
A scope guard object that automatically closes a timing section on destruction. More... | |
class | Kaskade::TaskTiming |
A class that gathers data on task timing and provides gnuplot visualization. More... | |
Functions | |
template<class LocalCoordinate > | |
LocalCoordinate::field_type | checkInside (Dune::GeometryType const >, LocalCoordinate const &x) |
Checks whether a point is inside or outside the reference element, and how much. More... | |
template<class Cell > | |
std::pair< Dune::FieldVector< double, Cell::dimension >, Dune::FieldVector< double, Cell::dimension > > | boundingBox (Cell const &cell) |
Computes the bounding box of a cell. More... | |
template<int d, class Real > | |
Real | Kaskade::bezier (int p, int i, Dune::FieldVector< Real, d > const &x) |
Computes the i-th Bezier function of order p at the given point in the unit simplex. More... | |
constexpr int | Kaskade::binomial (int n, int k) |
Computes the binomial coefficient \( \binom{n}{k} \). More... | |
template<size_t m> | |
size_t | Kaskade::multinomial (std::array< size_t, m > const &ks) |
Computes the multinomial coefficient \( \binom{n}{k_1,\dots,k_m}\), where \( n = \sum_{i=1}^m k_i \). More... | |
constexpr int | Kaskade::numberOfMultiindices (int d, int m) |
Computes the number of multiindices of order m and dimension d. More... | |
template<int d> | |
std::array< size_t, d > | Kaskade::multiIndex (size_t m, size_t n) |
Random access to multiindices. More... | |
template<class It > | |
void | Kaskade::next_multiindex (It first, It last, int const m) |
Computes the next nonnegative multiindex of order m. More... | |
std::vector< double > | Kaskade::lobattoNodes (int n) |
Computes the n Lobatto nodes on [-1,1] including the interval end points. More... | |
std::vector< double > | Kaskade::radauNodes (int n) |
Computes the n Radau nodes on [-1,1] including the right end point. More... | |
std::vector< double > | Kaskade::gaussNodes (int n) |
Computes the n Gauß nodes on [-1,1] excluding the interval end points. More... | |
bool | Kaskade::getKaskadeOptions (int argc, char **argv, Options const &options) |
Supports the extraction of program option values from command line arguments. More... | |
template<class A > | |
auto | Kaskade::moveUnique (A &&a) |
Creates a dynamically allocated object from an r-value reference. More... | |
template<class A > | |
A | Kaskade::duplicate (A const &a) |
Creates a copy of a given object. More... | |
template<class R > | |
R | Kaskade::power (R base, int exp) |
Computes integral powers of values in a multiplicative half-group. More... | |
template<class R > | |
R | Kaskade::square (R x) |
returns the square of its argument. More... | |
template<class R > | |
R | Kaskade::faculty (int n) |
Computes faculties. Is this somewhere in the standard? no. More... | |
constexpr int | Kaskade::sqrti (int n) |
Computes floor(sqrt(n)) for integers n at compile time. More... | |
template<class T > | |
T & | Kaskade::removeConst (T const &t) |
A convenience template for removing the const qualifier from references and pointers. More... | |
Miscellaneous support routines and classes.
Real Kaskade::bezier | ( | int | p, |
int | i, | ||
Dune::FieldVector< Real, d > const & | x | ||
) |
Computes the i-th Bezier function of order p at the given point in the unit simplex.
The Bezier functions are the summands in
\[ 1 = 1^p = \left(\sum_{j=0}^d \lambda_j\right)^p, \]
where \( \lambda_j \) are the barycentric coordinates.
p | the order of Bezier functions |
i | the number of the Bezier function of order p. 0 <= i < numberOfMultiindices(d,p) has to hold. |
x | the simplex coordinates, i.e. \( \lambda = (x,1-\|x\|_1) \). |
|
constexpr |
Computes the binomial coefficient \( \binom{n}{k} \).
Precondition: n >= k >= 0
Definition at line 29 of file combinatorics.hh.
Referenced by Kaskade::bezier(), Kaskade::HierarchicExtensionSimplexShapeFunctionSet< ctype, dimension, Scalar >::HierarchicExtensionSimplexShapeFunctionSet(), and Kaskade::HierarchicSimplexShapeFunctionSet< ctype, dimension, Scalar >::HierarchicSimplexShapeFunctionSet().
std::pair< Dune::FieldVector< double, Cell::dimension >, Dune::FieldVector< double, Cell::dimension > > boundingBox | ( | Cell const & | cell | ) |
Computes the bounding box of a cell.
Cell | the type of cell considered. Usually a typename GridView::template Codim<0>::Entity. |
Definition at line 766 of file fixdune.hh.
Referenced by GeometricObject::ImplementationDetail::DistanceImpl< Metric, Scalar, dim, Ball< Scalar, dim >, BoundingBox< Scalar, dim > >::apply(), GeometricObject::ImplementationDetail::IntersectionCheckImpl< Metric, Scalar, dim, Ball< Scalar, dim >, BoundingBox< Scalar, dim > >::apply(), GeometricObject::BoundingBox< Scalar, dimension >::BoundingBox(), BisectionPolicy< Cell, BoundingBox >::createChildren(), GeometricObject::distance(), GeometricObject::BoundingBoxWrapper< FastBoundingBox< Scalar, dim > >::init(), GeometricObject::BoundingBoxWrapper< FastBoundingBox< Scalar, dim > >::initCorners(), GeometricObject::BoundingBoxWrapper< FastBoundingBox< Scalar, dim > >::initEdges(), GeometricObject::BoundingBoxWrapper< FastBoundingBox< Scalar, dim > >::initFaces(), and GeometricObject::intersects().
LocalCoordinate::field_type checkInside | ( | Dune::GeometryType const & | gt, |
LocalCoordinate const & | x | ||
) |
Checks whether a point is inside or outside the reference element, and how much.
Computational geometry predicates are notoriously difficult to get right. In Dune, it may happen that a point inside a cell is reported not to be contained in any of its children. Often this is due to rounding errors. This function returns not just a binary decision subject to rounding errors, but a floating point value that gives the magnitude of being inside or outside, such that the result can be tested by arbitrary tolerances or sorted according to magnitude.
LocalCoordinate | the type of local coordinate vectors. Usually a Dune::FieldVector<ctype,dim>. |
Definition at line 741 of file fixdune.hh.
Referenced by Kaskade::findCell(), and Kaskade::findCellOnLevel().
A Kaskade::duplicate | ( | A const & | a | ) |
Creates a copy of a given object.
Use this if you need a temporary copy that can be moved from.
A | needs to be copy constructible and (preferably) moveable. |
Definition at line 48 of file memory.hh.
Referenced by Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeJacobiMultiGrid().
R Kaskade::faculty | ( | int | n | ) |
std::vector< double > Kaskade::gaussNodes | ( | int | n | ) |
Computes the n Gauß nodes on [-1,1] excluding the interval end points.
bool Kaskade::getKaskadeOptions | ( | int | argc, |
char ** | argv, | ||
Options const & | options | ||
) |
Supports the extraction of program option values from command line arguments.
Based on boost::program_options, command line arguments (or their default values) can be stored in variables as follows:
If an argument "--help" is given, the list of allowed options is printed.
std::vector< double > Kaskade::lobattoNodes | ( | int | n | ) |
Computes the n Lobatto nodes on [-1,1] including the interval end points.
auto Kaskade::moveUnique | ( | A && | a | ) |
Creates a dynamically allocated object from an r-value reference.
Use this if you get an object of a runtime-polymorphic class hierarchy only as an r-value, e.g. from a makeSomething() function, and you don't care for its actual type:
A | needs to be move (preferably) or copy constructible. |
Definition at line 33 of file memory.hh.
Referenced by Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeBlockJacobiPMultiGrid(), Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeJacobiMultiGrid(), Kaskade::makeMultigrid(), Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeMultiplicativePMultiGrid(), and Kaskade::makePMultiGridStack().
std::array< size_t, d > Kaskade::multiIndex | ( | size_t | m, |
size_t | n | ||
) |
Random access to multiindices.
This yields multiindices of length d that add up to at most m:
\[ \sum_{i=0}^{d-1} k_i \le m \]
d | the number of indices |
m | the maximum sum of the indices |
n | compute the n-th index (0 <= n < numberOfMultiindices(d,m) has to hold). |
Definition at line 118 of file combinatorics.hh.
Referenced by Kaskade::multiIndex().
size_t Kaskade::multinomial | ( | std::array< size_t, m > const & | ks | ) |
Computes the multinomial coefficient \( \binom{n}{k_1,\dots,k_m}\), where \( n = \sum_{i=1}^m k_i \).
m | the number of multinomials |
ks | an array of nonegative integers |
The multinomial coefficient is
\[ \binom{n}{k_1,\dots,k_m} = \frac{n!}{k_1! \cdots k_m!}. \]
Definition at line 61 of file combinatorics.hh.
Referenced by Kaskade::bezier().
void Kaskade::next_multiindex | ( | It | first, |
It | last, | ||
int const | m | ||
) |
Computes the next nonnegative multiindex of order m.
Using this function one can cycle through all integer tuples that sum up to at most m.
Definition at line 157 of file combinatorics.hh.
Referenced by Kaskade::HierarchicExtensionSimplexShapeFunctionSet< ctype, dimension, Scalar >::HierarchicExtensionSimplexShapeFunctionSet(), and Kaskade::HierarchicSimplexShapeFunctionSet< ctype, dimension, Scalar >::HierarchicSimplexShapeFunctionSet().
|
constexpr |
Computes the number of multiindices of order m and dimension d.
A multiindex of dimension d is a d-tuple \( (i_1,\dots,i_d)\) of nonnegative integers that sums up to at most m:
\[ \sum_{j=1}^d i_j \le m \]
Multiindices correspond directly to points on a cartesian grid covering the d-dimensional unit simplex and to Bezier functions.
Definition at line 89 of file combinatorics.hh.
Referenced by Kaskade::multiIndex(), and Kaskade::numberOfMultiindices().
R Kaskade::power | ( | R | base, |
int | exp | ||
) |
Computes integral powers of values in a multiplicative half-group.
Definition at line 30 of file power.hh.
Referenced by Kaskade::bezier(), Kaskade::coarsening(), Kaskade::Elastomechanics::MaterialLaws::MooneyRivlin< dimension >::d2(), estimateGeometricTruncationError(), and Kaskade::power().
std::vector< double > Kaskade::radauNodes | ( | int | n | ) |
Computes the n Radau nodes on [-1,1] including the right end point.
T & Kaskade::removeConst | ( | T const & | t | ) |
A convenience template for removing the const qualifier from references and pointers.
Use this instead of
in polymorphic lambda functions where the type T is not named (just inferred from auto) and hence more involved to obtain, e.g. using decltype.
Definition at line 30 of file typeTraits.hh.
|
constexpr |
Computes floor(sqrt(n)) for integers n at compile time.
Definition at line 88 of file power.hh.
Referenced by Kaskade::sqrti().
R Kaskade::square | ( | R | x | ) |
returns the square of its argument.
Definition at line 49 of file power.hh.
Referenced by Kaskade::Elastomechanics::MaterialLaws::MooneyRivlin< dimension >::d1(), Kaskade::EmbeddedErrorEstimator< VariableSetDescription, Scaling >::estimate(), and Kaskade::CoarseningDetail::GetLocalTransferProjection< Cell >::operator()().