NeoPZ
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
TPZSwelling Class Reference

The TPZSwelling class implements a numerical model of swelling material coupling flow through porous media with ionic transport. More...

#include <swelling.h>

Inheritance diagram for TPZSwelling:
[legend]
Collaboration diagram for TPZSwelling:
[legend]

Public Member Functions

 TPZSwelling (int matindex, STATE lambda, STATE shear, STATE alfa, STATE M, STATE Gamma, STATE Kperm, STATE DPlus, STATE DMinus, STATE rHinder, STATE Cfc, STATE Nf0, STATE NPlus0, STATE NMinus0)
 Constructor of the class, where the user needs to specify the most important parameters. More...
 
virtual ~TPZSwelling ()
 
virtual int Dimension () const override
 Dimension of the problem. More...
 
virtual int NStateVariables () const override
 Number of state variables, in this case:
3 displacements, 1 pressure, 3 eletrochemical potencials, 1 eletrical potencial. More...
 
virtual void Print (std::ostream &out) override
 Prints out the data associated with the material. More...
 
virtual std::string Name () override
 Returns the name of the material. More...
 
void SetComputationMode (int mode)
 
virtual int VariableIndex (const std::string &name) override
 
virtual int NSolutionVariables (int var) override
 Returns the number of solution variables associated with a variable index. More...
 
virtual void Solution (TPZMaterialData &data, int var, TPZVec< STATE > &Solout) override
 Returns the solution associated with the var index based on the finite element approximation. More...
 
virtual int ClassId () const override
 Define the class id associated with the class. More...
 
Contribute methods (weak formulation)
virtual void Contribute (TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override
 It computes a contribution to the stiffness matrix and load vector at one integration point. More...
 
virtual void Contribute (TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef) override
 It computes a contribution to the residual vector at one integration point. More...
 
virtual void ContributeBC (TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
 It computes a contribution to the stiffness matrix and load vector at one BC integration point. More...
 
virtual void ContributeBC (TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
 It computes a contribution to the stiffness matrix and load vector at one BC integration point. More...
 
- Public Member Functions inherited from TPZMaterial
 TPZMaterial (int id)
 Creates a material object and inserts it in the vector of material pointers of the mesh. More...
 
 TPZMaterial ()
 Default constructor. More...
 
 TPZMaterial (const TPZMaterial &mat)
 Creates a material object based on the referred object and inserts it in the vector of material pointers of the mesh. More...
 
TPZMaterialoperator= (const TPZMaterial &copy)
 operator = More...
 
virtual ~TPZMaterial ()
 Default destructor. More...
 
virtual void FillDataRequirements (TPZMaterialData &data)
 Fill material data parameter with necessary requirements for the. More...
 
virtual void FillDataRequirements (TPZVec< TPZMaterialData > &datavec)
 Fill material data parameter with necessary requirements for the Contribute method. Here, in base class, all requirements are considered as necessary. Each derived class may optimize performance by selecting only the necessary data. More...
 
virtual void FillBoundaryConditionDataRequirement (int type, TPZMaterialData &data)
 This method defines which parameters need to be initialized in order to compute the contribution of the boundary condition. More...
 
virtual void FillBoundaryConditionDataRequirement (int type, TPZVec< TPZMaterialData > &datavec)
 This method defines which parameters need to be initialized in order to compute the contribution of the boundary condition. More...
 
virtual void FillDataRequirementsInterface (TPZMaterialData &data)
 This method defines which parameters need to be initialized in order to compute the contribution of interface elements. More...
 
virtual void FillDataRequirementsInterface (TPZMaterialData &data, TPZVec< TPZMaterialData > &datavec_left, TPZVec< TPZMaterialData > &datavec_right)
 This method defines which parameters need to be initialized in order to compute the contribution of interface elements. More...
 
int Id () const
 
void SetId (int id)
 
virtual int NFluxes ()
 Returns the number of components which form the flux function. More...
 
int NumLoadCases ()
 returns the number of load cases for this material object More...
 
virtual int MinimumNumberofLoadCases ()
 returns the minimum number of load cases for this material More...
 
void SetNumLoadCases (int numloadcases)
 changes the number of load cases for this material More...
 
void SetPostProcessIndex (int index)
 indicates which variable should be post processed More...
 
virtual TPZBndCondCreateBC (TPZMaterial *reference, int id, int typ, TPZFMatrix< STATE > &val1, TPZFMatrix< STATE > &val2)
 Creates an object TPZBndCond derived of TPZMaterial. More...
 
void SetForcingFunction (TPZAutoPointer< TPZFunction< STATE > > fp)
 Sets a procedure as source function for the material. More...
 
void SetForcingFunction (void(*fp)(const TPZVec< REAL > &loc, TPZVec< STATE > &result), int porder)
 
void SetForcingFunction (void(*fp)(const TPZVec< REAL > &loc, TPZVec< STATE > &result, TPZFMatrix< STATE > &gradu), int porder)
 
TPZAutoPointer< TPZFunction< STATE > > & ForcingFunction ()
 Returns a procedure as source function for the material. More...
 
void SetForcingFunctionExact (TPZAutoPointer< TPZFunction< STATE > > fp)
 Sets a procedure as exact solution for the problem. More...
 
TPZAutoPointer< TPZFunction< STATE > > & ForcingFunctionExact ()
 Returns a procedure as exact solution for the problem. More...
 
void SetTimeDependentForcingFunction (TPZAutoPointer< TPZFunction< STATE > > fp)
 Sets a procedure as time dependent source function for the material. More...
 
TPZAutoPointer< TPZFunction< STATE > > & TimeDependentForcingFunction ()
 Returns a procedure as time dependent source function for the material. More...
 
void SetTimeDependentFunctionExact (TPZAutoPointer< TPZFunction< STATE > > fp)
 Sets a procedure as time dependent exact solution for the problem. More...
 
TPZAutoPointer< TPZFunction< STATE > > & TimedependentFunctionExact ()
 Returns a procedure as time dependent exact solution for the problem. More...
 
void SetBCForcingFunction (TPZAutoPointer< TPZFunction< STATE > > fp)
 Sets a procedure as variable boundary condition. More...
 
TPZAutoPointer< TPZFunction< STATE > > & BCForcingFunction ()
 Returns a procedure as variable boundary condition. More...
 
void SetTimedependentBCForcingFunction (TPZAutoPointer< TPZFunction< STATE > > fp)
 Sets a procedure as time variable boundary condition. More...
 
TPZAutoPointer< TPZFunction< STATE > > & TimedependentBCForcingFunction ()
 Returns a procedure as time variable boundary condition. More...
 
virtual int HasForcingFunction ()
 Directive that gives true if the material has a forcing function. More...
 
virtual int HasForcingFunctionExact ()
 Directive that gives true if the material has a function exact. More...
 
virtual int HasBCForcingFunction ()
 Directive that gives true if the material has a bc forcing function exact. More...
 
virtual int HasTimedependentFunctionExact ()
 Directive that gives true if the material has a time dependent function exact. More...
 
virtual int HasTimedependentForcingFunction ()
 Directive that gives true if the material has a time dependent forcing function. More...
 
virtual int HasTimedependentBCForcingFunction ()
 Directive that gives true if the material has a time dependent bc forcing function. More...
 
virtual int IntegrationRuleOrder (int elPMaxOrder) const
 Gets the order of the integration rule necessary to integrate an element with polinomial order p. More...
 
virtual int IntegrationRuleOrder (TPZVec< int > &elPMaxOrder) const
 Gets the order of the integration rule necessary to integrate an element multiphysic. More...
 
virtual void Errors (TPZMaterialData &data, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &errors)
 
virtual void Errors (TPZVec< TPZMaterialData > &data, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &errors)
 
virtual void Errors (TPZVec< REAL > &x, TPZVec< STATE > &sol, TPZFMatrix< STATE > &dsol, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux, TPZVec< STATE > &uexact, TPZFMatrix< STATE > &duexact, TPZVec< REAL > &val)
 Computes the error due to the difference between the interpolated flux
and the flux computed based on the derivative of the solution. More...
 
virtual void ErrorsHdiv (TPZMaterialData &data, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &values)
 
virtual int NEvalErrors ()
 Returns the number of norm errors. Default is 3: energy, L2 and H1. More...
 
virtual TPZMaterialNewMaterial ()
 To create another material of the same type. More...
 
virtual void SetData (std::istream &data)
 Reads data of the material from a istream (file data) More...
 
virtual void Clone (std::map< int, TPZMaterial * > &matvec)
 Creates a copy of the material object and put it in the vector which is passed on. More...
 
virtual int FluxType ()
 To return a numerical flux type to apply over the interfaces of the elements. More...
 
virtual void ContributeErrors (TPZMaterialData &data, REAL weight, TPZVec< REAL > &nk, int &errorid)
 
virtual REAL ComputeSquareResidual (TPZVec< REAL > &X, TPZVec< STATE > &sol, TPZFMatrix< STATE > &dsol)
 Computes square of residual of the differential equation at one integration point. More...
 
virtual int PushMemItem (int sourceIndex=-1)
 Pushes a new entry in the context of materials with memory, returning its index at the internal storage stack. More...
 
virtual void FreeMemItem (int index)
 Frees an entry in the material with memory internal history storage. More...
 
void SetLinearContext (bool IsLinear)
 Sets fLinearContext attribute. More...
 
bool GetLinearContext () const
 Returns fLinearContext attribute. More...
 
virtual void Contribute (TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)
 It computes a contribution to the stiffness matrix and load vector at one integration point to multiphysics simulation. More...
 
virtual void Contribute (TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ef)
 It computes a contribution to the stiffness matrix and load vector at one integration point to multiphysics simulation. More...
 
virtual void ContributeBC (TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
 It computes a contribution to the stiffness matrix and load vector at one BC integration point to multiphysics simulation. More...
 
virtual void ContributeBC (TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
 It computes a contribution to the stiffness matrix and load vector at one BC integration point to multiphysics simulation. More...
 
int ClassId () const override
 Unique identifier for serialization purposes. More...
 
void Write (TPZStream &buf, int withclassid) const override
 Saves the element data to a stream. More...
 
void Read (TPZStream &buf, void *context) override
 Reads the element data from a stream. More...
 
virtual void Solution (TPZVec< TPZMaterialData > &datavec, int var, TPZVec< STATE > &Solout)
 Returns the solution associated with the var index based on the finite element approximation. More...
 
virtual void Solution (TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleftvec, TPZVec< TPZMaterialData > &datarightvec, int var, TPZVec< STATE > &Solout)
 Returns the solution associated with the var index based on the finite element approximation around one interface element. More...
 
virtual void Solution (TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleftvec, TPZVec< TPZMaterialData > &datarightvec, int var, TPZVec< STATE > &Solout, TPZCompEl *left, TPZCompEl *ritgh)
 Returns the solution associated with the var index based on the finite element approximation around one interface element. More...
 
virtual void Flux (TPZVec< REAL > &x, TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux)
 Computes the value of the flux function to be used by ZZ error estimator. More...
 
- Public Member Functions inherited from TPZSavable
 TPZSavable ()
 
virtual ~TPZSavable ()
 
virtual std::list< std::map< std::string, uint64_t > > VersionHistory () const
 
virtual std::pair< std::string, uint64_t > Version () const
 
virtual bool Compare (TPZSavable *copy, bool override=false)
 Compares the object for identity with the object pointed to, eventually copy the object. More...
 
virtual bool Compare (TPZSavable *copy, bool override=false) const
 Compares the object for identity with the object pointed to, eventually copy the object. More...
 
- Public Member Functions inherited from TPZRegisterClassId
template<typename T >
 TPZRegisterClassId (int(T::*)() const)
 
 TPZRegisterClassId ()=default
 

Protected Member Functions

virtual void Solution (TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, int var, TPZVec< STATE > &Solout) override
 Computes a post-processed solution variable corresponding to the variable index. More...
 
- Protected Member Functions inherited from TPZMaterial

Private Member Functions

void ExactSolution (TPZVec< STATE > &mu, STATE ksi, STATE pres, TPZVec< STATE > &N)
 
void ComputeN (TPZVec< STATE > &N, TPZVec< STATE > &mu, STATE ksi)
 Computes the value of the N coeficients in function of ksi and mus, iterative method, inverting the Hessian of W. More...
 
void ComputeInitialGuess (TPZVec< STATE > &mu, STATE J, STATE &pres, STATE &ksi, TPZVec< STATE > &N)
 Computes the aproximate values of the pressure, ksi and N based on mu and J by direct inversion of the formulas. More...
 
void ComputeN (TPZVec< STATE > &mu, STATE ksi, STATE pressure, TPZVec< STATE > &N)
 Computes the value of the N coeficients in function of ksi and mus, iterative method, inverting the Hessian of W. More...
 
void NResidual (TPZVec< STATE > &mu, STATE ksi, STATE pressure, TPZVec< STATE > &N, TPZFMatrix< STATE > &res, TPZFMatrix< STATE > &tangent)
 Computes the residual and tangent vector of the system of equations which determines N. More...
 
int NumCases ()
 Methods needed to perform convergence checks. More...
 
void LoadState (TPZFMatrix< STATE > &state)
 Loads the state within the current object, to be used when computing the tangent matrix. More...
 
void ComputeTangent (TPZFMatrix< STATE > &tangent, TPZVec< REAL > &coefs, int cases)
 Computes the tangent matrix for a given loadcase. More...
 
void Residual (TPZFMatrix< STATE > &res, int cases)
 Computes the residual for the given state variable. More...
 

Private Attributes

int fComputationMode
 Computation mode:
0->residual wrt to time $ n $
1->residual and tangent wrt to time $ n+1 $. More...
 
TPZFMatrix< STATE > fKperm
 Hydraulic permeability $ [mm^4 / (Ns)] $. More...
 
STATE fLambda
 Compression modulus [N/mm^2]. More...
 
STATE fGamma
 Osmotic coeficient (no dimension) More...
 
STATE fShear
 Shear modulus [N/mm^2]. More...
 
STATE fAlfa
 Biot coupling coeficient (no dimension) More...
 
STATE fM
 Storage modulus [N/mm^2]. More...
 
STATE fDPlus
 Diffusion coeficient for cations [mm^2/s]. More...
 
STATE fDMinus
 Diffusion coeficient for anions [mm^2/s]. More...
 
STATE frHinder
 Hindrance factor (no dimension) More...
 
STATE fInitDeform
 Initial deformation (assuming everything occurs isotropically, a constant is suficient (no dimension) More...
 
STATE fCfc
 Fixed charge density [mmoleq/mm^3]. More...
 
STATE fNf0
 Initial fluid volume fraction (do dimension) More...
 
STATE fNPlus0
 Initial cation volume fraction (no dimension) More...
 
STATE fNMinus0
 Initial anion volume fraction (no dimension) More...
 
STATE fTheta
 Timestepping parameter theta (no dimension) More...
 
STATE fDelt
 Timestep [s]. More...
 

Static Private Attributes

static STATE gExtConc
 External concentration (used as reference value for pressure) [mmol/mm^3]. More...
 
static STATE gFaraday
 Faraday constant [C/mmol]. More...
 
static STATE gVPlus
 Molar volume cation [mm^3/mmol]. More...
 
static STATE gVMinus
 Molar volume anions [mm^3/mmol]. More...
 
static STATE gRGas
 gas constant [Nmm/(mmol K)] More...
 
static STATE gTemp
 Absolute temperature [K]. More...
 
static STATE gMuRef [3]
 Reference chemical potentials (order f,plus,minus) [mV]. More...
 
static TPZFMatrix< STATE > gState
 Variables which holds the state variables used in the check convergence procedure. More...
 
static TPZFMatrix< REAL > gphi
 Variables which holds the state variables used in the check convergence procedure. More...
 
static TPZFMatrix< REAL > gdphi
 

Additional Inherited Members

- Static Public Member Functions inherited from TPZSavable
static std::set< TPZRestoreClassBase * > & RestoreClassSet ()
 This static function guarantees that the gMap object is available when needed. More...
 
static std::map< int, TPZRestore_t > & ClassIdMap ()
 This static function guarantees that the gMap object is available when needed. More...
 
static std::pair< std::string, uint64_t > NeoPZVersion ()
 
static void Register (TPZRestoreClassBase *restore)
 
static void RegisterClassId (int classid, TPZRestore_t fun)
 
static TPZSavableCreateInstance (const int &classId)
 
- Static Public Attributes inherited from TPZMaterial
static REAL gBigNumber
 Big number to penalization method, used for Dirichlet conditions. More...
 
- Protected Attributes inherited from TPZMaterial
TPZAutoPointer< TPZFunction< STATE > > fForcingFunction
 Pointer to forcing function, it is the right member at differential equation. More...
 
TPZAutoPointer< TPZFunction< STATE > > fForcingFunctionExact
 Pointer to exact solution function, needed to calculate exact error. More...
 
TPZAutoPointer< TPZFunction< STATE > > fTimeDependentForcingFunction
 Pointer to time dependent forcing function, it is the right member at differential equation. More...
 
TPZAutoPointer< TPZFunction< STATE > > fTimedependentFunctionExact
 Pointer to time dependent exact solution function, needed to calculate exact error. More...
 
TPZAutoPointer< TPZFunction< STATE > > fBCForcingFunction
 Pointer to bc forcing function, it is a variable boundary condition at differential equation. More...
 
TPZAutoPointer< TPZFunction< STATE > > fTimedependentBCForcingFunction
 Pointer to time dependent bc forcing function, it is a variable boundary condition at differential equation. More...
 
bool fLinearContext
 Defines whether the equation context is linear solver or non linear. More...
 
int fNumLoadCases
 Defines the number of load cases generated by this material. More...
 
int fPostProcIndex
 indicates which solution should be used for post processing More...
 

Detailed Description

The TPZSwelling class implements a numerical model of swelling material coupling flow through porous media with ionic transport.

Definition at line 21 of file swelling.h.

Constructor & Destructor Documentation

◆ TPZSwelling()

TPZSwelling::TPZSwelling ( int  matindex,
STATE  lambda,
STATE  shear,
STATE  alfa,
STATE  M,
STATE  Gamma,
STATE  Kperm,
STATE  DPlus,
STATE  DMinus,
STATE  rHinder,
STATE  Cfc,
STATE  Nf0,
STATE  NPlus0,
STATE  NMinus0 
)

Constructor of the class, where the user needs to specify the most important parameters.

Parameters
matindexindex of material
lambdaCompression modulus
shearShear modulus
alfaBiot coupling coeficient
MStorage modulus
GammaOsmotic coeficient
Kpermhydraulic permeability (isotropic)
DPlusDiffusion coeficient for cations
DMinusDiffusion coeficient for anions
rHinderHindrance factor
CfcFixed charge density
Nf0Initial fluid volume fraction
NPlus0Initial cation volume fraction
NMinus0Initial anion volume fraction

Definition at line 45 of file swelling.cpp.

References fAlfa, fCfc, fComputationMode, fDelt, fDMinus, fDPlus, fGamma, fInitDeform, fKperm, fLambda, fM, fNf0, fNMinus0, fNPlus0, frHinder, fShear, fTheta, and TPZFMatrix< TVar >::Redim().

◆ ~TPZSwelling()

TPZSwelling::~TPZSwelling ( )
virtual

Definition at line 72 of file swelling.cpp.

Member Function Documentation

◆ ClassId()

int TPZSwelling::ClassId ( ) const
overridevirtual

Define the class id associated with the class.

This id has to be unique for all classes A non unique id is flagged at the startup of the program

Implements TPZSavable.

Definition at line 842 of file swelling.cpp.

References TPZMaterial::ClassId(), and Hash().

Referenced by Solution().

◆ ComputeInitialGuess()

void TPZSwelling::ComputeInitialGuess ( TPZVec< STATE > &  mu,
STATE  J,
STATE &  pres,
STATE &  ksi,
TPZVec< STATE > &  N 
)
private

Computes the aproximate values of the pressure, ksi and N based on mu and J by direct inversion of the formulas.

This method has been superseded by the direct computation ExactSolution

Definition at line 750 of file swelling.cpp.

References exp, fCfc, fGamma, fNf0, gFaraday, gMuRef, gRGas, gTemp, gVMinus, gVPlus, log, and sqrt.

Referenced by ContributeBC().

◆ ComputeN() [1/2]

void TPZSwelling::ComputeN ( TPZVec< STATE > &  N,
TPZVec< STATE > &  mu,
STATE  ksi 
)
private

Computes the value of the N coeficients in function of ksi and mus, iterative method, inverting the Hessian of W.

This method has been superseded by the direct computation ExactSolution

Referenced by ContributeBC(), and ExactSolution().

◆ ComputeN() [2/2]

void TPZSwelling::ComputeN ( TPZVec< STATE > &  mu,
STATE  ksi,
STATE  pressure,
TPZVec< STATE > &  N 
)
private

Computes the value of the N coeficients in function of ksi and mus, iterative method, inverting the Hessian of W.

This method has been superseded by the direct computation ExactSolution

◆ ComputeTangent()

void TPZSwelling::ComputeTangent ( TPZFMatrix< STATE > &  tangent,
TPZVec< REAL > &  coefs,
int  cases 
)
private

Computes the tangent matrix for a given loadcase.

Referenced by ContributeBC().

◆ Contribute() [1/2]

virtual void TPZSwelling::Contribute ( TPZMaterialData data,
REAL  weight,
TPZFMatrix< STATE > &  ek,
TPZFMatrix< STATE > &  ef 
)
inlineoverridevirtual

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 stiffness matrix
ef[out] is the load vector
Since
April 16, 2007

Implements TPZMaterial.

Definition at line 151 of file swelling.h.

◆ Contribute() [2/2]

virtual void TPZSwelling::Contribute ( TPZMaterialData data,
REAL  weight,
TPZFMatrix< STATE > &  ef 
)
inlineoverridevirtual

It computes a contribution to the residual vector at one integration point.

Parameters
data[in] stores all input data
weight[in] is the weight of the integration rule
ef[out] is the residual vector
Since
April 16, 2007

Reimplemented from TPZMaterial.

Definition at line 157 of file swelling.h.

References bc, and ContributeBC().

◆ ContributeBC() [1/2]

void TPZSwelling::ContributeBC ( TPZMaterialData data,
REAL  weight,
TPZFMatrix< STATE > &  ek,
TPZFMatrix< STATE > &  ef,
TPZBndCond bc 
)
overridevirtual

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 stiffness matrix
ef[out] is the load vector
bc[in] is the boundary condition material
Since
October 07, 2011

Implements TPZMaterial.

Definition at line 128 of file swelling.cpp.

References CheckConvergence(), ComputeInitialGuess(), ComputeN(), ComputeTangent(), DebugStop, ELU, ExactSolution(), fComputationMode, fInitDeform, TPZMaterial::gBigNumber, gdphi, gFaraday, gphi, gState, gVMinus, gVPlus, idf, LoadState(), main(), TPZBndCond::Material(), TPZVec< T >::NElements(), Norm(), NResidual(), NStateVariables(), NumCases(), TPZMaterialData::phi, Print(), TPZMatrix< TVar >::Print(), PZError, TPZFMatrix< TVar >::Redim(), test::res, Residual(), TPZMatrix< TVar >::Rows(), TPZVec< T >::size(), TPZMaterialData::sol, TPZMatrix< TVar >::SolveDirect(), sqrt, TPZBndCond::Type(), TPZBndCond::Val1(), TPZBndCond::Val2(), and rdt::values.

Referenced by Contribute().

◆ ContributeBC() [2/2]

virtual void TPZSwelling::ContributeBC ( TPZMaterialData data,
REAL  weight,
TPZFMatrix< STATE > &  ef,
TPZBndCond bc 
)
inlineoverridevirtual

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
ef[out] is the load vector
bc[in] is the boundary condition material
Since
April 16, 2007

Reimplemented from TPZMaterial.

Definition at line 168 of file swelling.h.

References ComputeInitialGuess(), ComputeN(), ComputeTangent(), TPZMaterial::ContributeBC(), ExactSolution(), LoadState(), main(), NResidual(), NumCases(), test::res, and Residual().

◆ Dimension()

virtual int TPZSwelling::Dimension ( ) const
inlineoverridevirtual

Dimension of the problem.

Implements TPZMaterial.

Definition at line 121 of file swelling.h.

◆ ExactSolution()

void TPZSwelling::ExactSolution ( TPZVec< STATE > &  mu,
STATE  ksi,
STATE  pres,
TPZVec< STATE > &  N 
)
private

Definition at line 779 of file swelling.cpp.

References ComputeN(), exp, fGamma, gFaraday, gRGas, gTemp, gVMinus, gVPlus, pow(), and val().

Referenced by ContributeBC().

◆ LoadState()

void TPZSwelling::LoadState ( TPZFMatrix< STATE > &  state)
private

Loads the state within the current object, to be used when computing the tangent matrix.

Referenced by ContributeBC().

◆ Name()

virtual std::string TPZSwelling::Name ( )
inlineoverridevirtual

Returns the name of the material.

Reimplemented from TPZMaterial.

Definition at line 131 of file swelling.h.

Referenced by Print().

◆ NResidual()

void TPZSwelling::NResidual ( TPZVec< STATE > &  mu,
STATE  ksi,
STATE  pressure,
TPZVec< STATE > &  N,
TPZFMatrix< STATE > &  res,
TPZFMatrix< STATE > &  tangent 
)
private

Computes the residual and tangent vector of the system of equations which determines N.

This method has been superseded by the direct computation ExactSolution

Referenced by ContributeBC().

◆ NSolutionVariables()

int TPZSwelling::NSolutionVariables ( int  var)
overridevirtual

Returns the number of solution variables associated with a variable index.

(e.g. 1 for scalar, 3 for vectorial, 9 for tensorial)

Reimplemented from TPZMaterial.

Definition at line 115 of file swelling.cpp.

References TPZMaterial::NSolutionVariables().

◆ NStateVariables()

virtual int TPZSwelling::NStateVariables ( ) const
inlineoverridevirtual

Number of state variables, in this case:
3 displacements, 1 pressure, 3 eletrochemical potencials, 1 eletrical potencial.

Implements TPZMaterial.

Definition at line 127 of file swelling.h.

References Print().

Referenced by ContributeBC().

◆ NumCases()

int TPZSwelling::NumCases ( )
private

Methods needed to perform convergence checks.

Number of cases which are considered for convergence checks

Referenced by ContributeBC().

◆ Print()

void TPZSwelling::Print ( std::ostream &  out)
overridevirtual

Prints out the data associated with the material.

Reimplemented from TPZMaterial.

Definition at line 76 of file swelling.cpp.

References fAlfa, fCfc, fComputationMode, fDelt, fDMinus, fDPlus, fGamma, fInitDeform, fKperm, fLambda, fM, fNf0, fNMinus0, fNPlus0, frHinder, fShear, fTheta, gFaraday, gMuRef, gRGas, gTemp, gVMinus, gVPlus, Name(), and TPZMaterial::Print().

Referenced by ContributeBC(), and NStateVariables().

◆ Residual()

void TPZSwelling::Residual ( TPZFMatrix< STATE > &  res,
int  cases 
)
private

Computes the residual for the given state variable.

Referenced by ContributeBC().

◆ SetComputationMode()

void TPZSwelling::SetComputationMode ( int  mode)
inline

Definition at line 133 of file swelling.h.

◆ Solution() [1/2]

void TPZSwelling::Solution ( TPZVec< STATE > &  Sol,
TPZFMatrix< STATE > &  DSol,
TPZFMatrix< REAL > &  axes,
int  var,
TPZVec< STATE > &  Solout 
)
overrideprotectedvirtual

Computes a post-processed solution variable corresponding to the variable index.

Reimplemented from TPZMaterial.

Definition at line 120 of file swelling.cpp.

References TPZMaterial::Solution().

Referenced by Solution().

◆ Solution() [2/2]

virtual void TPZSwelling::Solution ( TPZMaterialData data,
int  var,
TPZVec< STATE > &  Solout 
)
inlineoverridevirtual

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

Reimplemented from TPZMaterial.

Definition at line 294 of file swelling.h.

References TPZMaterialData::axes, ClassId(), DebugStop, TPZMaterialData::dsol, TPZVec< T >::size(), TPZMaterialData::sol, and Solution().

◆ VariableIndex()

int TPZSwelling::VariableIndex ( const std::string &  name)
overridevirtual

returns the variable index associated with the name

Reimplemented from TPZMaterial.

Definition at line 111 of file swelling.cpp.

References TPZMaterial::VariableIndex().

Member Data Documentation

◆ fAlfa

STATE TPZSwelling::fAlfa
private

Biot coupling coeficient (no dimension)

Definition at line 38 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fCfc

STATE TPZSwelling::fCfc
private

Fixed charge density [mmoleq/mm^3].

Definition at line 50 of file swelling.h.

Referenced by ComputeInitialGuess(), Print(), and TPZSwelling().

◆ fComputationMode

int TPZSwelling::fComputationMode
private

Computation mode:
0->residual wrt to time $ n $
1->residual and tangent wrt to time $ n+1 $.

Definition at line 28 of file swelling.h.

Referenced by ContributeBC(), Print(), and TPZSwelling().

◆ fDelt

STATE TPZSwelling::fDelt
private

Timestep [s].

Definition at line 61 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fDMinus

STATE TPZSwelling::fDMinus
private

Diffusion coeficient for anions [mm^2/s].

Definition at line 44 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fDPlus

STATE TPZSwelling::fDPlus
private

Diffusion coeficient for cations [mm^2/s].

Definition at line 42 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fGamma

STATE TPZSwelling::fGamma
private

Osmotic coeficient (no dimension)

Definition at line 34 of file swelling.h.

Referenced by ComputeInitialGuess(), ExactSolution(), Print(), and TPZSwelling().

◆ fInitDeform

STATE TPZSwelling::fInitDeform
private

Initial deformation (assuming everything occurs isotropically, a constant is suficient (no dimension)

Definition at line 48 of file swelling.h.

Referenced by ContributeBC(), Print(), and TPZSwelling().

◆ fKperm

TPZFMatrix<STATE> TPZSwelling::fKperm
private

Hydraulic permeability $ [mm^4 / (Ns)] $.

Definition at line 30 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fLambda

STATE TPZSwelling::fLambda
private

Compression modulus [N/mm^2].

Definition at line 32 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fM

STATE TPZSwelling::fM
private

Storage modulus [N/mm^2].

Definition at line 40 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fNf0

STATE TPZSwelling::fNf0
private

Initial fluid volume fraction (do dimension)

Definition at line 52 of file swelling.h.

Referenced by ComputeInitialGuess(), Print(), and TPZSwelling().

◆ fNMinus0

STATE TPZSwelling::fNMinus0
private

Initial anion volume fraction (no dimension)

Definition at line 56 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fNPlus0

STATE TPZSwelling::fNPlus0
private

Initial cation volume fraction (no dimension)

Definition at line 54 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ frHinder

STATE TPZSwelling::frHinder
private

Hindrance factor (no dimension)

Definition at line 46 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fShear

STATE TPZSwelling::fShear
private

Shear modulus [N/mm^2].

Definition at line 36 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ fTheta

STATE TPZSwelling::fTheta
private

Timestepping parameter theta (no dimension)

Definition at line 59 of file swelling.h.

Referenced by Print(), and TPZSwelling().

◆ gdphi

TPZFMatrix< REAL > TPZSwelling::gdphi
staticprivate

Definition at line 281 of file swelling.h.

Referenced by ContributeBC().

◆ gExtConc

STATE TPZSwelling::gExtConc
staticprivate

External concentration (used as reference value for pressure) [mmol/mm^3].

Definition at line 63 of file swelling.h.

◆ gFaraday

STATE TPZSwelling::gFaraday
staticprivate

Faraday constant [C/mmol].

Definition at line 81 of file swelling.h.

Referenced by ComputeInitialGuess(), ContributeBC(), ExactSolution(), and Print().

◆ gMuRef

STATE TPZSwelling::gMuRef
staticprivate

Reference chemical potentials (order f,plus,minus) [mV].

Definition at line 91 of file swelling.h.

Referenced by ComputeInitialGuess(), and Print().

◆ gphi

TPZFMatrix< REAL > TPZSwelling::gphi
staticprivate

Variables which holds the state variables used in the check convergence procedure.

Definition at line 281 of file swelling.h.

Referenced by ContributeBC().

◆ gRGas

STATE TPZSwelling::gRGas
staticprivate

gas constant [Nmm/(mmol K)]

Definition at line 87 of file swelling.h.

Referenced by ComputeInitialGuess(), ExactSolution(), and Print().

◆ gState

TPZFMatrix< STATE > TPZSwelling::gState
staticprivate

Variables which holds the state variables used in the check convergence procedure.

Definition at line 279 of file swelling.h.

Referenced by ContributeBC().

◆ gTemp

STATE TPZSwelling::gTemp
staticprivate

Absolute temperature [K].

Definition at line 89 of file swelling.h.

Referenced by ComputeInitialGuess(), ExactSolution(), and Print().

◆ gVMinus

STATE TPZSwelling::gVMinus
staticprivate

Molar volume anions [mm^3/mmol].

Definition at line 85 of file swelling.h.

Referenced by ComputeInitialGuess(), ContributeBC(), ExactSolution(), and Print().

◆ gVPlus

STATE TPZSwelling::gVPlus
staticprivate

Molar volume cation [mm^3/mmol].

Definition at line 83 of file swelling.h.

Referenced by ComputeInitialGuess(), ContributeBC(), ExactSolution(), and Print().


The documentation for this class was generated from the following files: