Projection class that allows to project an ElementFormula (i.e. More...

#include <formulafrmEformula.hh>

Inheritance diagram for concepts::FormulaFromElementFormula< dim, F, G >:
concepts::Formula< F >

Public Types

typedef Realtype< F >::type G
 
enum  PNF { Exc = 0, NaN, Zero }
 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. More...
 
typedef F value_type
 

Public Member Functions

void allowPeriodicity (const Point< Real, dim > &O, const Point< Real, dim > &eL)
 This method allows for periodic destination meshes. More...
 
virtual FormulaFromElementFormula< dim, F, G > * clone () const
 
 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)
 Constructor of the projection. More...
 
virtual F operator() (const Connector &cntr, const Real p, const Real t=0.0) const
 Convenience implementation, that by default ignores its elm param. More...
 
virtual F operator() (const Connector &cntr, const Real2d &p, const Real t=0.0) const
 Convenience implementation, that by default ignores its elm param. More...
 
virtual F operator() (const Connector &cntr, const Real3d &p, const Real t=0.0) const
 Convenience implementation, that by default ignores its elm param. More...
 
virtual F operator() (const Real p, const Real t=0.0) const
 Application operator. More...
 
virtual F operator() (const Real2d &p, const Real t=0.0) const
 Application operator. More...
 
virtual F operator() (const Real3d &p, const Real t=0.0) const
 Application operator. More...
 
void setPNF (PNF pnf)
 Sets the behaviour of the class in the case when a requested point is not found, e.g. More...
 
virtual ~FormulaFromElementFormula ()
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 

Private Member Functions

void addElements_ (const uint cKey, const typename JacobianCell< dim >::cell &cell, const HashMap< const concepts::ElementWithCell< G > * > &mapToElm, HashMap< HashMap< const concepts::ElementWithCell< G > * > > &cCtE)
 
void build_directMap_ (const typename JacobianCell< dim >::cell *cell, const concepts::HashMap< const concepts::ElementWithCell< G > * > &mapToElm)
 Method to build the O(1) access mapping defined via maxDoI . More...
 
void buildBoxCandidates_ (const Point< Real, dim > P, concepts::Set< CellStripeElement< dim > > &box) const
 Build Boxes of Vertex Cells that contain the point P. More...
 
bool checkCCells_ (const std::set< typename CellType< dim >::cell * > &allCell, typename std::set< typename CellType< dim >::cell * >::const_iterator iter)
 Method checks for the preassumed of ordered Cells, starting with coarsest cells, then finer ones etc. More...
 
bool compute0_ (const typename JacobianCell< dim >::cell &cCell, const Point< Real, dim > &P, const Point< Real, dim > &cEta, const Real t, F &val) const
 Evaluates the ElementFormula in the case of coarse Cell with key cKey has a local mapping. More...
 
bool compute1_ (const typename JacobianCell< dim >::cell &cCell, const Point< Real, dim > &P, const Point< Real, dim > &cEta, const Real t, F &val) const
 Evaluates the ElementFormula in the case of coarse Cell using quasi-binary search in the dichotomic tree. More...
 
compute_ (const Point< Real, dim > &p, const Real t) const
 Internal evaluation of local coordinate in the source space and its corresponding element. More...
 
 FormulaFromElementFormula (const concepts::ElementFormulaContainer< F, G > frm, const concepts::RCP< const concepts::Set< uint > > s_attrb, const concepts::RCP< const concepts::Set< uint > > d_attrb, const Point< Real, dim > scale, const Point< Real, dim > shift, const uint N, const uint cardF_2, const Real tol, const Real interior, const concepts::Set< CellBox< dim > > cBoxes, const concepts::Set< CCell_F< dim > > rCC, const bool allowPer, const Point< Real, dim > O, const Point< Real, dim > eL, const Set< uint > mapCCells, const Set< const typename JacobianCell< dim >::cell * > coarsestCell, const Set< const typename JacobianCell< dim >::cell * > irrCCells, const concepts::HashMap< CellMap< dim, G > > mapToFinestElm, const uint maxDoI, const PNF pnf, const HashMap< HashMap< const concepts::ElementWithCell< G > * > > coarseToElm)
 
void getExtremalLevel (const typename JacobianCell< dim >::cell *cell, uint(&maxL)[dim], uint(&minL)[dim])
 
bool getOrigPoint_ (const typename JacobianCell< dim >::cell &cell, const Point< Real, dim > P, Point< Real, dim > &eta, const uint N, const Real interior, const Point< Real, dim > &x0=Real3d(0.49, 0.51, 0.49)) const
 Projected Newton scheme based on the parameter interior. More...
 
handlePNF_ () const
 
bool isCoarse_ (const typename JacobianCell< dim >::cell &cell)
 
template<uint dim1 = dim, typename std::enable_if< dim1==1u, uint >::type = 0>
void setLocalMap_ (const typename JacobianCell< dim >::cell *coarseCell, const HashMap< const ElementWithCell< G > * > &mapToElm, CellMap< dim, G > &map, const std::array< uint, dim > idxL, const std::array< uint, dim > idxR) const
 
template<uint dim2 = dim, typename std::enable_if< dim2==2u, uint >::type = 0>
void setLocalMap_ (const typename JacobianCell< dim >::cell *coarseCell, const HashMap< const ElementWithCell< G > * > &mapToElm, CellMap< dim, G > &map, const std::array< uint, dim > idxL, const std::array< uint, dim > idxR) const
 
template<uint dim3 = dim, typename std::enable_if< dim3==3u, uint >::type = 0>
void setLocalMap_ (const typename JacobianCell< dim >::cell *coarseCell, const HashMap< const ElementWithCell< G > * > &mapToElm, CellMap< dim, G > &map, const std::array< uint, dim > idxL, const std::array< uint, dim > idxR) const
 
void shift_DBox_ (const Point< Real, dim > &p, Point< Real, dim > &P) const
 shifts the physical point to the default periodic box out of the destination mesh, then applies shift and scale transformation More...
 

Private Attributes

bool allowPer_
 Flag defining whether the destination mesh is considered periodic for the formula projection. More...
 
const uint cardF_2_
 Square root of amount of interior points in a curved cell for distance evaluation. More...
 
concepts::Set< CellBox< dim > > cBoxes_
 Coarse boxes of vertex coarse cells. More...
 
Z2 cCellType_
 type of last Coarse Cell 0 = Coarse Cell with mapping => O(1) search 1 = Coarse cell has no boxing => O(log) search More...
 
Set< const typename JacobianCell< dim >::cell * > coarsestCell_
 Set of coarsest cells. More...
 
HashMap< HashMap< const concepts::ElementWithCell< G > * > > coarseToElm_
 maps from coarse cell key to all finest element children belonging to that coarse cell; this is kind of a bucket More...
 
const concepts::RCP< const concepts::Set< uint > > d_attrb_
 Destination attributes. More...
 
Point< Real, dim > eL_
 Size of the periodic box, in a case of a periodic formula projection. More...
 
const Real interior_
 Interior in $(0,0.5)$. Defines projection rule, meant to be very small. More...
 
Set< const typename JacobianCell< dim >::cell * > irrCCells_
 Set of coarsest cells that have an irregular children refinement pattern with irregular degree bigger as requested. More...
 
const JacobianCell< dim >::cell * lastCell_
 last cell to avoid dynamic cast More...
 
const JacobianCell< dim >::cell * lastCoarseCell_
 last coarse cell to avoid dynamic cast More...
 
const ElementWithCell< G > * lastElm_
 
Set< uint > mapCCells_
 Keys of all coarse cells that provide local mapping. More...
 
concepts::HashMap< CellMap< dim, G > > mapToFinestElm_
 
const uint maxDoI_
 Determines the internal map for resolving meshes with hanging nodes. More...
 
const uint N_
 Number of Newton steps. More...
 
Point< Real, dim > O_
 Origin of the periodic box, in a case of a periodic formula projection. More...
 
PNF pnf_
 
concepts::Set< CCell_F< dim > > rCC_
 Remaining coarse cells that do not provide box design, i.e. non vertex. More...
 
const concepts::RCP< const concepts::Set< uint > > s_attrb_
 Source attributes. More...
 
const Point< Real, dim > scale_
 Scaling box, from destination mesh to source mesh. More...
 
const Point< Real, dim > shift_
 Shift box, from destination mesh to source mesh. More...
 
const Real tol_
 Stop criterion for interior Newton scheme. More...
 
const concepts::ElementFormulaContainer< F, Gu_
 The original projected formula. More...
 
bool wasCurved_
 flag is last finecell was a curved one More...
 
bool wasCurved_c_
 flag if last coarse cell was curved More...
 

Static Private Attributes

static uint max_lvl
 

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.

Enumerator
Exc 
NaN 
Zero 

Definition at line 531 of file formulafrmEformula.hh.

Constructor & Destructor Documentation

◆ FormulaFromElementFormula() [1/2]

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.

◆ FormulaFromElementFormula() [2/2]

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::RCP< const concepts::Set< uint > >  s_attrb,
const concepts::RCP< const concepts::Set< uint > >  d_attrb,
const Point< Real, dim >  scale,
const Point< Real, dim >  shift,
const uint  N,
const uint  cardF_2,
const Real  tol,
const Real  interior,
const concepts::Set< CellBox< dim > >  cBoxes,
const concepts::Set< CCell_F< dim > >  rCC,
const bool  allowPer,
const Point< Real, dim >  O,
const Point< Real, dim >  eL,
const Set< uint >  mapCCells,
const Set< const typename JacobianCell< dim >::cell * >  coarsestCell,
const Set< const typename JacobianCell< dim >::cell * >  irrCCells,
const concepts::HashMap< CellMap< dim, G > >  mapToFinestElm,
const uint  maxDoI,
const PNF  pnf,
const HashMap< HashMap< const concepts::ElementWithCell< G > * > >  coarseToElm 
)
private

Member Function Documentation

◆ addElements_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
void concepts::FormulaFromElementFormula< dim, F, G >::addElements_ ( const uint  cKey,
const typename JacobianCell< dim >::cell &  cell,
const HashMap< const concepts::ElementWithCell< G > * > &  mapToElm,
HashMap< HashMap< const concepts::ElementWithCell< G > * > > &  cCtE 
)
private
Parameters
cKeykey of the coarse (grand)father cell
cellto look up for children @mapToElm map from key to all finest element @cCtE coarse cell to Element mapping

◆ 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.

◆ build_directMap_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
void concepts::FormulaFromElementFormula< dim, F, G >::build_directMap_ ( const typename JacobianCell< dim >::cell *  cell,
const concepts::HashMap< const concepts::ElementWithCell< G > * > &  mapToElm 
)
private

Method to build the O(1) access mapping defined via maxDoI .

Based on the level of irregularity on the input coarse Cell cell , the mapping mapToFinestElm_ is built in case the submesh has lower irregularity degree as maxDoI .

Parameters
cellCurrent coarse Cell
mapToElmMapping from key to all finest Elements

◆ buildBoxCandidates_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
void concepts::FormulaFromElementFormula< dim, F, G >::buildBoxCandidates_ ( const Point< Real, dim >  P,
concepts::Set< CellStripeElement< dim > > &  box 
) const
private

Build Boxes of Vertex Cells that contain the point P.

Parameters
PPhysical Point to look for.
boxBox will be builded within the routine

◆ checkCCells_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
bool concepts::FormulaFromElementFormula< dim, F, G >::checkCCells_ ( const std::set< typename CellType< dim >::cell * > &  allCell,
typename std::set< typename CellType< dim >::cell * >::const_iterator  iter 
)
private

Method checks for the preassumed of ordered Cells, starting with coarsest cells, then finer ones etc.

Retruns true if the order is correct.

◆ 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

◆ compute0_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
bool concepts::FormulaFromElementFormula< dim, F, G >::compute0_ ( const typename JacobianCell< dim >::cell &  cCell,
const Point< Real, dim > &  P,
const Point< Real, dim > &  cEta,
const Real  t,
F &  val 
) const
private

Evaluates the ElementFormula in the case of coarse Cell with key cKey has a local mapping.

Parameters
cEtais the computed local coordinate in the coarse cell

◆ compute1_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
bool concepts::FormulaFromElementFormula< dim, F, G >::compute1_ ( const typename JacobianCell< dim >::cell &  cCell,
const Point< Real, dim > &  P,
const Point< Real, dim > &  cEta,
const Real  t,
F &  val 
) const
private

Evaluates the ElementFormula in the case of coarse Cell using quasi-binary search in the dichotomic tree.

Parameters
cEtais the computed local coordinate in the coarse cell

◆ compute_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
F concepts::FormulaFromElementFormula< dim, F, G >::compute_ ( const Point< Real, dim > &  p,
const Real  t 
) const
private

Internal evaluation of local coordinate in the source space and its corresponding element.

General strategy:

Search on lastElement first to be the containing element. If that fails, search on the lastCoarsecell, and update last element if succeed.

If that also fails, then we will look in all coarse cells for the correct one containing the point and update the last coarse Cell, and last Elements.

◆ getExtremalLevel()

template<uint dim, class F , typename G = typename Realtype<F>::type>
void concepts::FormulaFromElementFormula< dim, F, G >::getExtremalLevel ( const typename JacobianCell< dim >::cell *  cell,
uint(&)  maxL[dim],
uint(&)  minL[dim] 
)
private

◆ getOrigPoint_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
bool concepts::FormulaFromElementFormula< dim, F, G >::getOrigPoint_ ( const typename JacobianCell< dim >::cell &  cell,
const Point< Real, dim >  P,
Point< Real, dim > &  eta,
const uint  N,
const Real  interior,
const Point< Real, dim > &  x0 = Real3d(0.49, 0.51, 0.49) 
) const
private

Projected Newton scheme based on the parameter interior.

Parameters
elmpossible source space element to look in
cellcorresponding elm's cell that supports jacobian functionality
Pphysical Point
etalocal point that if exists in elm is mapped by underlying mapping to P
interiorheuritic parameter to avoid vertex locks
x0start value
Note
We do not take _i = 0.5$ to avoir vertex oscillation

◆ handlePNF_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
F concepts::FormulaFromElementFormula< dim, F, G >::handlePNF_ ( ) const
private

◆ 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

Definition at line 623 of file formulafrmEformula.hh.

◆ isCoarse_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
bool concepts::FormulaFromElementFormula< dim, F, G >::isCoarse_ ( const typename JacobianCell< dim >::cell &  cell)
private

◆ operator()() [1/6]

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/6]

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

Convenience implementation, that by default ignores its elm param.

Reimplemented from concepts::Formula< F >.

◆ operator()() [3/6]

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

Convenience implementation, that by default ignores its elm param.

Reimplemented from concepts::Formula< F >.

◆ operator()() [4/6]

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()() [5/6]

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

Application operator.

Evaluates the formula.

Parameters
pPoint in space
tPoint in time

Implements concepts::Formula< F >.

◆ operator()() [6/6]

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

Application operator.

Evaluates the formula.

Parameters
pPoint in space
tPoint in time

Implements concepts::Formula< F >.

◆ setLocalMap_() [1/3]

template<uint dim, class F , typename G = typename Realtype<F>::type>
template<uint dim1 = dim, typename std::enable_if< dim1==1u, uint >::type = 0>
void concepts::FormulaFromElementFormula< dim, F, G >::setLocalMap_ ( const typename JacobianCell< dim >::cell *  coarseCell,
const HashMap< const ElementWithCell< G > * > &  mapToElm,
CellMap< dim, G > &  map,
const std::array< uint, dim >  idxL,
const std::array< uint, dim >  idxR 
) const
private

◆ setLocalMap_() [2/3]

template<uint dim, class F , typename G = typename Realtype<F>::type>
template<uint dim2 = dim, typename std::enable_if< dim2==2u, uint >::type = 0>
void concepts::FormulaFromElementFormula< dim, F, G >::setLocalMap_ ( const typename JacobianCell< dim >::cell *  coarseCell,
const HashMap< const ElementWithCell< G > * > &  mapToElm,
CellMap< dim, G > &  map,
const std::array< uint, dim >  idxL,
const std::array< uint, dim >  idxR 
) const
private

◆ setLocalMap_() [3/3]

template<uint dim, class F , typename G = typename Realtype<F>::type>
template<uint dim3 = dim, typename std::enable_if< dim3==3u, uint >::type = 0>
void concepts::FormulaFromElementFormula< dim, F, G >::setLocalMap_ ( const typename JacobianCell< dim >::cell *  coarseCell,
const HashMap< const ElementWithCell< G > * > &  mapToElm,
CellMap< dim, G > &  map,
const std::array< uint, dim >  idxL,
const std::array< uint, dim >  idxR 
) const
private

◆ 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.

◆ shift_DBox_()

template<uint dim, class F , typename G = typename Realtype<F>::type>
void concepts::FormulaFromElementFormula< dim, F, G >::shift_DBox_ ( const Point< Real, dim > &  p,
Point< Real, dim > &  P 
) const
private

shifts the physical point to the default periodic box out of the destination mesh, then applies shift and scale transformation

Parameters
pPhysical point of the destination mesh
PPhysical point of the source mesh

Member Data Documentation

◆ allowPer_

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

Flag defining whether the destination mesh is considered periodic for the formula projection.

Definition at line 678 of file formulafrmEformula.hh.

◆ cardF_2_

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

Square root of amount of interior points in a curved cell for distance evaluation.

Definition at line 648 of file formulafrmEformula.hh.

◆ cBoxes_

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

Coarse boxes of vertex coarse cells.

Definition at line 745 of file formulafrmEformula.hh.

◆ cCellType_

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

type of last Coarse Cell 0 = Coarse Cell with mapping => O(1) search 1 = Coarse cell has no boxing => O(log) search

Definition at line 674 of file formulafrmEformula.hh.

◆ coarsestCell_

template<uint dim, class F , typename G = typename Realtype<F>::type>
Set<const typename JacobianCell<dim>::cell* > concepts::FormulaFromElementFormula< dim, F, G >::coarsestCell_
private

Set of coarsest cells.

Definition at line 696 of file formulafrmEformula.hh.

◆ coarseToElm_

template<uint dim, class F , typename G = typename Realtype<F>::type>
HashMap<HashMap<const concepts::ElementWithCell<G>* > > concepts::FormulaFromElementFormula< dim, F, G >::coarseToElm_
mutableprivate

maps from coarse cell key to all finest element children belonging to that coarse cell; this is kind of a bucket

Definition at line 741 of file formulafrmEformula.hh.

◆ d_attrb_

template<uint dim, class F , typename G = typename Realtype<F>::type>
const concepts::RCP<const concepts::Set<uint> > concepts::FormulaFromElementFormula< dim, F, G >::d_attrb_
private

Destination attributes.

Definition at line 636 of file formulafrmEformula.hh.

◆ eL_

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

Size of the periodic box, in a case of a periodic formula projection.

Definition at line 683 of file formulafrmEformula.hh.

◆ interior_

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

Interior in $(0,0.5)$. Defines projection rule, meant to be very small.

Definition at line 650 of file formulafrmEformula.hh.

◆ irrCCells_

template<uint dim, class F , typename G = typename Realtype<F>::type>
Set<const typename JacobianCell<dim>::cell* > concepts::FormulaFromElementFormula< dim, F, G >::irrCCells_
private

Set of coarsest cells that have an irregular children refinement pattern with irregular degree bigger as requested.

Definition at line 700 of file formulafrmEformula.hh.

◆ lastCell_

template<uint dim, class F , typename G = typename Realtype<F>::type>
const JacobianCell<dim>::cell* concepts::FormulaFromElementFormula< dim, F, G >::lastCell_
mutableprivate

last cell to avoid dynamic cast

Definition at line 657 of file formulafrmEformula.hh.

◆ lastCoarseCell_

template<uint dim, class F , typename G = typename Realtype<F>::type>
const JacobianCell<dim>::cell* concepts::FormulaFromElementFormula< dim, F, G >::lastCoarseCell_
mutableprivate

last coarse cell to avoid dynamic cast

Definition at line 660 of file formulafrmEformula.hh.

◆ lastElm_

template<uint dim, class F , typename G = typename Realtype<F>::type>
const ElementWithCell<G>* concepts::FormulaFromElementFormula< dim, F, G >::lastElm_
mutableprivate

Definition at line 654 of file formulafrmEformula.hh.

◆ mapCCells_

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

Keys of all coarse cells that provide local mapping.

Definition at line 750 of file formulafrmEformula.hh.

◆ mapToFinestElm_

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

Definition at line 735 of file formulafrmEformula.hh.

◆ max_lvl

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

Definition at line 522 of file formulafrmEformula.hh.

◆ maxDoI_

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

Determines the internal map for resolving meshes with hanging nodes.

In particular this is the maximal degree of Irregularity taken into account with the O(1) map.

Definition at line 688 of file formulafrmEformula.hh.

◆ N_

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

Number of Newton steps.

Definition at line 646 of file formulafrmEformula.hh.

◆ O_

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

Origin of the periodic box, in a case of a periodic formula projection.

Definition at line 681 of file formulafrmEformula.hh.

◆ pnf_

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

Definition at line 691 of file formulafrmEformula.hh.

◆ rCC_

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

Remaining coarse cells that do not provide box design, i.e. non vertex.

Definition at line 748 of file formulafrmEformula.hh.

◆ s_attrb_

template<uint dim, class F , typename G = typename Realtype<F>::type>
const concepts::RCP<const concepts::Set<uint> > concepts::FormulaFromElementFormula< dim, F, G >::s_attrb_
private

Source attributes.

Definition at line 634 of file formulafrmEformula.hh.

◆ scale_

template<uint dim, class F , typename G = typename Realtype<F>::type>
const Point<Real, dim> concepts::FormulaFromElementFormula< dim, F, G >::scale_
private

Scaling box, from destination mesh to source mesh.

Definition at line 639 of file formulafrmEformula.hh.

◆ shift_

template<uint dim, class F , typename G = typename Realtype<F>::type>
const Point<Real, dim> concepts::FormulaFromElementFormula< dim, F, G >::shift_
private

Shift box, from destination mesh to source mesh.

Definition at line 641 of file formulafrmEformula.hh.

◆ tol_

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

Stop criterion for interior Newton scheme.

Definition at line 644 of file formulafrmEformula.hh.

◆ u_

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

The original projected formula.

Definition at line 631 of file formulafrmEformula.hh.

◆ wasCurved_

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

flag is last finecell was a curved one

Definition at line 663 of file formulafrmEformula.hh.

◆ wasCurved_c_

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

flag if last coarse cell was curved

Definition at line 666 of file formulafrmEformula.hh.


The documentation for this class was generated from the following file:
Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich