Element on an edge representing the normal derivatives of neighbouring elements, especially their mean or jump. More...

#include <function.hh>

Classes

struct  mapPart
 

Public Types

typedef concepts::ElementAndFacette< hp3D::Element< F > > UnderlyingElement
 

Public Member Functions

void addElement (const hp3D::Hexahedron &hex, uint k, Real weight=1.0)
 Adds the contribution to the Neumann trace from the element on one of the (at most two) sides. More...
 
void buildT ()
 After adding the Elements this routine builds the TMatrix and polynomial informations on the element. More...
 
void computeShpfctVals (const concepts::Real2d &xP, concepts::Array< Real > &val) const
 Evaluation method to get basisfunctions evaluated on requested points this should be used for functions, i.e. More...
 
 NeumannTraceElement3d (concepts::QuadNd &cell, const ushort *p)
 Constructor. More...
 
virtual const concepts::TMatrix< F > & T () const
 
const concepts::Sequence< UnderlyingElement > & uelm () const
 Returns the Underlying Element(s) More...
 
virtual ~NeumannTraceElement3d ()
 

Protected Member Functions

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

Private Member Functions

void compute_ (const concepts::Real2d &q, concepts::Array< Real > &val) const
 Routine calculating the normal derivative of the basis functions (of underlying quads), with normalvector pointing outwards the respective underlying Quad on a given Set of Points. More...
 

Private Attributes

mapPart coarseCell_
 
bool isComputed_
 
std::unique_ptr< concepts::ShapeFunction1D< Real > > shpfct_
 
concepts::TMatrix< Real > T_
 
concepts::Sequence< UnderlyingElementuelm_
 
concepts::Array< Real > values_
 The actual storage of the values of the shape functions along the edge on integration points. More...
 
concepts::Sequence< Real > weights_
 Weight for each side. This defines the TraceType on the Element. More...
 

Detailed Description

template<class F = Real>
class hp3D::NeumannTraceElement3d< F >

Element on an edge representing the normal derivatives of neighbouring elements, especially their mean or jump.

Also a PLUS and MINUS setting is involded, that is a FIRST trace kind, with respect to the K+ and K- introduced elements via the setted Normalvector.

The number of shape functions are the sum of the number of shape functions of the two neighbouring 2D elements.

Author
Robert Gruhlke 2013

Definition at line 20 of file function.hh.

Member Typedef Documentation

◆ UnderlyingElement

Definition at line 44 of file neumannTraceElement3d.hh.

Constructor & Destructor Documentation

◆ NeumannTraceElement3d()

template<class F = Real>
hp3D::NeumannTraceElement3d< F >::NeumannTraceElement3d ( concepts::QuadNd cell,
const ushort *  p 
)
inline

Constructor.

For irregular meshes the underlying elements of a neumanntraceelement is a fine Quad and a neighboured coarse Quad. The coarseCell describes the edgemapping of the coarseedge of which the cell is a subset (part). The coarseCellpart also marks if a neumanntraceElement is build upon regular or irregular mesh

Parameters
cellgeometric cell of this element
phighest Polynomial degree of underlying Elements
coarseCellfor irregular meshes where cell's underlying

Definition at line 120 of file neumannTraceElement3d.hh.

◆ ~NeumannTraceElement3d()

template<class F = Real>
virtual hp3D::NeumannTraceElement3d< F >::~NeumannTraceElement3d ( )
virtual

Member Function Documentation

◆ addElement()

template<class F = Real>
void hp3D::NeumannTraceElement3d< F >::addElement ( const hp3D::Hexahedron hex,
uint  k,
Real  weight = 1.0 
)

Adds the contribution to the Neumann trace from the element on one of the (at most two) sides.

Basis functions between quads do not have to be continous. In this case the highest polynomial degree along edge is set, i.e. for integration rule

Parameters
hexHexahedron with k-th edge is NeumannTraceElement's topological cell.

◆ buildT()

template<class F = Real>
void hp3D::NeumannTraceElement3d< F >::buildT ( )

After adding the Elements this routine builds the TMatrix and polynomial informations on the element.

◆ compute_()

template<class F = Real>
void hp3D::NeumannTraceElement3d< F >::compute_ ( const concepts::Real2d q,
concepts::Array< Real > &  val 
) const
private

Routine calculating the normal derivative of the basis functions (of underlying quads), with normalvector pointing outwards the respective underlying Quad on a given Set of Points.

Parameters
qevaluation domain points
valof evaluations of shapefunctions on q values

◆ computeShpfctVals()

template<class F = Real>
void hp3D::NeumannTraceElement3d< F >::computeShpfctVals ( const concepts::Real2d xP,
concepts::Array< Real > &  val 
) const

Evaluation method to get basisfunctions evaluated on requested points this should be used for functions, i.e.

hp2D::Value

Parameters
xPrequested Points to evaluate on
valevaluated shpfct, method resizes array.

◆ info()

template<class F = Real>
virtual std::ostream& hp3D::NeumannTraceElement3d< F >::info ( std::ostream &  os) const
protectedvirtual

◆ T()

template<class F = Real>
virtual const concepts::TMatrix<F>& hp3D::NeumannTraceElement3d< F >::T ( ) const
inlinevirtual

Definition at line 126 of file neumannTraceElement3d.hh.

◆ uelm()

template<class F = Real>
const concepts::Sequence<UnderlyingElement>& hp3D::NeumannTraceElement3d< F >::uelm ( ) const
inline

Returns the Underlying Element(s)

Definition at line 160 of file neumannTraceElement3d.hh.

Member Data Documentation

◆ coarseCell_

template<class F = Real>
mapPart hp3D::NeumannTraceElement3d< F >::coarseCell_
private

Definition at line 205 of file neumannTraceElement3d.hh.

◆ isComputed_

template<class F = Real>
bool hp3D::NeumannTraceElement3d< F >::isComputed_
private

Definition at line 188 of file neumannTraceElement3d.hh.

◆ shpfct_

template<class F = Real>
std::unique_ptr<concepts::ShapeFunction1D<Real> > hp3D::NeumannTraceElement3d< F >::shpfct_
private

Definition at line 185 of file neumannTraceElement3d.hh.

◆ T_

template<class F = Real>
concepts::TMatrix<Real> hp3D::NeumannTraceElement3d< F >::T_
private

Definition at line 182 of file neumannTraceElement3d.hh.

◆ uelm_

template<class F = Real>
concepts::Sequence<UnderlyingElement> hp3D::NeumannTraceElement3d< F >::uelm_
private

Definition at line 190 of file neumannTraceElement3d.hh.

◆ values_

template<class F = Real>
concepts::Array<Real> hp3D::NeumannTraceElement3d< F >::values_
mutableprivate

The actual storage of the values of the shape functions along the edge on integration points.

The values corresponding to the local shapefunctions of its underlying quad(s) are computed with the normalvector pointing outwards the quad respectivly.

This have to be noticed for defining the weights for trace types. mutable to work with const compute_ as needed for const

Definition at line 199 of file neumannTraceElement3d.hh.

◆ weights_

template<class F = Real>
concepts::Sequence<Real> hp3D::NeumannTraceElement3d< F >::weights_
private

Weight for each side. This defines the TraceType on the Element.

Definition at line 202 of file neumannTraceElement3d.hh.


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