Class documentation of Concepts

Loading...
Searching...
No Matches
concepts::FormulaFromElementFormula< dim, F, G > Class Template Referenceabstract

#include <formulafrmEformula.hh>

Inheritance diagram for concepts::FormulaFromElementFormula< dim, F, G >:
concepts::Formula< F > concepts::PiecewiseFormulaBase< F > concepts::ElementFormula< F, concepts::Realtype< F >::type > concepts::Cloneable concepts::OutputOperator

Public Types

enum  PNF { Exc = 0 , NaN , Zero }
 
typedefvalue_type
 
typedef Realtype< F >::type G
 

Public Member Functions

 FormulaFromElementFormula (const concepts::ElementFormulaContainer< F, G > frm, const concepts::SpaceOnCoarseCells< dim, G > &sSpc, const Set< uint > *s_attrb=0, const Set< uint > *d_attrb=0, const Point< Real, dim > scale=Point< Real, dim >(1), const Point< Real, dim > shift=Point< Real, dim >(0), const uint N=10, const uint cardF_2=2, const Real tol=EPS, const Real interior=0.00000001, const uint maxDoI=3)
 
virtualoperator() (const Real p, const Real t=0.0) const
 
virtualoperator() (const Real2d &p, const Real t=0.0) const
 
virtualoperator() (const Real3d &p, const Real t=0.0) const
 
virtualoperator() (const Connector &cntr, const Real p, const Real t=0.0) const
 
virtualoperator() (const Connector &cntr, const Real2d &p, const Real t=0.0) const
 
virtualoperator() (const Connector &cntr, const Real3d &p, const Real t=0.0) const
 
void allowPeriodicity (const Point< Real, dim > &O, const Point< Real, dim > &eL)
 
void setPNF (PNF pnf)
 
virtual FormulaFromElementFormula< dim, F, G > * clone () const
 Virtual copy constructor.
 
virtualoperator() (const ElementWithCell< G > &elm, const Real p, const Real t=0.0) const
 
virtualoperator() (const ElementWithCell< G > &elm, const Real2d &p, const Real t=0.0) const
 
virtualoperator() (const ElementWithCell< G > &elm, const Real3d &p, const Real t=0.0) const
 
virtualoperator() (const ElementWithCell< concepts::Realtype< F >::type > &elm, const Real p, const Real t=0.0) const=0
 
virtualoperator() (const ElementWithCell< concepts::Realtype< F >::type > &elm, const Real2d &p, const Real t=0.0) const=0
 
virtualoperator() (const ElementWithCell< concepts::Realtype< F >::type > &elm, const Real3d &p, const Real t=0.0) const=0
 
virtual const F & dflt_value () const
 Gives default value.
 
virtual F & dflt_value ()
 Gives default value.
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream.
 

Detailed Description

template<uint dim, class F, typename G = typename Realtype<F>::type>
class concepts::FormulaFromElementFormula< dim, F, G >

Projection class that allows to project an ElementFormula (i.e. a ElementFormulaContainer) from a Space sSpc with a submesh sMsh defined by attributes onto a different Mesh that can be modified by given attributes.

It is possible that the destination Mesh dMesh is shifted or scaled. That is let the source mesh form a domain Omega_1. Let Omega_2 be a subset of Omega_1, then the destination mesh is allowed to be of the form scale * Omega_2 + shift. The scaling is allowed in each coordinate.

Furthermore the class supports periodic extensions of Elementformulas. In this case the destination mesh is bigger than the source mesh. For this we assume that the destination mesh is a periodic layer of a default box possible in various directions. The default box can be of the form scale * B + shift, where B is a subset of Omega_1. This is done with the help of the allowPeriodicity routine.

Functions build on vectorial spaces can be used as well, but the space components have to be built on the same mesh (including refinement). Then as sSpc one simply can add one particular space component and its mesh will be used. Note: If the mesh is not the same and one inserts one particular (scalar) source (space)mesh, then Newton subroutine will fail to find the cells.

This class works for tensored cells(Edge1d, Quad2d, Hexaedron3d) only in dimension 1, 2 and 3 that can be convex and non-convex as well. For an extension to more cells, the stripe algorithm and newtons start values may be adapted.

This class is not meant to provide a fast solution in terms of integration speed etc. This is clear as this class needs to compute for each point requested out of the destination space the corresponding element in the original space followed by evaluating the elementformula itsself.

Important notes:

Assume the source space msh and the destination space msh (up to scaling shiftiny and attribute modifiying) are the same, then the class does not simplify to the identity on general ElementFormulas.

However if the Elementformula is continuous over the mesh cell interfaces (points, edges, faces) then the latter holds true.

If the Elementformula is discontinuous across an interface, make sure that physical Points to evaluate do not intersect such interfaces. Then the (discrete) identity mapping remains.

Example of non identity: 1) The ElementFormla is an ElementFormulaVector<dim>(sSpc, Vector<F> sol, Grad<F>()), then it is discontinuous across interfaces. If one applies static overall trapeze integration rule, the physical Points can touch the interfaces. The interfaces has at least two adjoining elements K and L. Since the ElementFormulaVector is evaluated per element it gives different values on K and L for the same physical (and therefore adapted corresponding local ) point. Since the internal element search routine is deterministic, it is possible that the ElementformulaVector is evaluated on the Element K only for the same physical point coming from two different cells. The difference Projection - Elementformula cannot be zero in this case.

2) The projection for obvious reasons depends on the source and destination space and more important on the underlying quadrature points. Think of periodic projection of an elementformla defined on $(0,1)^2$ to $(0,N) \times (0,1)$, $N$ is natural. Both meshes consist of one cell only. Due to the underlying quadrature rule, the point in the destination space get stretched. The evaluation happens on different points in $(0,1)^2$ now, since physical strechted quadrature points after shifted to the default box do not coincide with original quadrature points. In this case the pointwise identity cannot hold true in general. In order to still get $||proj||_{(0,1)x[0,1|} = N *
||frm||_{(0,1)}^2$, what may is called generalized identity, one needs finer mesh and finer quadrature.

It is also possible only to project just parts of the Elementformula. This is described in the rule of attributes.

The rule of attributes:

The parameter s_attrib, shrinks the Set of all cells within the sourceSpace to look up. This applicates for example if we only want to project a part of a original solution, e.g. close to a corner etc. By default all cells are used.

The parameter d_attrib, shrinks the Set of all cells on which the projection is evaluated. If the destination cell is not requested we extend the projection by zero on these cells. By default all cells are used.

Note that (without applying scaling and shifting) the cells avaiable after shrinking by d_attrib form a domain that must be a subset of the domain formed by cells after shrinking with s_attrib. If this is not the case the underlying newton scheme cannot find a original element, and therefore an error occurs. On remaining cells not requested by d_attrib the formula is extended by zero.

Short description of underlying algorithm:

For a given physical point in the destination mesh we look for the original local coordinate and the corresponding element in the original space. We search within vertex cells by an stripe intersection that we call box search, that is a p=\infty norm idea. For non vertex cells search for the correct cell in the order of minimal distance of P and requested Set of Point evaluations of the underlying map, controlled by cardF_2. For all cells we apply newtons routine that if succeeding returns original local coordinate. Vertexcells are found in O(log(N)) where N is the number of vertex cells in the source space. Non VertexCells are found in O(N), where N is the number of non vertex Cells in the source space.

Note that the newton routine which is applied withing [0,1]^dim can have different behaviour, that is

  • convergence to the loca coordinate
  • leave [0,1]^dim, then we project back to the boundary
  • visits vertex points (e.g by projection), then we apply interior rule
  • no convergence at all, either increase Newton steps N, the element does not contain the point, or newton scheme oscilation between same points

Note that oscilation is not handled at the moment.

Projection rule defined by interior:

Attention:

When the newton scheme ends in a Vertex of [0,1]^2, we project it to the interior of [0,1]^2, in the case of the vertex is not a solution already. Therefore the parameter interior is only allowed to be within (0,0.5).

For example the point (1, 0)^T will be projected to (1-interior, interior)^T

This applicates since in the case of vertex, the jacobian is not invertible anymore. The deriviation dependent newton scheme fails.

One workaround of the latter one could be to use derivation free root finder in the vertex endup case. Since this is supposed to be very slow, at the moment the above heuristic approach is made.

In actual practical testing this behaves nice for vertex cells, non and convex non vertex cells in 2d.

Furthermore interior should not be choosen too small, i.e. so that 1-interior = 1 in machine precision. But interior should be chosen smaller than the smallest non zero quadrature points' coordinate.

With this idea we heuristically ensure that the newton step goes away from vertex in the next step. This corresponds with an educated guess of a new startpoint.

Precondition
The input SpaceOnCoarseCell provides all cells strictly beginning with the coarse cells then finer cells, else a exception is thrown.
Author
Robert Gruhlke, 2015, 2016
Adrien Semin, 2016

Definition at line 518 of file formulafrmEformula.hh.

Member Typedef Documentation

◆ G

template<typename F >
typedef Realtype<F>::type concepts::Formula< F >::G
inherited

Definition at line 37 of file formula.hh.

◆ value_type

template<typename F >
typedef F concepts::Formula< F >::value_type
inherited

Definition at line 36 of file formula.hh.

Member Enumeration Documentation

◆ PNF

template<uint dim, class F , typename G = typename Realtype<F>::type>
enum concepts::FormulaFromElementFormula::PNF

Exc = 0, by default exception is thrown NaN = 1, if point is not found return a NaN Zero = 2, if point is not found return zero.

Definition at line 531 of file formulafrmEformula.hh.

Constructor & Destructor Documentation

◆ FormulaFromElementFormula()

template<uint dim, class F , typename G = typename Realtype<F>::type>
concepts::FormulaFromElementFormula< dim, F, G >::FormulaFromElementFormula ( const concepts::ElementFormulaContainer< F, G frm,
const concepts::SpaceOnCoarseCells< dim, G > &  sSpc,
const Set< uint > *  s_attrb = 0,
const Set< uint > *  d_attrb = 0,
const Point< Real, dim >  scale = PointReal, dim >(1),
const Point< Real, dim >  shift = PointReal, dim >(0),
const uint  N = 10,
const uint  cardF_2 = 2,
const Real  tol = EPS,
const Real  interior = 0.00000001,
const uint  maxDoI = 3 
)

Constructor of the projection

Parameters
uFormula that needs to be projected
sSpcSource Mesh onto which u is defined
s_attrbSource attributes
d_attrbDestination attributes
scalePossible scale parameter
shiftPossible shift parameter
tolStop criterion for interior Newton scheme
NNumber of Newton steps
cardF_2Square root of amount of interior points in a curved cell for distance evaluation.
interiorinterior in (0,0.5). defines projection rule, meant to be very small
maxDoINumber of maximal degree of irregularity of mesh, such that O(1) map is build from coarse cell to its submesh, if submesh exceeds this number no map is build and default search is performed in O(log(N))
Warning
Do not overuse the maxDoI in terms of size to avoid big memory allocation, e.g. caused by geometric meshes.

◆ ~FormulaFromElementFormula()

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual concepts::FormulaFromElementFormula< dim, F, G >::~FormulaFromElementFormula ( )
inlinevirtual

Definition at line 567 of file formulafrmEformula.hh.

Member Function Documentation

◆ allowPeriodicity()

template<uint dim, class F , typename G = typename Realtype<F>::type>
void concepts::FormulaFromElementFormula< dim, F, G >::allowPeriodicity ( const Point< Real, dim > &  O,
const Point< Real, dim > &  eL 
)

This method allows for periodic destination meshes. In this case the destination domain is bigger then the original domain.

We assume the default mesh that is shift copied to periodic mesh be described within a cartesian box that is line/rectangle/hexaedron. The latter one is described by origin position and by its edge lengths. The orientation of the box is to be assumed by default cartesian orientation.

Example: in two dimensions, a rectangle with origin $O$ and with rectangle size $eL$ will produce the domain $(O(1),O(1)+eL(1))
\times (O(2),O(2)+eL(2))$.

Parameters
OOrigin of the default box
eLSize of the default box
Note
The default box lies in the destination domain.

◆ clone()

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual FormulaFromElementFormula< dim, F, G > * concepts::FormulaFromElementFormula< dim, F, G >::clone ( ) const
virtual

Virtual copy constructor.

Implements concepts::Formula< F >.

◆ dflt_value() [1/2]

template<typename F >
virtual F & concepts::PiecewiseFormulaBase< F >::dflt_value ( )
inlinevirtualinherited

Gives default value.

Definition at line 83 of file piecewiseFormula.hh.

◆ dflt_value() [2/2]

template<typename F >
virtual const F & concepts::PiecewiseFormulaBase< F >::dflt_value ( ) const
inlinevirtualinherited

Gives default value.

Definition at line 81 of file piecewiseFormula.hh.

◆ info()

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual std::ostream & concepts::FormulaFromElementFormula< dim, F, G >::info ( std::ostream &  os) const
inlineprotectedvirtual

Returns information in an output stream.

Reimplemented from concepts::PiecewiseFormulaBase< F >.

Definition at line 623 of file formulafrmEformula.hh.

◆ operator()() [1/12]

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual F concepts::FormulaFromElementFormula< dim, F, G >::operator() ( const Connector cntr,
const Real  p,
const Real  t = 0.0 
) const
virtual

Convenience implementation, that by default ignores its elm param.

Reimplemented from concepts::Formula< F >.

◆ operator()() [2/12]

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual F concepts::FormulaFromElementFormula< dim, F, G >::operator() ( const Connector cntr,
const Real2d p,
const Real  t = 0.0 
) const
virtual

Reimplemented from concepts::Formula< F >.

◆ operator()() [3/12]

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual F concepts::FormulaFromElementFormula< dim, F, G >::operator() ( const Connector cntr,
const Real3d p,
const Real  t = 0.0 
) const
virtual

Reimplemented from concepts::Formula< F >.

◆ operator()() [4/12]

virtual F concepts::ElementFormula< F, concepts::Realtype< F >::type >::operator() ( const ElementWithCell< concepts::Realtype< F >::type > &  elm,
const Real  p,
const Real  t = 0.0 
) const
pure virtualinherited

Evaluates the formula.

Parameters
elmElement
pPoint in space in local element coordinates
tPoint in time

◆ operator()() [5/12]

virtual F concepts::ElementFormula< F, concepts::Realtype< F >::type >::operator() ( const ElementWithCell< concepts::Realtype< F >::type > &  elm,
const Real2d p,
const Real  t = 0.0 
) const
pure virtualinherited

Evaluates the formula.

Parameters
elmElement
pPoint in space in local element coordinates
tPoint in time

◆ operator()() [6/12]

virtual F concepts::ElementFormula< F, concepts::Realtype< F >::type >::operator() ( const ElementWithCell< concepts::Realtype< F >::type > &  elm,
const Real3d p,
const Real  t = 0.0 
) const
pure virtualinherited

Evaluates the formula.

Parameters
elmElement
pPoint in space in local element coordinates
tPoint in time

◆ operator()() [7/12]

template<typename F >
virtual F concepts::PiecewiseFormulaBase< F >::operator() ( const ElementWithCell< G > &  elm,
const Real  p,
const Real  t = 0.0 
) const
inlinevirtualinherited

Definition at line 52 of file piecewiseFormula.hh.

◆ operator()() [8/12]

template<typename F >
virtual F concepts::PiecewiseFormulaBase< F >::operator() ( const ElementWithCell< G > &  elm,
const Real2d p,
const Real  t = 0.0 
) const
inlinevirtualinherited

Definition at line 57 of file piecewiseFormula.hh.

◆ operator()() [9/12]

template<typename F >
virtual F concepts::PiecewiseFormulaBase< F >::operator() ( const ElementWithCell< G > &  elm,
const Real3d p,
const Real  t = 0.0 
) const
inlinevirtualinherited

Definition at line 63 of file piecewiseFormula.hh.

◆ operator()() [10/12]

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual F concepts::FormulaFromElementFormula< dim, F, G >::operator() ( const Real  p,
const Real  t = 0.0 
) const
virtual

Application operator. Evaluates the formula.

Parameters
pPoint in space
tPoint in time

Implements concepts::Formula< F >.

◆ operator()() [11/12]

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual F concepts::FormulaFromElementFormula< dim, F, G >::operator() ( const Real2d p,
const Real  t = 0.0 
) const
virtual

◆ operator()() [12/12]

template<uint dim, class F , typename G = typename Realtype<F>::type>
virtual F concepts::FormulaFromElementFormula< dim, F, G >::operator() ( const Real3d p,
const Real  t = 0.0 
) const
virtual

◆ setPNF()

template<uint dim, class F , typename G = typename Realtype<F>::type>
void concepts::FormulaFromElementFormula< dim, F, G >::setPNF ( PNF  pnf)
inline

Sets the behaviour of the class in the case when a requested point is not found, e.g. cause the cell is not there or the Newton Scheme fails. For details see documentation on the enum struct.

Parameters
pnfPoint not found handling. Choose between three behaviours Exc, NaN, Zero.

Definition at line 617 of file formulafrmEformula.hh.


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