The Shape hierarchy

Approximation spaces

An essential ingredient of Galerkin finite element approximations if the definition of a finite dimensional subspace whose basis functions will be used to build the algebraic system of equations. The NeoPZ environment allows to user to define a variety of approximation spaces : \(H^1, H(div), H(curl)\) or discontinuous spaces

In NeoPZ the different flavors of approximation spaces are based on polynomial spaces associated with a master element. These basic “building blocks” of the polynomial spaces are each implemented in a class dedicated to each element topology. The shape template classes are :

The full (doxygen generated) documentation of the classes can be found at TPZShapeData.

Based on these classes a generic classes define \(H^1,\;H(div)\;H(curl)\) functions associated with each of the shape building blocks.

TPZShapeH1

the TPZShapeH1 class computes H1 compatible shape functions. The method Initialize takes as argument the polynomial order associated with each side of the element, an integer identifier associated with each node and initializes datastructures in the TPZShapeData data structure.

The method Shape fills in the data structure fPhi and fDPhi of the TPZShapeData object.

template<class TSHAPE>
struct TPZShapeH1

Class that implements the computation of H1 shape functions with variable connect order.

Subclassed by TPZShapeHDivKernel2D< TSHAPE >

TPZShapeHDiv

template<class TSHAPE>
struct TPZShapeHDiv

Traditional HDiv spaces, data structures that do not depend on the geometric map.

Subclassed by TPZShapeHDivCollapsed< TSHAPE >

TPZShapeHDivKernel

Warning

doxygenstruct: Cannot find class “TPZShapeHDivKernel” in doxygen xml output for project “NeoPZ” from directory: /Users/runner/work/neopz/neopz/build/docs/doxygen/xml

TPZShapeHDivConstant

template<class TSHAPE>
struct TPZShapeHDivConstant : public TPZShapeHCurlNoGrads<TSHAPE>

Traditional HDiv spaces, data structures that do not depend on the geometric map.

Public Static Functions

static void Initialize(TPZVec<int64_t> &ids, TPZVec<int> &connectorders, const TPZVec<int> &sideorient, TPZShapeData &data)

Should be called once per element. Initializes the data structure.

TPZShapeHCurl

template<class TSHAPE>
struct TPZShapeHCurl

Generates HCurl spaces in the reference element using constant vector fields and H1 shape functions.

Public Static Functions

static void Initialize(const TPZVec<int64_t> &ids, TPZVec<int> &connectorders, TPZShapeData &data)

Should be called once per element. Initializes the data structure.

static void ComputeVecandShape(TPZShapeData &data)

Computes the pair (vec, h1 shape) that are used for creating the HCurl shape functions.

static int ComputeNConnectShapeF(const int connect, const int order)

Computes the number of shape functions for a given connect.

Parameters
  • connect[in] which connect

  • order[in] connect ordder

static int NHCurlShapeF(const TPZShapeData &data)

Number of HCurl shape functions for a given (filled) TPZShapeData.

static void Shape(TPZVec<REAL> &pt, TPZShapeData &data, TPZFMatrix<REAL> &phi, TPZFMatrix<REAL> &curlphi)

Calculates phi and curl phi at a given integration point.

Parameters
  • pt[in] integration point in the reference element

  • shapedata[in]

static void StaticIndexShapeToVec(const TPZVec<int> &connectOrder, const TPZVec<int64_t> &firstH1ShapeFunc, const TPZVec<int> &sidesH1Ord, const TPZVec<int64_t> &nodeIds, TPZVec<unsigned int> &shapeCountVec, TPZVec<std::pair<int, int64_t>> &indexVecShape)

Returns a matrix index of the shape and vector associate to element.

Parameters
  • connectOrder[in] Order of the HCurl connects

  • sirstH1ShapeFunc[in] Index of the first H1 shape function associated with a given side

  • sidesH1Ord[in] H1 order of each side

  • nodeIds[in] Ids of the nodes. Used for adjusting order of HCurl functions.

  • shapeCountVec[out] Stores the number of shape function for each connect

  • indexVecShape[out] Indicates the pair vector/H1 shape function that will be used for constructing the HCurl functions

static void CalcH1ShapeOrders(const TPZVec<int> &ordHCurl, TPZVec<int> &ord)

Internal method for calculating the needed order of H1 functions for each side given the HCurl connect orders listed in ordHCurl param.

Note

since the H1 vertex functions are always needed and their order does not change, ord has size NSides-NCornerNodes.

Parameters
  • ordHCurl[in] The effective order of each HCurl connect associated with a the side of shape TSIDESHAPE

  • the[out] order of the h1 connects needed for computing the hcurl functions

static int MaxOrder(const int ordh1)

Maximum (actual) polynomial order. Useful for choosing integration rule.

TPZShapeData

class TPZShapeData : public virtual TPZSavable

This class implements a structure to store data regarding shape functions in the master element coordinate system.

Subclassed by TPZMaterialData

Public Types

enum MShapeFunctionType

Type of the shape function associated with an element

Values:

enumerator EEmpty
enumerator EScalarShape

Scalar shape functions (H1, L2)

enumerator EVecandShape

Composite shape function (scalar function and vector field, HDiv and HCurl spaces)

enumerator EVecShape

Vector shape function (e.g. vector H1 space)

Public Functions

TPZShapeData()

Default constructor.

TPZShapeData(const TPZShapeData&) = default

Copy constructor.

TPZShapeData(TPZShapeData&&) = default

Move constructor.

virtual ~TPZShapeData() = default

Destructor.

TPZShapeData &operator=(const TPZShapeData&) = default

Copy assignment operator.

TPZShapeData &operator=(TPZShapeData&&) = default

Move assignment operator.

virtual int ClassId() const override

Read and Write methods.

virtual void Write(TPZStream &buf, int withclassid) const override

Writes this object to the TPZStream buffer.

Include the classid if withclassid = true

virtual void Read(TPZStream &buf, void *context) override

Reads an objects from the TPZStream buffer.

std::string ShapeFunctionType() const

Shape function type as a string.

void Print(std::ostream &out) const

Prints the data.

void PrintMathematica(std::ostream &out) const

Prints the data in a format suitable for Mathematica.

Public Members

MShapeFunctionType fShapeType = {EEmpty}

Type of shape function.

TPZManVector<int64_t, 8> fCornerNodeIds

Corner node ids determine the parameter transformations to the sides.

TPZManVector<int, 27> fH1ConnectOrders

Connect orders determine the number of shape functions.

TPZManVector<int, 27> fH1NumConnectShape

Number of shape functions by connect.

TPZManVector<TPZTransform<REAL>, 20> fSideTransforms

Parametric transforms from the interior to the side.

TPZManVector<int, 20> fSideTransformationId

Transformation ids for each side.

TPZFNMatrix<MatDataNumPhi, REAL> fPhi

Vector of shapefunctions (format is dependent on the value of shapetype) over the master element.

TPZFNMatrix<3 * MatDataNumPhi, REAL> fDPhi

Values of the derivative of the shape functions over the master element.

TPZManVector<int, 27> fHDivConnectOrders

Connect orders determine the number of shape functions.

TPZManVector<int, 27> fHDivNumConnectShape

Number of shape functions by connect.

TPZManVector<int, 20> fSideOrient

The orientation of the sides (either -1 or +1)

TPZFNMatrix<3 * MatDataNumDir> fMasterDirections

Directions on the master element.

TPZManVector<std::pair<int, int64_t>, MatDataNumPhi> fSDVecShapeIndex

Correspondence between direction vector index and index of the shape functions. Used for H(div) and H(curl) approximation spaces.

The list of shape template classes

class TPZShapePoint : public pztopology::TPZPoint

Compute the single shape function associated with a point.

Shape

class TPZShapeLinear : public pztopology::TPZLine

Implements the shape functions of a linear (1D) element.

Shape

The linear shape functions form the basis of all other shape function computations

The range of the master element is -1,1

The orthogonal function which generates the linear shape functions can be modified by changing the function pointer fOrthogonal

all static tables and functions concerning one-d elements will be grouped in this

class TPZShapeTriang : public pztopology::TPZTriangle

Implements the shape functions of a triangular (2D) element.

Shape

The triangular shape functions are also used in 3D elements

The range of the master element is 0,1

class TPZShapeQuad : public pztopology::TPZQuadrilateral

Implements the shape functions of a quadrilateral (2D) element.

Shape

The quadrilateral shape functions are also used in 3D elements

The range of the master element is -1,1

class TPZShapeCube : public pztopology::TPZCube

Implements the shape functions of a hexahedral (3D) element.

Shape The range of the master element is [-1 ,1]

class TPZShapeTetra : public pztopology::TPZTetrahedron

Implements the shape functions of a tetrahedral (3D) element.

Shape The range of the master element is 0,1

class TPZShapePrism : public pztopology::TPZPrism

Implements the shape functions of a prism (3D) element.

Shape The range of the master element is \( [-1,1] \) for the \( z \) direction and \( [0,1] \) for the triangular plane

class TPZShapePiram : public pztopology::TPZPyramid

Implements the shape functions of a pyramid (3D) element.

Shape The range of the master element is -1,1 for the base and [0,1] for the height

Subclassed by pzshape::TPZShapePiramHdiv