Available Materials in NeoPZ

Note

There are many deprecated weak formulations in the folder Material/needrefactor. In case you need one of them, please feel free to contact us at <neopz@googlegroups.com> or to adapt them and submit a pull request.

Note

If you write a material that you think it should be available in NeoPZ, please feel free to submit a pull request.

Projections Group

This group has materials related projecting a given solution into a discrete approximation space. These classes can be used to see how to interact with H1(Ω), H(div,Ω) and H(curl,Ω) approximation spaces.

TPZL2Projection

template<class TVar = STATE>
class TPZL2Projection : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatErrorSingleSpace<STATE>>

Implements an L2 projection of a given solution on a scalar approximation space.

This material uses a scalar approximation space in a geometric domain of dimension defined by the user (1, 2 or 3D). It can project several solutions at once, but for calculating the error, it will only consider one solution at a time. The solutions to be projected are given by the Forcing Function of the TPZMaterialT. It is expected to be a TPZVec<TVar> of size nsol. For the error analysis, the solution is set through the ExactSol of TPZMatErrorSingleSpace.

Public Types

enum ESolutionVars

Solution indices of post-processing.

Values:

enumerator ENone
enumerator ESolution
enumerator EDerivative

Public Functions

TPZL2Projection() = default

Default constructor.

TPZL2Projection(int id, int dim, int nstate = 1)

Class constructor.

Parameters
  • id – material id

  • dim – problem dimension

  • nstate – number of state variables

TPZL2Projection(int id, int dim, int nstate, const TPZVec<TVar> &sol)

Class constructor setting a default constant solution.

Parameters
  • id – material id

  • dim – problem dimension

  • nstate – number of state variables

  • sol – constant solution. Ignored if a forcing function is set.

inline virtual std::string Name() const override

Returns the name of the material.

inline void SetScaleFactor(REAL scale)

Set a scale factor for the stiffness matrix and right hand side the default value of the scale factor is 1.

inline virtual int Dimension() const override

Returns the integrable dimension of the material.

inline virtual void SetDimension(int dim)

Sets problem dimension.

inline virtual int NStateVariables() const override

Returns the number of state variables associated with the material.

virtual TPZMaterial *NewMaterial() const override

To create another material of the same type.

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

It returns the variable index associated with the name.

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

virtual int ClassId() const override

Unique identifier for each class.

TPZHDivProjection

template<class TVar = STATE>
class TPZHDivProjection : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatErrorSingleSpace<STATE>>

Implements a projection of a given solution on a HDiv approximation space.

This material uses a HDiv approximation space in a geometric domain of dimension defined by the user (2 or 3D) and projects a given solution using the HDiv norm. The solution to be projected is set by the Forcing Function of the TPZMaterialT. It is expected to be a TPZVec<TVar> of size 4, as the solution is always in a 3D space and the last position is the divergence. For the error analysis, the solution is set through the ExactSol of TPZMatErrorSingleSpace.

Public Types

enum ESolutionVars

Solution indices of post-processing.

Values:

enumerator ENone
enumerator ESolution
enumerator EDivergence

Public Functions

TPZHDivProjection()

Default constructor.

TPZHDivProjection(int id, int dim)

Creates the TPZHDivProjection material.

Parameters
  • id – material identifier

  • dim – dimension

inline virtual std::string Name() const override

Returns the name of the material.

inline void SetScaleFactor(REAL scale)

Set a scale factor for the stiffness matrix and right hand side the default value of the scale factor is 1.

inline virtual int Dimension() const override

Returns the integrable dimension of the material.

inline virtual void SetDimension(int dim)

Sets problem dimension.

inline virtual int NStateVariables() const override

Returns the number of state variables associated with the material.

virtual TPZMaterial *NewMaterial() const override

To create another material of the same type.

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

It returns the variable index associated with the name.

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

virtual int ClassId() const override

Unique identifier for each class.

TPZHCurlProjection

template<class TVar = STATE>
class TPZHCurlProjection : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatErrorSingleSpace<STATE>>

Implements a projection of a given solution on a HCurl approximation space.

This material uses a HCurl approximation space in a geometric domain of dimension defined by the user (2 or 3D) and projects a given solution using the HCurl norm. The solution to be projected is set by the Forcing Function of the TPZMaterialT. It is expected to be a TPZVec<TVar> of size 3+ dim curl, as the solution is always in a 3D space and the last positions are the curl. For the error analysis, the solution is set through the ExactSol of TPZMatErrorSingleSpace.

Public Types

enum ESolutionVars

Solution indices of post-processing.

Values:

enumerator ENone
enumerator ESolution
enumerator ECurl

Public Functions

TPZHCurlProjection()

Default constructor.

TPZHCurlProjection(int id, int dim)

Creates the TPZHCurlProjection material.

Parameters
  • id – material identifier

  • dim – dimension

inline virtual std::string Name() const override

Returns the name of the material.

inline void SetScaleFactor(REAL scale)

Set a scale factor for the stiffness matrix and right hand side the default value of the scale factor is 1.

inline virtual int Dimension() const override

Returns the integrable dimension of the material.

inline void SetDimension(int dim)

Sets problem dimension.

inline virtual int NStateVariables() const override

Returns the number of state variables associated with the material.

virtual TPZMaterial *NewMaterial() const override

To create another material of the same type.

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

It returns the variable index associated with the name.

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

virtual int ClassId() const override

Unique identifier for each class.

The Lagrange Multiplier Group

TPZLagrangeMultiplier

template<class TVar = STATE>
class TPZLagrangeMultiplier : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatInterfaceSingleSpace<STATE>>, public TPZLagrangeMultiplierBase

Material which implements a Lagrange Multiplier for single space materials.

Lagrange

inline virtual void SetMultiplier(const TVar mult)

Sets the multiplier.

inline TVar Multiplier()

Gets the multiplier.

TPZLagrangeMultiplierCS

template<class TVar = STATE>
class TPZLagrangeMultiplierCS : public TPZMatBase<STATE, TPZMatCombinedSpacesT<STATE>, TPZMatInterfaceCombinedSpaces<STATE>>, public TPZLagrangeMultiplierBase

Material which implements a Lagrange Multiplier for combined spaces.

Lagrange

inline virtual void SetMultiplier(const TVar mult)

Sets the multiplier.

inline TVar Multiplier()

Gets the multiplier.

TPZNullMaterial

When dealing with combined approximation spaces, there are atomic FEM meshes associated with a multiphysic mesh. In these cases, one needs a dummy material in the atomic meshes, since no computation needs to be done. This can be easily done by using the following class

template<class TVar = STATE>
class TPZNullMaterial : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>>

Poisson

Table of Contents

Material for solving the Poisson Equation

divu=f

for one, two and three-dimensional domains. The right hand side f can be a function of the spatial coordinates, and the class can solve for multiple right hand sides at once, and it can perform error analysis based on a given solution.

TPZMatPoisson

template<class TVar = STATE>
class TPZMatPoisson : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatErrorSingleSpace<STATE>, TPZMatLoadCases<STATE>>

Implements a H1 formulation of the poisson equation.

This material uses a scalar approximation space in a geometric domain of dimension defined by the user (1, 2 or 3D). It can solve for multiple values of rhs at once, but for calculating the error, it will only consider one solution at a time.

Public Types

enum ESolutionVars

Solution indices of post-processing.

Values:

enumerator ENone
enumerator ESolution
enumerator EDerivative

Public Functions

TPZMatPoisson() = default

Default constructor.

TPZMatPoisson(int id, int dim)

Class constructor.

Parameters
  • id – material id

  • dim – problem dimension

  • nstate – number of state variables

inline virtual std::string Name() const override

Returns the name of the material.

inline void SetScaleFactor(REAL scale)

Set a scale factor for the stiffness matrix and right hand side the default value of the scale factor is 1.

inline virtual int Dimension() const override

Returns the integrable dimension of the material.

inline virtual int NStateVariables() const override

Returns the number of state variables associated with the material.

inline virtual void SetDimension(int dim)

Sets problem dimension.

virtual TPZMaterial *NewMaterial() const override

To create another material of the same type.

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

It returns the variable index associated with the name.

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

virtual int ClassId() const override

Unique identifier for each class.

Elasticity Group

Table of Contents

This group has materials related to problems involving elasticity.

TPZElasticity2D

class TPZElasticity2D : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatErrorSingleSpace<STATE>, TPZMatLoadCases<STATE>>

This class implements a two dimensional elastic material in plane stress or strain.

Elasticity

inline void SetElasticity(STATE E, STATE nu)

Set elasticity parameters.

inline void SetElasticityFunction(ElasticityFunctionType func)

Set a variable elasticity and poisson coefficient.

inline void SetPlaneStrain()

Set the material configuration to plane strain.

inline void SetPlaneStress()

Set the material configuration to plane stress.

inline void SetBodyForce(STATE fx, STATE fy)

Set forcing function.

inline STATE E()

Returns the elasticity modulus E.

inline STATE Nu()

Returns the poison coefficient modulus E.

void SetPreStress(STATE Sigxx, STATE Sigyy, STATE Sigxy, STATE Sigzz)

Set PresStress Tensor.

Contribute methods

void Contribute(const TPZMaterialDataT<STATE> &data, STATE weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override

Calculates the element stiffness matrix.

void ContributeBC(const TPZMaterialDataT<STATE> &data, STATE weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCondT<STATE> &bc) override

Applies the element boundary conditions.

Solution

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

Returns the variable index associated with the name.

virtual int NSolutionVariables(int var) const override

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

void Solution(const TPZMaterialDataT<STATE> &data, int var, TPZVec<STATE> &Solout) override

Returns the solution associated with the var index based on the finite element approximation.

Errors

inline int NEvalErrors() const override

Returns the number of norm errors.

Default is 3: energy, L2 and H1.

ReadWrite

virtual int ClassId() const override

Read and Write methods.

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

TPZElasticity2D(int id, STATE E, STATE nu, STATE fx, STATE fy, int planestress = 1)

Creates an elastic material with:

Parameters
  • id – material id

  • E – elasticity modulus

  • nu – poisson coefficient

  • fx – forcing function x=fx

  • fy – forcing function y=fy

  • plainstressplainstress=1 indicates use of plainstress

inline virtual int Dimension() const override

Returns the model dimension.

virtual int NStateVariables() const override

Returns the number of state variables associated with the material.

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

Print the material data.

inline virtual std::string Name() const override

Returns the material name.

virtual TPZMaterial *NewMaterial() const override

Creates a new material from the current object ??

TPZElasticity3D

class TPZElasticity3D : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatErrorSingleSpace<STATE>, TPZMatLoadCases<STATE>>

This class implements a 3D isotropic elasticity material.

Since

Aug 31, 2005.

Public Functions

TPZElasticity3D(int nummat, STATE E, STATE poisson, TPZVec<STATE> &force, STATE preStressXX = 0., STATE preStressYY = 0., STATE preStressZZ = 0.)

Class constructor.

Parameters
  • nummat – - material ID.

  • E – - Young’s modulus.

  • poisson – - poisson’s ratio

  • force – - external forces

TPZElasticity3D(int nummat)

Constructor.

Parameters

nummat – - material ID.

TPZElasticity3D()

Default constructor.

TPZElasticity3D(const TPZElasticity3D &cp)

Copy constructor.

virtual ~TPZElasticity3D()

Default destructor.

inline virtual int Dimension() const override

Returns model dimension.

inline virtual int NStateVariables() const override

Number of state variables.

inline virtual int IntegrationRuleOrder(int elPMaxOrder) const override

Gets the order of the integration rule necessary to integrate an element with polinomial order p.

virtual void Print(std::ostream &out) const override

Print material report.

inline void SetPostProcessingDirection(TPZVec<REAL> &Direction)

Direction to post process stress and strain.

Result of post processing is (Stress.Direction) or (Strain.Direction)

inline virtual std::string Name() const override

Material name.

void ContributeBC(const TPZMaterialDataT<STATE> &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCondT<STATE> &bc) override

Implements Dirichlet and Neumann boundary conditions.

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

Returns index of post-processing variable.

virtual int NSolutionVariables(int var) const override

Number of data of variable var.

void Solution(const TPZMaterialDataT<STATE> &data, int var, TPZVec<STATE> &Solout) override

Returns the solution associated with the var index based on the finite element approximation.

inline virtual int NEvalErrors() const override

Returns the number of norm errors: 3 (Semi H1, L2 and H1)

void FillDataRequirements(TPZMaterialData &data) const override

Fill material data parameter with necessary requirements for the Contribute method.

inline void FillBoundaryConditionDataRequirements(int type, TPZMaterialData &data) const override

This method defines which parameters need to be initialized in order to compute the contribution of the boundary condition.

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

Saves the element data to a stream.

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

Reads the element data from a stream.

virtual int ClassId() const override

Unique identifier for each class.

inline virtual TPZMaterial *NewMaterial() const override

Creates a new material from the current object ??

Electromagnetics Group

This group has materials related to electromagnetic problems. For usage of these materials, the reader is referred to the wgma library

Scattering Analysis of Planar Waveguides

The following classes implement a H1(Ω) formulation for the scattering analysis of planar waveguides.

Assuming that

  • the computational domain is in the yz plane

  • the direction of no variation of the waveguide is x

  • at the incidence plane z=z, denoted by Γ there is a source fs:fs(y) that propagates with ejβ(zz) in the +z direction and with ejβ(zz) in the z direction (such as a mode field obtained by a previous one-dimensional modal analysis, ffor instance),

this material implements:

Find non-trivial {ϕ}[C]N such that:

[A]{ϕ}=[f],

where

[A]ij=Ω[czϕizϕjz+cyϕiyϕjyk0csϕiϕj]dΩ[f]i=2Γczϕijβfs(y)dΓ

where we have, for both TE and TM modes, the following meanings:

TE

TM

ϕ

Ex

Hx

cz

μr,yy1

ϵr,yy1

cy

μr,zz1

ϵr,zz1

cs

μr,xx1

ϵr,xx1

Note that the media are isotropic everywhere except in the PML regions, which justify the subscripts in the table above.

This formulation was presented originally in this work 3, and it is specially interesting when analysing waveguide discontinuities.

Note

For simplicity, in the implementation, the computational domain is assumed to be in the xy plane (source propagates in the x direction).

TPZPlanarWGScattering

Implements the scattering analysis of a planar waveguide.

class TPZPlanarWGScattering : public TPZMatBase<CSTATE, TPZMatSingleSpaceT<CSTATE>>

This class implements the weak statement for the scattering analysis of planar waveguides using H1 elements.

It can either approximate TE or TM modes of a 2D waveguide with a 1D cross-section. Currently, only sources aligned with the y direction of a 2D domain can be used. This source is implemented through a ForcingFunctionBCType. The matrix-type parameter is not used, and the vector-type is expected to contain source-value in the first position and beta-value in the second position.

Note

Formulation taken from: Yasuhide Tsuji and Masanori Koshiba, “Finite Element Method Using Port Truncation by Perfectly Matched Layer Boundary Conditions for Optical Waveguide Discontinuity Problems,” J. Lightwave Technol. 20, 463- (2002)

Subclassed by TPZPlanarWGScatteringPML

ParamMethods

void SetWavelength(STATE lambda)

Sets the wavelength being analysed.

inline STATE GetWavelength() const

Gets the current wavelength.

void SetPermeability(const CSTATE ur)

Sets the permeability of the material.

virtual TPZVec<CSTATE> GetPermeability(const TPZVec<REAL> &x) const

Gets the permeability of the material.

void SetPermittivity(const CSTATE er)

Sets the permittivity of the material.

virtual TPZVec<CSTATE> GetPermittivity(const TPZVec<REAL> &x) const

Gets the permittivity of the material.

SolutionMethods

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

Variable index of a given solution.

virtual int NSolutionVariables(int var) const override

Number of variables associated with a given solution.

void Solution(const TPZMaterialDataT<CSTATE> &data, int var, TPZVec<CSTATE> &solout) override

Computes the solution at an integration point.

Public Types

enum ModeType

Sets for which mode type the problem will be solved.

Values:

enumerator TE
enumerator TM

Public Functions

TPZPlanarWGScattering(int id, const CSTATE ur, const CSTATE er, const STATE lambda, const ModeType mode, const REAL &scale = 1.)

Constructor taking a few material parameters.

Note

the scale param might help with floating point arithmetics on really small domains.

Parameters
  • id[in] Material identifier.

  • ur[in] Relative permeability.

  • er[in] Relative permittivity.

  • scale[in] Scale for geometric domain.

virtual TPZPlanarWGScattering *NewMaterial() const override

Create another material of the same type.

inline virtual std::string Name() const override

Returns the name of the material.

inline virtual int Dimension() const override

Returns the integrable dimension of the material.

inline virtual int NStateVariables() const override

Returns the number of state variables associated with the material.

TPZPlanarWGScatteringPML

Implements a PML for the planar waveguide scattering analysis.

class TPZPlanarWGScatteringPML : public TPZPlanarWGScattering

Public Functions

TPZPlanarWGScatteringPML(const int id, const TPZPlanarWGScattering &mat)

Creates PML based on another domain region.

void SetAttX(const REAL pmlBegin, const STATE alpha, const REAL d)

Sets information regarding the attenuation of the PML in the x-direction.

void SetAttY(const REAL pmlBegin, const STATE alpha, const REAL d)

Sets information regarding the attenuation of the PML in the y-direction.

virtual TPZVec<CSTATE> GetPermeability(const TPZVec<REAL> &x) const override

Gets the permeability of the material.

virtual TPZVec<CSTATE> GetPermittivity(const TPZVec<REAL> &x) const override

Gets the permittivity of the material.

virtual TPZPlanarWGScatteringPML *NewMaterial() const override

Create another material of the same type.

inline virtual std::string Name() const override

Returns the name of the material.

1

J.-F. Lee, D.-K. Sun, and Z.J. Cendes. Full-wave analysis of dielectric waveguides using tangential vector finite elements. IEEE Transactions on Microwave Theory and Techniques, 39(8):1262–1271, 1991.

2

Francisco T. Orlandini. A study on the construction of $H(curl,\Omega )$- elements and their application on waveguide problem. Master’s thesis, Universidade Estadual de Campinas, 2018.

3

Yasuhide Tsuji and Masanori Koshiba. Finite element method using port truncation by perfectly matched layer boundary conditions for optical waveguide discontinuity problems. J. Lightwave Technol., 20(3):463, Mar 2002. URL: http://opg.optica.org/jlt/abstract.cfm?URI=jlt-20-3-463.

Darcy flow

Materials for solving the Darcy flow equation:

σ=fσ=Ku

TPZDarcyFlow

class TPZDarcyFlow : public TPZMatBase<STATE, TPZMatSingleSpaceT<STATE>, TPZMatErrorSingleSpace<STATE>, TPZIsotropicPermeability>

This class implements an H1-conforming approximation for the Darcy flow equation for isotropic materials.

The Darcy flow equation is given by:

Ku=f,
where u is the pressure field to be solved, K is the permeability tensor and f is the source term.

See

TPZMixedDarcyFlow For an approximation using the mixed method.

See

TPZIsotropicPermeability For setting the permeability field.

Subclassed by TPZHybridDarcyFlow

Public Functions

TPZDarcyFlow()

Default constructor.

TPZDarcyFlow(int id, int dim)

Class constructor.

Parameters
  • id[in] material id

  • dim[in] problem dimension

inline virtual std::string Name() const override

Returns a ‘std::string’ with the name of the material.

inline virtual int Dimension() const override

Returns the problem dimension.

inline virtual int NStateVariables() const override

Returns the number of state variables.

inline int NEvalErrors() const override

Returns the number of errors to be evaluated.

Returns the number of errors to be evaluated, that is, the number of error norms associated with the problem.

virtual void SetDimension(int dim)

Sets problem dimension.

void Contribute(const TPZMaterialDataT<STATE> &data, STATE weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override

It computes a contribution to the stiffness matrix and load vector at one integration point.

Parameters
  • data[in] stores all input data

  • weight[in] is the weight of the integration rule

  • ek[out] is the element matrix

  • ef[out] is the rhs vector

void ContributeBC(const TPZMaterialDataT<STATE> &data, STATE weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCondT<STATE> &bc) override

It computes a contribution to the stiffness matrix and load vector at one BC integration point.

Parameters
  • data[in] stores all input data

  • weight[in] is the weight of the integration rule

  • ek[out] is the element matrix

  • ef[out] is the rhs vector

  • bc[in] is the boundary condition material

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

Returns an integer associated with a post-processing variable name.

Parameters

name[in] string containing the name of the post-processing variable. Ex: “Pressure”.

virtual int NSolutionVariables(int var) const override

Returns an integer with the dimension of a post-processing variable.

Parameters

var[in] index of the post-processing variable, according to TPZDarcyFlow::VariableIndex method.

void Solution(const TPZMaterialDataT<STATE> &data, int var, TPZVec<STATE> &solOut) override

Returns the solution associated with the var index based on the finite element approximation at a point.

Parameters
  • data[in] material data associated with a given integration point

  • var[in] index of the variable to be calculated

  • solOut[out] vector to store the solution

void GetSolDimensions(uint64_t &u_len, uint64_t &du_row, uint64_t &du_col) const override

Get the dimensions of the solution for each state variable.

This will be used for initializing the corresponding TPZMaterialData

Parameters
  • u_len[out] solution dimension

  • du_row[out] number of rows of the derivative

  • du_col[out] number of columns of the derivative

void Errors(const TPZMaterialDataT<STATE> &data, TPZVec<REAL> &errors) override

Calculates the approximation error at a point.

Parameters
  • data[in] material data of the integration point

  • errors[out] calculated errors

virtual int ClassId() const override

Returns an unique class identifier.

virtual TPZMaterial *NewMaterial() const override

Creates another material of the same type.

virtual void Print(std::ostream &out) const override

Prints data associated with the material.

TPZMixedDarcyFlow

class TPZMixedDarcyFlow : public TPZMatBase<STATE, TPZMatCombinedSpacesT<STATE>, TPZMatErrorCombinedSpaces<STATE>, TPZIsotropicPermeability>

This class implements a mixed approximation for the Darcy flow equation for isotropic materials.

The Darcy flow equation is given by:

σ=f,
where σ=Ku and u are the flux and pressure field to be solved respectively, K is the permeability tensor and f is the source term.

See

TPZDarcyFlow For an H1-conforming approximation.

See

TPZIsotropicPermeability For setting the permeability field.

Subclassed by TPZMixedDarcyFractureFlow

Public Functions

TPZMixedDarcyFlow()

Default constructor.

TPZMixedDarcyFlow(int id, int dim)

Class constructor.

Parameters
  • id[in] material id

  • dim[in] problem dimension

TPZMixedDarcyFlow(const TPZMixedDarcyFlow &copy)

copy constructor

TPZMixedDarcyFlow &operator=(const TPZMixedDarcyFlow &copy)

copy constructor

inline virtual std::string Name() const override

Returns a ‘std::string’ with the name of the material.

inline virtual int Dimension() const override

Returns the problem dimension.

inline virtual int NStateVariables() const override

Returns the number of state variables.

inline int NEvalErrors() const override

Returns the number of errors to be evaluated.

Returns the number of errors to be evaluated, that is, the number of error norms associated with the problem.

virtual void SetDimension(int dim)

Sets problem dimension.

void Contribute(const TPZVec<TPZMaterialDataT<STATE>> &datavec, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override

It computes a contribution to the stiffness matrix and load vector at one integration point.

Parameters
  • datavec[in] stores all input data

  • weight[in] is the weight of the integration rule

  • ek[out] is the element matrix

  • ef[out] is the rhs vector

void ContributeBC(const TPZVec<TPZMaterialDataT<STATE>> &datavec, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCondT<STATE> &bc) override

It computes a contribution to the stiffness matrix and load vector at one BC integration point.

Parameters
  • datavec[in] stores all input data

  • weight[in] is the weight of the integration rule

  • ek[out] is the element matrix

  • ef[out] is the rhs vector

  • bc[in] is the boundary condition material

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

Returns an integer associated with a post-processing variable name.

Parameters

name[in] string containing the name of the post-processing variable. Ex: “Pressure”.

virtual int NSolutionVariables(int var) const override

Returns an integer with the dimension of a post-processing variable.

Parameters

var[in] index of the post-processing variable, according to TPZDarcyFlow::VariableIndex method.

void Solution(const TPZVec<TPZMaterialDataT<STATE>> &datavec, int var, TPZVec<STATE> &solOut) override

Returns the solution associated with the var index based on the finite element approximation at a point.

Parameters
  • datavec[in] material data associated with a given integration point

  • var[in] index of the variable to be calculated

  • solOut[out] vector to store the solution

void Errors(const TPZVec<TPZMaterialDataT<STATE>> &data, TPZVec<REAL> &errors) override

Calculates the approximation error at a point.

Parameters
  • data[in] material data of the integration point

  • errors[out] calculated errors

virtual int ClassId() const override

Returns an unique class identifier.

virtual TPZMaterial *NewMaterial() const override

Creates another material of the same type.

virtual void Print(std::ostream &out) const override

Prints data associated with the material.

TPZIsotropicPermeability

class TPZIsotropicPermeability : public TPZSavable

This class implements the interface with the methods required to handle the permeability field of an isotropic material.

Subclassed by TPZIsotropicPermeabilityBC

Public Functions

void SetConstantPermeability(STATE constant)

Set a constant permeability to the material.

Parameters

constant[in] permeability value

void SetPermeabilityFunction(PermeabilityFunctionType &perm_function)

Set a varying permeability field to the material.

Parameters

perm_function[in] function that describes the permeability field

STATE GetPermeability(const TPZVec<REAL> &coord)

Return the permeability value at a coordinate.

Parameters

coord[in] coordinate of interest

virtual int ClassId() const override

Unique identifier for each class.

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

Reads an objects from the TPZStream buffer.

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

Writes this object to the TPZStream buffer.

Include the classid if withclassid = true