The Material hierarchy

The Material Group

The TPZMaterial hierarchy is where the weak statement of a differential equation is defined in NeoPZ.

For a single integration point, the main role of the material is to compute the contribution to the FEM matrix and to the right hand side vector and, after the linear system is solved, the post-processing variables.

Different regions of the domain can have different materials, or maybe multiple instances of a same material with different properties (e.g. in scenarios where the domain is homogeneous by parts).

For a list of materials (weak statements) available in NeoPZ, please check Available Materials in NeoPZ.

The Material Hierarchy Explained

During a simulation, the computational element knows which approximation space (or which combination of approximation spaces) is being used. The element will then call its associated material for computations on the integration points. These calls vary wildly depending on the simulation being performed: using discontinuous approximation spaces one often has to perform calculations on the interface between elements, which is not needed for continuous spaces, for instance.

In order not to pollute the base class of the TPZMaterial hierarchy with every possible interface, the following inheritance scheme was devised:

class TPZMaterial : public virtual TPZSavable

This class has all the common type-agnostic interface that is available in a given material.

template<class TVar>
class TPZMaterialT : public TPZMaterial

This class complements the interfaces of TPZMaterial with some type-related interface, which depends if the material is solving a problem resulting in a real matrix TVar = STATE or a complex matrix TVar = CSTATE. It also introduces the Forcing Function (more about that later).

template<class TVar, class ...Interfaces>
class TPZMatBase : public TPZMaterialT<TVar>, public Interfaces

The TPZMatBase is a variadic class template, meaning that it can be instantiated with an arbitrary number of the Interfaces parameters. This will be used so one can create interfaces according to their needs. The class TPZMatBase is a subclass of TPZMaterialT<TVar> and of each Interfaces parameter. You can read more about these interfaces in the section Material Interfaces.

Note

TPZMatBase is a variadic class template. Therefore

auto *mat = dynamic_cast<TPZMatBase<STATE>*>(someMat);

is expected to fail, since TPZMatBase<STATE> and TPZMatBase<STATE, Interface1, Interface2> are different types. Cast to TPZMaterial is always possible, to TPZMaterialT<TVar> with the correct TVar type or to any of the interface classes.

There are two special interfaces that have their own detailed documentation, TPZMatSingleSpaceT and TPZMatCombinedSpacesT. These are the interfaces that define, respectively, the behaviour of materials implementing a formulation with only one approximation space and with combined approximation spaces. See The Single Space Material Group and The Combined Spaces Material Group for further details on these interfaces.

Forcing Function

The TPZMaterialT::ForcingFunction is a mechanism for position-dependent functions to be evaluated at an integration point. Its type is defined as:

template<class TVar>
using ForcingFunctionType = std::function<void(const TPZVec<REAL> &loc, TPZVec<TVar> &result)>

Alias for forcing function type.

Lambdas can be used for easily creating forcing functions as

/*let us assume that mat inherits from a TPZMaterialT<STATE>*/


auto forcingFunction = [](const TPZVec<REAL>&x, TPZVec<STATE>&u){
  u[0] = x[0]*x[0];
  u[1] = x[1]*x[1];
  u[2] = x[2]*x[2];
};

//pOrder is the suggested integration rule for the forcing function
constexpr int pOrder{2};
mat->SetForcingFunction(forcingFunction,pOrder);

The class TPZMaterialT contains the forcing function-related methods.

Further documentation on TPZMaterial Hierarchy

class TPZMaterial : public virtual TPZSavable

This abstract class is the basis of the TPZMaterial hierarchy.

It contains the type-agnostic interface that materials should implement, and it applies for materials using both single approximation space and combined approximation spaces (multiphysics materials).

Note

Actual materials should derive from TPZMatBase

Subclassed by TPZMaterialT< STATE >, TPZMaterialT< CSTATE >, TPZMaterialT< TVar >

Constructors

TPZMaterial() = default

Default constructor.

inline explicit TPZMaterial(int id)

Constructor taking material identifier.

TPZMaterial(const TPZMaterial&) = default

Copy constructor.

TPZMaterial(TPZMaterial&&) = default

Move constructor.

Utilities

inline virtual std::string Name() const

Returns the name of the material.

virtual int Dimension() const = 0

Returns the integrable dimension of the material.

inline virtual int Id() const

Returns the material identifier.

inline void SetId(int id)

Sets the material identifier.

virtual bool HasForcingFunction() const = 0

Whether a forcing function has been set.

virtual void Print(std::ostream &out = std::cout) const

Prints data associated with the material.

virtual TPZMaterial *NewMaterial() const

Create another material of the same type.

inline void SetBigNumber(REAL bignumber)

Change the value of the penalty constant.

inline REAL BigNumber()

Access method to the penaly constant.

PostProcess

virtual int NStateVariables() const = 0

Returns the number of state variables associated with the material.

virtual int VariableIndex(const std::string &name) const

Returns the variable index associated with a given name.

virtual int NSolutionVariables(int var) const

Returns the number of variables associated with the variable indexed by var.

Parameters

var – Index variable into the solution, is obtained by calling VariableIndex

template<class TVar>
class TPZMaterialT : public TPZMaterial

This class provides additional interfaces for materials.

While the TPZMaterial class had only the type-agnostic interface, this class has type information.

Note

Actual materials should derive from TPZMatBase

Subclassed by TPZBndCondBase< TVar, Interfaces >, TPZMatBase< TVar, Interfaces >

ForcingFunction

virtual void SetForcingFunction(ForcingFunctionType<TVar> f, int pOrder)

Sets a forcing function to the material.

This method is used for setting functions used in the weak formulation.

Parameters
  • f[in] Forcing function

  • pOrder[in] Suggested integration rule order.

inline virtual bool HasForcingFunction() const override

Whether a forcing function has been set.

inline int ForcingFunctionPOrder() const

Get forcing function polynomial order.

inline virtual const ForcingFunctionType<TVar> &ForcingFunction() const

Gets reference to the forcing function.

Public Functions

TPZMaterialT() = default

Default constructor.

inline explicit TPZMaterialT(int id)

Constructor taking material identifier.

virtual TPZBndCondT<TVar> *CreateBC(TPZMaterial *reference, int id, int type, const TPZFMatrix<TVar> &val1, const TPZVec<TVar> &val2)

Creates an associated boundary condition.

Parameters
  • reference[in] The volumetric material associated with the BC.

  • id[in] Boundary condition identifier.

  • type[in] Type of the boundary condition.

  • val1[in] Value to be set at the element matrix.

  • val2[in] Value to be set at the rhs vector.

virtual int ClassId() const override

Unique identifier for each class.

template<class TVar, class ...Interfaces>
class TPZMatBase : public TPZMaterialT<TVar>, public Interfaces

The class that Materials should inherit from. One mandatory interface is either TPZMatSingleSpaceT, for formulations with only one approximation space, or TPZMatCombinedSpacesT, for formulations that use multiple approximation spaces. The other interfaces will enhance the functionality of the material.

MandatoryOverride

virtual int Dimension() const override = 0

Returns the integrable dimension of the material.

virtual int NStateVariables() const override = 0

Returns the number of state variables associated with the material.

OptionalOverride

inline virtual std::string Name() const override

Returns the name of the material.

inline virtual int NSolutionVariables(int var) const override

Returns the number of variables associated with the variable indexed by var.

Parameters

var – Index variable into the solution, is obtained by calling VariableIndex

inline virtual int VariableIndex(const std::string &name) const override

Returns the variable index associated with a given name.

Material Data

At an integration point, the information available to the material is given by the TPZMaterialData hierarchy.

Note

See TPZShapeData for information regarding the reference element

These classes relate to the information that is needed to compute the element matrices. They contain the basis functions and its derivatives in the deformed element, as well as informations regarding the geometric transform that relates the master and the deformed elements.

class TPZMaterialData : public TPZShapeData, public virtual TPZSavable

This class implements a type-agnostic interface between TPZCompEl::CalcStiff and Contribute methods of the materials.

It requests to the material which attributes must be computed by the computational element and trigger their computation.

Subclassed by TPZMaterialDataT< TVar >

Flags

Flags indicating which attributes should be available when computing the FEM matrix.

bool fNeedsSol = {false}

Whether the solution is needed for computing the FEM matrix.

bool fNeedsNeighborSol = {false}

Whether the neighbour’s solution is needed for computing the FEM matrix.

bool fNeedsHSize = {false}

Whether the element size is needed for computing the FEM matrix.

bool fNeedsNeighborCenter = {false}

Whether the neighbour’s center coordinate is needed for computing the FEM matrix.

bool fNeedsDeformedDirectionsFad = {false}

Whether FAD directions are needed for computing the FEM matrix.

Can be used for HDiv and HCurl approximation spaces

bool fNeedsNormal = {false}

Whether the normal vector is needed for computing the FEM matrix.

bool fActiveApproxSpace = {true}

Whether the material has an active approximation space.

Data

Attributes to be computed in CalcStiff

TPZFNMatrix<MatDataNumPhi, REAL> phi

Vector of shapefunctions (format is dependent on the value of shapetype)

TPZFNMatrix<MatDataNumDPhi, REAL> dphix

Values of the derivative of the shape functions in the deformed element.

Note

This derivative is calculated in the element axes. For obtaining the derivative in the xyz-coordinate system, one should call TPZAxesTools<REAL>::Axes2XYZ(dphix, dphix_xyz, data.axes) in the contribute method.

TPZFNMatrix<MatDataNumPhi, REAL> divphi

Values of the divergence of the shape functions in the mapped element (only applicable to H(div) spaces)

TPZFNMatrix<MatDataNumPhi, REAL> curlphi

Values of the curl of the shape functions in the mapped element (only applicable to H(curl) spaces)

TPZFNMatrix<9, REAL> axes

Axes between the R3 space and the element’s coordinates. Used for the derivatives of the shape functions.

TPZFNMatrix<9, REAL> jacobian

Value of the jacobian at the integration point.

TPZFNMatrix<9, REAL> jacinv

Value of the inverse of the jacobian at the integration point.

TPZManVector<REAL, 3> normal

Normal to the element at the integration point.

TPZManVector<REAL, 3> x

Value of the coordinate at the integration point.

TPZManVector<REAL, 3> xParametric

Value of the coordinate at the integration point corresponding to the x-parametric coordinate (master element)

int p = {-1}

Maximum polinomial order of the shape functions.

REAL HSize

Measure of the size of the element.

REAL detjac

Determinant of the jacobian.

TPZManVector<REAL, 3> XCenter

Value of the coordinate at the center of the element.

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

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

TPZFNMatrix<MatDataNumDir> fDeformedDirections

Directions on the deformed element. Used for H(div) and H(curl) approximation spaces.

template<class TVar>
class TPZMaterialDataT : public TPZMaterialData

This class implements an interface between TPZCompEl::CalcStiff and Contribute methods of the materials.

Since it takes the solution type into account, it can store the solution at an integration point.

Data

Attributes to be computed in CalcStiff.

TPZSolVec<TVar> sol

Vector of the solutions at the integration point.

TPZGradSolVec<TVar> dsol

Vector of the derivatives of the solution at the integration point.

TPZSolVec<TVar> divsol

Vector of the divergence of the solution at the integration point (only of H(div) spaces)

TPZSolVec<TVar> curlsol

Vector of the curl of the solution at the integration point (only of H(curl) spaces).

Boundary Conditions

The Boundary Condition Hierarchy explained

Usually there is a different contribution to be made by the materials associated with a boundary. For this purpose, there is a special hierarchy for boundary conditions, the TPZBndCond hierarchy. The boundary conditions are usually created from a given TPZMaterial and it will be of the type TPZBndCondBase, a variadic class template with as many InterfaceBC parameters as are present in original material. Using the CreateBC interface therefore ensures that the boundary condition will contain all the interfaces present in the original material.

The programmer will note that there is a unique boundary interface associated with each interface. The interface type is declared by a “using InterfaceBC” statement declared in each material interface

class TPZBndCond : public virtual TPZSavable
template<class TVar>
class TPZBndCondT : public TPZBndCond
template<class TVar, class ...Interfaces>
class TPZBndCondBase : public TPZMaterialT<TVar>, public TPZBndCondT<TVar>, public Interfaces

Creating Boundary Conditions

The following snippet illustrates how to create a boundary condition

/*
* mat: pointer to a given material
*      representing a weak form using only real variables
* cmesh: pointer to the computational mesh
*/

/*
val1 is the matrix of values that goes on the FEM matrix.
In this example, it is a 1x1 matrix
*/
TPZFMatrix<STATE> v1(1,1,0.);
//val2 is the vector of values that goes on the rhs
TPZVec<STATE> v2 = {0};
/*
boundary identifier: each identifier should be unique in a given mesh
*/
const int bcIdentifier{-1};
/*
boundary condition type: In the available NeoPZ materials,
                         Dirichlet=0,
                         Neumann=1,
                         Robin=2
                         but in your own materials any value can be used
*/
const int bcType = 0;

//Creates a boundary condition associated with mat
TPZBndCond * bnd = mat->CreateBC(mat,bcIdentifier,boundType,val1,val2);

//inserts into the computational mesh
cmesh->InsertMaterialObject(bnd)

Read the next section for a better comprehension of the val1 and val2 variables.

Note

The bcIdentifier should correspond to a (geometrical) mesh region (each element in the mesh are associated with a region through an identifier).

Note

The created boundary condition will automatically inherit any interface that mat has.

Implementing Boundary Conditions

Examples of boundary conditions can be seen in the Available Materials in NeoPZ, but for clarity, a quick example follows.

Taking the TPZMatSingleSpaceT class as an example, the boundary condition contribution is performed on

virtual void ContributeBC(const TPZMaterialDataT<TVar> &data, REAL weight, TPZFMatrix<TVar> &ek, TPZFMatrix<TVar> &ef, TPZBndCondT<TVar> &bc) = 0

The function parameter bc has all the information regarding the boundary condition. We will give special attention to four methods that are specially useful for computing contributions to boundaries:

inline virtual bool HasForcingFunctionBC() const final
inline const TPZFMatrix<TVar> &Val1() const
inline const TPZVec<TVar> &Val2() const
inline REAL BigNumber()

With HasForcingFunctionBC, one can check if the boundary condition has an associated Forcing Function, which is a position-dependent function that might be used to specify boundary conditions. Its signature is:

template<class TVar>
using ForcingFunctionBCType = std::function<void(const TPZVec<REAL> &loc, TPZVec<TVar> &rhsVal, TPZFMatrix<TVar> &matVal)>

Alias for boundary condition forcing function type.

When forcing functions are not being used, than Val1 and Val2 should be used for extracting the user-prescribed values for the boundary conditions. A simple example to illustrate their usage follows:

 /*
 * This is inside the ContributeBC method
 * we suppose that the material is expecting
 * a vec of size 3 for val 2
 * and a 3x3 matrix for val1
 */

 //TPZManVector has static and dynamic storage
 TPZManVector<TVar,3> v2(3,0.);
 //TPZFNMatrix is a matrix equivalent of TPZManVector
 TFNMatrix<9,TVar> v1 (3,3,0.);

 if(bc.HasForcingFunctionBC()){
   bc.ForcingFunctionBC()(data.x,v2,v1);
 } else {
   v2 = bc.Val2();
   v1 = bc.Val1();
 }

 /*example of dirichlet BC using big number*/
 const auto &big = TPZMaterial::BigNumber();
 switch (bc.Type()){
 // Dirichlet condition
 case 0 : {
   for(auto in = 0 ; in < phr; in++) {
     ef(in,0) += big * v2[0] * phi.GetVal(in,0) * weight;
       for (auto jn = 0 ; jn < phr; jn++) {
         ek(in,jb) += big * phi(in,0) * phi(jn,0) * weight;
       }//jn
     }//in
   break;
}
//more code...

Note

Homogeneous Dirichlet boundary conditions can also be implemented by making use of the method

inline void FilterEquations(TPZVec<int64_t> &origindex, TPZVec<int64_t> &destindex) const

This is a slightly more advanced usage. As soon as there is a tutorial or quick example on that, this section will be updated.

Further documentation on TPZBndCond hierarchy

class TPZBndCond : public virtual TPZSavable

This class defines a type-agnostic interface for boundary condition in TPZMaterial objects.

Note

Type-related methods are declared and defined in TPZBndCondT

Subclassed by TPZBndCondT< TVar >

Public Functions

virtual void SetMaterial(TPZMaterial *mat) = 0

Set associated material. It also gets its attributes.

virtual TPZMaterial *Material()

Get associated material.

inline int Type() const

Get type of BC.

inline void SetType(int type)

Set type of BC.

virtual int Id() const = 0

Material identifier.

virtual int ClassId() const override

Unique identifier for each class.

virtual bool HasForcingFunctionBC() const = 0

Whether the boundary condition has a forcing function.

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

Reads an objects from the TPZStream buffer.

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

Writes this object to the TPZStream buffer.

Include the classid if withclassid = true

template<class TVar>
class TPZBndCondT : public TPZBndCond

This class defines an interface for boundary conditions.

It complements the interfaces defined in TPZBndCond with type information.

Note

Boundary conditions will actually be created as TPZBndCondBase

Subclassed by TPZBndCondBase< TVar, Interfaces >

ReadWrite

virtual int ClassId() const override

Unique identifier for each class.

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

Reads an objects from the TPZStream buffer.

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

Writes this object to the TPZStream buffer.

Include the classid if withclassid = true

Public Functions

inline void SetForcingFunctionBC(ForcingFunctionBCType<TVar> f)

Sets a forcing function to the boundary material.

This method is used for setting functions used in the weak formulation.

Parameters
  • f[in] Forcing function

  • pOrder[in] Suggested integration rule order.

inline virtual bool HasForcingFunctionBC() const final

Whether a forcing function has been set.

inline const ForcingFunctionBCType<TVar> &ForcingFunctionBC() const

Gets reference to the boundary condition forcing function.

inline const TPZFMatrix<TVar> &Val1() const

Gets the FEM matrix value for the boundary condition.

inline const TPZVec<TVar> &Val2() const

Gets the right hand side value for the boundary condition.

inline void SetVal1(const TPZFMatrix<TVar> &val1)

Sets the FEM matrix value for the boundary condition.

inline void SetVal2(const TPZVec<TVar> &val2)

Sets the right hand side value for the boundary condition.

template<class TVar, class ...Interfaces>
class TPZBndCondBase : public TPZMaterialT<TVar>, public TPZBndCondT<TVar>, public Interfaces

This class ensures that the boundary condition material is compatible with any given interfaces.

Public Functions

TPZBndCondBase(TPZMaterial *material, int id, int type, const TPZFMatrix<TVar> &val1, const TPZVec<TVar> &val2)

Constructor from an existing material setting fVal1 and fVal2

virtual int ClassId() const override

Unique identifier for each class.

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

Reads an objects from the TPZStream buffer.

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

Writes this object to the TPZStream buffer.

Include the classid if withclassid = true

virtual void SetMaterial(TPZMaterial*) final

Set associated material. It also gets its attributes.

inline virtual int Dimension() const final

Returns the integrable dimension of the material.

inline virtual int NStateVariables() const final

Returns the number of state variables associated with the material.

inline virtual int Id() const override

Returns the material identifier.

inline virtual TPZMaterial *NewMaterial() const override

Create another material of the same type.

Material Interfaces

Additional functionality may be incorporated to a material through an interface.

The most important interfaces are, as already mentioned, TPZMatSingleSpaceT and TPZMatCombinedSpacesT. These classes provide the interface for materials to perform computations using only one approximation space or a combination of approximation spaces, respectively. These classes share many similarities, as can be seen on The Single Space Material Group and The Combined Spaces Material Group.

However, a material can have an arbitrary number of interfaces. Examples of available interfaces are:

Interfaces are specific to the number of spaces of the material. For details on the interfaces available for each material group, see The Single Space Material Group and The Combined Spaces Material Group. The common interfaces are in Common Interfaces.

Creating custom material interfaces

Material interfaces can be easily created to adapt a material to one’s needs.

However, special care must be taken due to the materials associated with boundary conditions: In some cases the boundary condition material might need to forward function calls to its associated material, or sometimes it needs to be cast to the interface type.

For this to happen, an auxiliary class for the interface is needed, as it can be seen in the following example:

//forward declare bc auxiliary class
class TPZMyInterfaceBC;

class TPZMyInterface : public TPZSavable{
public:
  using TInterfaceBC = TPZMyInterfaceBC;
  //declare all the virtual methods and etc as usual
  //@{
  //! Read and Write methods
  int ClassId() const override;
  void Read(TPZStream&,int) override;
  void Write(TPZStream&,void*) override;
  //@}
};

//maybe this forward declaration is needed
class TPZMaterial;

class TPZMyInterfaceBC : public TPZMyInterface{
protected:
  TPZMyInterface *fMatMyInterface{nullptr};
  void SetMaterialImpl(TPZMaterial* mat){
    fMatMyInterface = dynamic_cast<TPZMyInterface*>(mat);
  }
public:
//overload the methods as needed.
};

In this example, an interface is given in the class TPZMyInterface. The rules for TPZMyInterface to be a valid interface are, as follows:

  • It must derive from TPZSavable.

  • It must have an associated boundary condition class, even if there is nothing to be done in this class.

  • The associated boundary condition class must be publicly avaliable under the alias TPZMyInterface::TInterfaceBC.

  • The associated boundary condition class must have a (protected) method void SetMaterialImpl(TPZMaterial *)(). This method can be useful if a casted pointer to the associated material is needed when forwarding calls.

  • The read and write methods must be overriden.