hp3D Namespace Reference

3D hp-FEM for H1-conforming elements. More...

Classes

class  Advection
 A function class to calculate element matrices of the bilinear form. More...
 
class  APrioriRefinement
 Carries out a-priori given refinements. More...
 
class  ArrayElementFormula
 Array of formula values on quadrature points. More...
 
class  ArrayHexaWeights
 Class to represent the quadrature weights on all quadrature points. More...
 
class  BilinearFormTwoPartDeriv
 A class to calculate the element matrices for partial derivatives. More...
 
class  BuildDofsBase
 Responsible to build the degrees of freedom in a space. More...
 
class  BuildEdgeDofs
 Responsible to build the edge degrees of freedom in a space. More...
 
class  BuildFaceDofs
 Responsible to build the face degrees of freedom in a space. More...
 
class  BuildFaceDofsHypTrunk
 Responsible to build the inner degrees of freedom in a (trunk) space. More...
 
class  BuildFaceDofsLinTrunk
 Responsible to build the face degrees of freedom in a (trunk) space. More...
 
class  BuildInnerDofs
 Responsible to build the inner degrees of freedom in a space. More...
 
class  BuildInnerDofsHypTrunk
 Responsible to build the inner degrees of freedom in a (trunk) space. More...
 
class  BuildInnerDofsLinTrunk
 Responsible to build the inner degrees of freedom in a (trunk) space. More...
 
class  BuildVertexDofs
 Responsible to build the vertex degrees of freedom in a space. More...
 
class  DaugeWeight
 Weight implemented by Monique Dauge in Melina. More...
 
class  DistancePost
 A function class for weighted regularization. More...
 
class  DivDiv
 A function class to calculate element matrices for the Div u*Div v Bilinearform. More...
 
class  Element
 Abstract class for a 3D FEM element. More...
 
class  ElementFormulaVectorOnTrace
 Class for evaluation of solutions on the TraceSpace with an ElementFunction that is only given on the whole space. More...
 
class  FormulaFromWeight
 Makes it possible to plot a given Weight. More...
 
class  Grad
 The gradient of an approximated function in a FE space. More...
 
class  Hexahedron
 A 3D FEM element: a hexahedron. More...
 
class  HexahedronFaceBase
 Class to construct an quadrilateral element out of an hexahedron. More...
 
class  HexahedronFaceFirst
 
class  HexahedronGraphics
 Handles graphics for hexahedral 3D hp FEM elements. More...
 
class  HexFunctions
 Auxiliary functions for hexahedra. More...
 
class  Hook
 A function class to calculate element matrices for the Bilinearform of linear Elasticity. More...
 
class  Identity
 A function class to calculate element matrices for the mass matrix. More...
 
class  Laplace
 A function class to calculate element matrices for the Laplacian. More...
 
class  Laplacian
 The Laplacian of an approximated function in a FE space. More...
 
class  LinearElasticity
 A class to calculate the element matrices for problems of linear elasticity in 3D. More...
 
class  LinearFormHelper_1
 Helper class for linearforms l(v), where v is a one form. More...
 
class  MaxwellSharedData
 Shared data for RotRot and DivDiv. More...
 
class  NeumannTraceElement3d
 Element on an edge representing the normal derivatives of neighbouring elements, especially their mean or jump. More...
 
class  NeumannTraceSpace3d
 A the NeumannTrace space of a given 3D - Finite Element space. More...
 
class  Postprocess14
 Computes x to the power of 1/4. More...
 
class  Postprocess3
 Computes x to the power of 3. More...
 
class  Postprocess34
 Computes x to the power of 3/4. More...
 
class  Postprocess4
 Computes x to the power of 4. More...
 
class  PostprocessRoot
 Returns the square root of x. More...
 
class  PostprocessRoot4
 Returns the forth root of x. More...
 
class  PostprocessSquare
 Squares the given x. More...
 
class  ProductOfAll
 Computes the product of all distances in the singular set singularities. More...
 
class  RefineOrRaise
 Refines element or raises its polynomial degree. More...
 
class  Riesz
 Linear form in 3D. More...
 
class  RotRot
 A function class to calculate element matrices for the Rot u*Rot v Bilinearform. More...
 
class  ShapeFunction3D
 Collecting the data of a 3D shape function in one class. More...
 
class  ShortestDist
 A function class for weighted regularization, which returns the square of the shortest distance of a point to the singular edges and vertices. More...
 
class  ShortestDistLimited
 A function class for weighted regularization, which returns the minimum of a value and the square of the shortest distance of a point to the singular edges and vertices. More...
 
class  SingularEdge
 Class for storing a singular edge with coordinates of its corners. More...
 
class  SingularSet
 Class for handling a set of singular edges and vertices. More...
 
class  SingularVertex
 Class for storing a singular vertex with the coordinates. More...
 
class  Space
 A 3D hp FEM space with continuous, picewise polynomial basis functions. More...
 
class  SpaceTransition
 Maps a solution vector from one space to another. More...
 
class  SumFactorization
 Sum factorization for an element matrix. More...
 
class  TraceSpace
 Builds the trace space of a FE space consisting of hexahedral elements. More...
 
class  TransitionPair
 Element pair for the SpaceTransition. More...
 
class  TrivialWeight
 A function class for trivial (constant equal 1.0) weight function. More...
 
class  Value
 The approximated function in a FE space. More...
 
class  ZeroTangentialValue
 Sets the outer product of the coefficients and the normal vector to zero in every node on an edge with chosen attribute. More...
 

Enumerations

enum  partDerivType { NO_DERIV = 0, X_DERIV = 1, Y_DERIV = 2, Z_DERIV = 3 }
 

Functions

template<class F >
ShapeFunction3D< F > makeShapeFunction3D (const Hexahedron &elm)
 
std::ostream & operator<< (std::ostream &os, const DaugeWeight &p)
 
template<typename DistClass , typename Function >
std::ostream & operator<< (std::ostream &os, const DistancePost< DistClass, Function > &p)
 
std::ostream & operator<< (std::ostream &os, const Postprocess14 &p)
 
std::ostream & operator<< (std::ostream &os, const Postprocess3 &p)
 
std::ostream & operator<< (std::ostream &os, const Postprocess34 &p)
 
std::ostream & operator<< (std::ostream &os, const Postprocess4 &p)
 
std::ostream & operator<< (std::ostream &os, const PostprocessRoot &p)
 
std::ostream & operator<< (std::ostream &os, const PostprocessRoot4 &p)
 
std::ostream & operator<< (std::ostream &os, const PostprocessSquare &p)
 
std::ostream & operator<< (std::ostream &os, const ProductOfAll &p)
 
std::ostream & operator<< (std::ostream &os, const ShortestDist &p)
 
std::ostream & operator<< (std::ostream &os, const ShortestDistLimited &p)
 
std::ostream & operator<< (std::ostream &os, const TrivialWeight &p)
 
template<class F >
void setupAdvection (vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf, const concepts::ElementFormulaContainer< concepts::Point< F, 3 > > frm)
 Function to setup a bilinear form related to the vector Advection, namely. More...
 
template<class F >
void setupIdentity (vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf)
 Function to setup a bilinear form related to the vector Identity, namely. More...
 
template<class F >
void setupIdentity (vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &bf, const concepts::ElementFormulaContainer< concepts::Mapping< F, 3 > > frm, bool transposed=false)
 Function to setup a bilinear form related to the vector Identity, namely. More...
 

Detailed Description

3D hp-FEM for H1-conforming elements.

The Space can be built using full tensor product polynomial spaces in the elements or trunk spaces (and much more) by changing the way the degrees of freedom are built in the BuildDofsBase class (and its specializations).

Author
Philipp Frauenfelder, 2001
Todo:
Three dimensional hp-DGFEM in hexahedra for locally refined meshes.

Enumeration Type Documentation

◆ partDerivType

Enumerator
NO_DERIV 
X_DERIV 
Y_DERIV 
Z_DERIV 

Definition at line 31 of file bf3D_partialDeriv.hh.

Function Documentation

◆ makeShapeFunction3D()

template<class F >
ShapeFunction3D<F> hp3D::makeShapeFunction3D ( const Hexahedron elm)

◆ operator<<() [1/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const DaugeWeight p 
)

◆ operator<<() [2/13]

template<typename DistClass , typename Function >
std::ostream& hp3D::operator<< ( std::ostream &  os,
const DistancePost< DistClass, Function > &  p 
)
inline

Definition at line 175 of file shortestDist.hh.

◆ operator<<() [3/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const Postprocess14 p 
)

◆ operator<<() [4/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const Postprocess3 p 
)

◆ operator<<() [5/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const Postprocess34 p 
)

◆ operator<<() [6/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const Postprocess4 p 
)

◆ operator<<() [7/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const PostprocessRoot p 
)

◆ operator<<() [8/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const PostprocessRoot4 p 
)

◆ operator<<() [9/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const PostprocessSquare p 
)

◆ operator<<() [10/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const ProductOfAll p 
)

◆ operator<<() [11/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const ShortestDist p 
)

◆ operator<<() [12/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const ShortestDistLimited p 
)

◆ operator<<() [13/13]

std::ostream& hp3D::operator<< ( std::ostream &  os,
const TrivialWeight p 
)

◆ setupAdvection()

template<class F >
void hp3D::setupAdvection ( vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &  bf,
const concepts::ElementFormulaContainer< concepts::Point< F, 3 > >  frm 
)

Function to setup a bilinear form related to the vector Advection, namely.

\[\sum_{i,j=1}^3 \int_K w_j u_i\frac{\partial v_i}{\partial_{x_j}}\mathrm{d}x = \sum_{i=1}^3 \int_K \underline{w}\cdot\nabla v_i u_i \mathrm{d}x = \int_K \underline{w}^\top\nabla\underline{v}\underline{u}\mathrm{d}x\]

on vectorial spaces, where $\underline{u}$ and $\underline{v}$ are trial and test functions respectively, $\underline{w}$ is a given vectorial formula, and

\[ \nabla \underline{v} = \left(\begin{array}{ccc} \frac{\partial u_1}{\partial x_1} & \frac{\partial u_1}{\partial x_2} & \frac{\partial u_1}{\partial x_3} \\ \frac{\partial u_2}{\partial x_1} & \frac{\partial u_2}{\partial x_2} & \frac{\partial u_2}{\partial x_3} \\ \frac{\partial u_3}{\partial x_1} & \frac{\partial u_3}{\partial x_2} & \frac{\partial u_3}{\partial x_3}. \end{array}\right) \]

◆ setupIdentity() [1/2]

template<class F >
void hp3D::setupIdentity ( vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &  bf)

Function to setup a bilinear form related to the vector Identity, namely.

\[\int_K \underline{u}\cdot\underline{v}\mathrm{d}x\]

on vectorial spaces.

◆ setupIdentity() [2/2]

template<class F >
void hp3D::setupIdentity ( vectorial::BilinearForm< F, typename concepts::Realtype< F >::type > &  bf,
const concepts::ElementFormulaContainer< concepts::Mapping< F, 3 > >  frm,
bool  transposed = false 
)

Function to setup a bilinear form related to the vector Identity, namely.

\[\int_K \underline{u}^\top A\underline{v}\mathrm{d}x\]

wih a matrix function A on vectorial spaces.

If transposed is true then the transposed of A is taken instead.

Page URL: http://wiki.math.ethz.ch/bin/view/Concepts/WebHome
21 August 2020
© 2020 Eidgenössische Technische Hochschule Zürich