NeoPZ
|
Material to solve a mixed poisson problem 3D by multiphysics simulation. More...
#include <pzmatmixedpoisson3d.h>
Public Member Functions | |
TPZMatMixedPoisson3D () | |
TPZMatMixedPoisson3D (int matid, int dim) | |
virtual | ~TPZMatMixedPoisson3D () |
TPZMatMixedPoisson3D (const TPZMatMixedPoisson3D ©) | |
copy constructor More... | |
TPZMatMixedPoisson3D & | operator= (const TPZMatMixedPoisson3D ©) |
virtual std::string | Name () override |
Returns the name of the material. More... | |
int | Dimension () const override |
Returns the integrable dimension of the material. More... | |
void | SetDimension (int dim) |
int | MatId () |
virtual TPZMaterial * | NewMaterial () override |
To create another material of the same type. More... | |
virtual int | NStateVariables () const override |
Returns the number of state variables associated with the material. More... | |
void | SetPermeability (REAL perm) |
void | SetPermeabilityTensor (TPZFMatrix< REAL > K, TPZFMatrix< REAL > invK) |
void | SetReactionTerm (REAL alpha) |
void | SetViscosity (REAL visc) |
void | SetInternalFlux (REAL flux) |
void | SetPermeabilityFunction (TPZAutoPointer< TPZFunction< STATE > > fp) |
TPZAutoPointer< TPZFunction< STATE > > | PermeabilityFunction () |
void | SetfReactionTermFunction (TPZAutoPointer< TPZFunction< STATE > > fp) |
TPZAutoPointer< TPZFunction< STATE > > | ReactionTermFunction () |
virtual int | IntegrationRuleOrder (int elPMaxOrder) const override |
Gets the order of the integration rule necessary to integrate an element with polinomial order p. More... | |
virtual int | IntegrationRuleOrder (TPZVec< int > &elPMaxOrder) const override |
Gets the order of the integration rule necessary to integrate an element multiphysic. More... | |
void | Print (std::ostream &out) override |
Prints out the data associated with the material. More... | |
Contribute methods | |
virtual void | Contribute (TPZVec< TPZMaterialData > &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. This method use normalized piola contravariant mapping for nonlinear mappings. With second integration by parts. More... | |
void | ContributeWithoutSecondIntegration (TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) |
This method use piola contravariant mapping for nonlinear mappings. More... | |
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 | ContributeBC (TPZVec< TPZMaterialData > &datavec, 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 > &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 | FillDataRequirements (TPZVec< TPZMaterialData > &datavec) override |
It computes a contribution to stiffness matrix and load vector at one BC integration point. More... | |
virtual void | FillBoundaryConditionDataRequirement (int type, TPZVec< TPZMaterialData > &datavec) override |
This method defines which parameters need to be initialized in order to compute the contribution of the boundary condition. More... | |
int | VariableIndex (const std::string &name) override |
int | NSolutionVariables (int var) override |
Returns the number of variables associated with the variable indexed by var. More... | |
void | Solution (TPZVec< TPZMaterialData > &datavec, int var, TPZVec< STATE > &Solout) override |
Returns the solution associated with the var index based on the finite element approximation. More... | |
void | Solution (TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, int var, TPZVec< STATE > &Solout) override |
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... | |
void | Errors (TPZVec< REAL > &x, TPZVec< STATE > &u, TPZFMatrix< STATE > &dudx, TPZFMatrix< REAL > &axes, TPZVec< STATE > &, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &values) override |
Computes the error due to the difference between the interpolated flux and the flux computed based on the derivative of the solution. More... | |
void | ErrorsHdiv (TPZMaterialData &data, TPZVec< STATE > &u_exact, TPZFMatrix< STATE > &du_exact, TPZVec< REAL > &values) override |
void | UseSecondIntegrationByParts () |
bool | IsUsedSecondIntegration () |
virtual int | ClassId () const override |
Define the class id associated with the class. 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... | |
TPZMaterial & | operator= (const TPZMaterial ©) |
operator = More... | |
virtual | ~TPZMaterial () |
Default destructor. More... | |
virtual void | FillDataRequirements (TPZMaterialData &data) |
Fill material data parameter with necessary requirements for the. 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 | 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 TPZBndCond * | CreateBC (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 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 int | NEvalErrors () |
Returns the number of norm errors. Default is 3: energy, L2 and H1. 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 (TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef) |
It computes a contribution to the residual vector at one integration point. 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 (TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc) |
It computes a contribution to the stiffness matrix and load vector at one BC integration point. 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 (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 Attributes | |
int | fMatId |
REAL | fF |
int | fDim |
REAL | falpha |
REAL | fvisc |
fluid viscosity More... | |
TPZFMatrix< REAL > | fTensorK |
permeability tensor. Coeficient which multiplies the gradient operator More... | |
TPZFMatrix< REAL > | fInvK |
inverse of the permeability tensor. More... | |
TPZAutoPointer< TPZFunction< STATE > > | fPermeabilityFunction |
Pointer to forcing function, it is the Permeability and its inverse. More... | |
TPZAutoPointer< TPZFunction< STATE > > | fReactionTermFunction |
Pointer to forcing function, it is the reaction term. More... | |
TPZLagrangeMultiplier * | fmatLagr |
bool | fSecondIntegration |
Parameter to choose the second integration by parts in the variational formulation. More... | |
bool | fReactionTerm |
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... | |
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 TPZSavable * | CreateInstance (const int &classId) |
Static Public Attributes inherited from TPZMaterial | |
static REAL | gBigNumber |
Big number to penalization method, used for Dirichlet conditions. More... | |
Protected Member Functions inherited from TPZMaterial |
Material to solve a mixed poisson problem 3D by multiphysics simulation.
Definition at line 33 of file pzmatmixedpoisson3d.h.
TPZMatMixedPoisson3D::TPZMatMixedPoisson3D | ( | ) |
Valor da funcao de carga
Dimensao do dominio
Material id not initialized
Definition at line 26 of file pzmatmixedpoisson3d.cpp.
References falpha, fDim, fF, fInvK, fMatId, fmatLagr, fPermeabilityFunction, fReactionTerm, fReactionTermFunction, fSecondIntegration, fTensorK, fvisc, TPZMatrix< TVar >::Identity(), and TPZFMatrix< TVar >::Resize().
Referenced by NewMaterial().
TPZMatMixedPoisson3D::TPZMatMixedPoisson3D | ( | int | matid, |
int | dim | ||
) |
Valor da funcao de carga
Dimensao do dominio
Material id no initialized
Definition at line 53 of file pzmatmixedpoisson3d.cpp.
References falpha, fDim, fF, fInvK, fMatId, fmatLagr, fPermeabilityFunction, fReactionTerm, fReactionTermFunction, fSecondIntegration, fTensorK, fvisc, TPZMatrix< TVar >::Identity(), TPZFMatrix< TVar >::Redim(), and TPZFMatrix< TVar >::Resize().
|
virtual |
Definition at line 80 of file pzmatmixedpoisson3d.cpp.
TPZMatMixedPoisson3D::TPZMatMixedPoisson3D | ( | const TPZMatMixedPoisson3D & | copy | ) |
copy constructor
Definition at line 83 of file pzmatmixedpoisson3d.cpp.
References falpha, fDim, fF, fInvK, fMatId, fmatLagr, fReactionTerm, fSecondIntegration, fTensorK, and operator=().
|
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 805 of file pzmatmixedpoisson3d.cpp.
References TPZMaterial::ClassId(), and Hash().
Referenced by IsUsedSecondIntegration().
|
overridevirtual |
It computes a contribution to the stiffness matrix and load vector at one integration point. This method use normalized piola contravariant mapping for nonlinear mappings. With second integration by parts.
datavec | [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 |
Reimplemented from TPZMaterial.
Definition at line 132 of file pzmatmixedpoisson3d.cpp.
References TPZAxesTools< TVar >::Axes2XYZ(), TPZMatrix< TVar >::Cols(), TPZLagrangeMultiplier::Contribute(), ContributeWithoutSecondIntegration(), DebugStop, TPZFunction< TVar >::Execute(), falpha, fDim, fF, TPZMaterial::fForcingFunction, fInvK, fmatLagr, fPermeabilityFunction, fReactionTerm, fReactionTermFunction, fSecondIntegration, fTensorK, fvisc, TPZVec< T >::NElements(), test::res, TPZMatrix< TVar >::Rows(), TPZLagrangeMultiplier::SetMultiplier(), TPZVec< T >::size(), and TPZFMatrix< TVar >::Zero().
Referenced by IntegrationRuleOrder().
|
inlineoverridevirtual |
It computes a contribution to the stiffness matrix and load vector at one integration point.
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 |
Implements TPZMaterial.
Definition at line 209 of file pzmatmixedpoisson3d.h.
References bc, ContributeBC(), and DebugStop.
|
overridevirtual |
It computes a contribution to the stiffness matrix and load vector at one BC integration point.
datavec | [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 |
Reimplemented from TPZMaterial.
Definition at line 459 of file pzmatmixedpoisson3d.cpp.
References DebugStop, TPZFunction< TVar >::Execute(), TPZMaterial::ForcingFunction(), TPZMaterial::gBigNumber, TPZMaterial::HasForcingFunction(), test::res, TPZMatrix< TVar >::Rows(), TPZVec< T >::size(), TPZBndCond::Type(), TPZBndCond::Val1(), and TPZBndCond::Val2().
Referenced by Contribute().
|
inlineoverridevirtual |
It computes a contribution to the stiffness matrix and load vector at one BC integration point.
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 |
Implements TPZMaterial.
Definition at line 233 of file pzmatmixedpoisson3d.h.
References DebugStop, Errors(), ErrorsHdiv(), FillBoundaryConditionDataRequirement(), FillDataRequirements(), NSolutionVariables(), Solution(), rdt::values, and VariableIndex().
void TPZMatMixedPoisson3D::ContributeWithoutSecondIntegration | ( | TPZVec< TPZMaterialData > & | datavec, |
REAL | weight, | ||
TPZFMatrix< STATE > & | ek, | ||
TPZFMatrix< STATE > & | ef | ||
) |
This method use piola contravariant mapping for nonlinear mappings.
Definition at line 316 of file pzmatmixedpoisson3d.cpp.
References DebugStop, TPZFunction< TVar >::Execute(), falpha, fDim, fF, TPZMaterial::fForcingFunction, fInvK, fPermeabilityFunction, fReactionTerm, fReactionTermFunction, fTensorK, fvisc, TPZVec< T >::NElements(), test::res, TPZMatrix< TVar >::Rows(), TPZVec< T >::size(), and TPZFMatrix< TVar >::Zero().
Referenced by Contribute(), and IntegrationRuleOrder().
|
inlineoverridevirtual |
Returns the integrable dimension of the material.
Implements TPZMaterial.
Definition at line 86 of file pzmatmixedpoisson3d.h.
References fDim.
Referenced by Solution().
|
overridevirtual |
Computes the error due to the difference between the interpolated flux
and the flux computed based on the derivative of the solution.
:: need to be corrected
Reimplemented from TPZMaterial.
Definition at line 780 of file pzmatmixedpoisson3d.cpp.
References DebugStop, fabs, fDim, TPZVec< T >::Fill(), TPZMaterial::NEvalErrors(), TPZVec< T >::Resize(), and Solution().
Referenced by ContributeBC().
|
overridevirtual |
Reimplemented from TPZMaterial.
Definition at line 763 of file pzmatmixedpoisson3d.cpp.
References abs(), TPZVec< T >::Fill(), and Solution().
Referenced by ContributeBC().
|
overridevirtual |
This method defines which parameters need to be initialized in order to compute the contribution of the boundary condition.
Reimplemented from TPZMaterial.
Definition at line 542 of file pzmatmixedpoisson3d.cpp.
References TPZVec< T >::size().
Referenced by ContributeBC().
|
overridevirtual |
It computes a contribution to stiffness matrix and load vector at one BC integration point.
datavec | [in] |
dataleft | [in] |
weight | [in] |
ek | [out] is the stiffness matrix |
ef | [out] is the load vector |
bc | [in] is the boundary condition object |
data | [in] |
dataleft | [in] |
weight | [in] |
ek | [out] is the stiffness matrix |
ef | [out] is the load vector |
bc | [in] is the boundary condition object |
datavec | [in] |
dataleft | [in] |
dataright | [in] |
weight | [in] |
ek | [out] is the stiffness matrix |
ef | [out] is the load vector |
data | [in] |
dataleft | [in] |
dataright | [in] |
weight | [in] |
ek | [out] is the stiffness matrix |
ef | [out] is the load vector |
datavec | [in] Data material vector |
var | [in] number of solution variables. See NSolutionVariables() method |
Solout | [out] is the solution vector |
Reimplemented from TPZMaterial.
Definition at line 750 of file pzmatmixedpoisson3d.cpp.
References TPZVec< T >::size().
Referenced by ContributeBC().
|
inlineoverridevirtual |
Gets the order of the integration rule necessary to integrate an element with polinomial order p.
Reimplemented from TPZMaterial.
Definition at line 162 of file pzmatmixedpoisson3d.h.
|
inlineoverridevirtual |
Gets the order of the integration rule necessary to integrate an element multiphysic.
Reimplemented from TPZMaterial.
Definition at line 168 of file pzmatmixedpoisson3d.h.
References Contribute(), ContributeWithoutSecondIntegration(), TPZMaterial::fForcingFunction, TPZFunction< TVar >::PolynomialOrder(), and Print().
|
inline |
Definition at line 320 of file pzmatmixedpoisson3d.h.
References ClassId(), and fSecondIntegration.
|
inline |
Definition at line 93 of file pzmatmixedpoisson3d.h.
References fMatId.
|
inlineoverridevirtual |
Returns the name of the material.
Reimplemented from TPZMaterial.
Definition at line 84 of file pzmatmixedpoisson3d.h.
Referenced by Print().
|
inlineoverridevirtual |
To create another material of the same type.
Reimplemented from TPZMaterial.
Definition at line 98 of file pzmatmixedpoisson3d.h.
References NStateVariables(), and TPZMatMixedPoisson3D().
|
overridevirtual |
Returns the number of variables associated with the variable indexed by var.
var | Index variable into the solution, is obtained by calling VariableIndex |
Reimplemented from TPZMaterial.
Definition at line 570 of file pzmatmixedpoisson3d.cpp.
References fDim, and TPZMaterial::NSolutionVariables().
Referenced by ContributeBC(), and Solution().
|
overridevirtual |
Returns the number of state variables associated with the material.
Implements TPZMaterial.
Definition at line 114 of file pzmatmixedpoisson3d.cpp.
Referenced by Hdiv2dPaper201504::CMeshFlux(), hdivCurvedJCompAppMath::CMeshFlux(), Hdiv3dPaper201504::CMeshFlux(), Hdiv2dPaper201504::CMeshPressure(), hdivCurvedJCompAppMath::CMeshPressure(), Hdiv3dPaper201504::CMeshPressure(), and NewMaterial().
TPZMatMixedPoisson3D & TPZMatMixedPoisson3D::operator= | ( | const TPZMatMixedPoisson3D & | copy | ) |
Definition at line 98 of file pzmatmixedpoisson3d.cpp.
References falpha, fDim, fF, fInvK, fMatId, fmatLagr, fReactionTerm, fSecondIntegration, fTensorK, and TPZMaterial::operator=().
Referenced by TPZMatMixedPoisson3D().
|
inline |
Definition at line 144 of file pzmatmixedpoisson3d.h.
References fPermeabilityFunction.
|
overridevirtual |
Prints out the data associated with the material.
Reimplemented from TPZMaterial.
Definition at line 118 of file pzmatmixedpoisson3d.cpp.
References fDim, fF, fMatId, Name(), and TPZMaterial::Print().
Referenced by IntegrationRuleOrder().
|
inline |
Definition at line 154 of file pzmatmixedpoisson3d.h.
References fReactionTermFunction.
|
inline |
Definition at line 88 of file pzmatmixedpoisson3d.h.
|
inline |
Definition at line 149 of file pzmatmixedpoisson3d.h.
|
inline |
Definition at line 135 of file pzmatmixedpoisson3d.h.
|
inline |
Definition at line 104 of file pzmatmixedpoisson3d.h.
References fDim, fInvK, fTensorK, and TPZFMatrix< TVar >::Resize().
|
inline |
Definition at line 139 of file pzmatmixedpoisson3d.h.
|
inline |
Definition at line 115 of file pzmatmixedpoisson3d.h.
References TPZMatrix< TVar >::Cols(), DebugStop, fDim, and TPZMatrix< TVar >::Rows().
Referenced by Hdiv2dPaper201504::CMeshMixed(), hdivCurvedJCompAppMath::CMeshMixed(), and Hdiv3dPaper201504::CMeshMixed().
|
inline |
Definition at line 125 of file pzmatmixedpoisson3d.h.
|
inline |
Definition at line 131 of file pzmatmixedpoisson3d.h.
|
overridevirtual |
Returns the solution associated with the var index based on the finite element approximation.
Reimplemented from TPZMaterial.
Definition at line 584 of file pzmatmixedpoisson3d.cpp.
References TPZAxesTools< TVar >::Axes2XYZ(), Dimension(), TPZFunction< TVar >::Execute(), fDim, fF, TPZMaterial::fForcingFunction, TPZMaterial::fForcingFunctionExact, NSolutionVariables(), test::res, TPZVec< T >::Resize(), and val().
Referenced by ContributeBC(), Errors(), and ErrorsHdiv().
|
overridevirtual |
Reimplemented from TPZMaterial.
Definition at line 726 of file pzmatmixedpoisson3d.cpp.
References TPZAxesTools< TVar >::Axes2XYZ(), fDim, NSolutionVariables(), TPZVec< T >::Resize(), and TPZMaterial::Solution().
|
overridevirtual |
Returns the solution associated with the var index based on the finite element approximation.
Reimplemented from TPZMaterial.
Definition at line 699 of file pzmatmixedpoisson3d.cpp.
References TPZVec< T >::Resize(), and TPZMaterialData::sol.
|
inline |
Definition at line 316 of file pzmatmixedpoisson3d.h.
|
overridevirtual |
Returns the variable index associated with the name
Reimplemented from TPZMaterial.
Definition at line 553 of file pzmatmixedpoisson3d.cpp.
References TPZMaterial::VariableIndex().
Referenced by ContributeBC().
|
protected |
Coeficiente que multiplica a pressão: termo de reacao
Definition at line 46 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), operator=(), and TPZMatMixedPoisson3D().
|
protected |
Dimensao do dominio
Definition at line 43 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), Dimension(), Errors(), NSolutionVariables(), operator=(), Print(), SetPermeability(), SetPermeabilityTensor(), Solution(), and TPZMatMixedPoisson3D().
|
protected |
Valor da funcao de carga
Definition at line 40 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), operator=(), Print(), Solution(), and TPZMatMixedPoisson3D().
|
protected |
inverse of the permeability tensor.
Definition at line 55 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), operator=(), SetPermeability(), and TPZMatMixedPoisson3D().
|
protected |
Material Id
Definition at line 37 of file pzmatmixedpoisson3d.h.
Referenced by MatId(), operator=(), Print(), and TPZMatMixedPoisson3D().
|
protected |
Definition at line 64 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), operator=(), and TPZMatMixedPoisson3D().
|
protected |
Pointer to forcing function, it is the Permeability and its inverse.
Definition at line 58 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), PermeabilityFunction(), and TPZMatMixedPoisson3D().
|
protected |
Definition at line 69 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), operator=(), and TPZMatMixedPoisson3D().
|
protected |
Pointer to forcing function, it is the reaction term.
Definition at line 61 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), ReactionTermFunction(), and TPZMatMixedPoisson3D().
|
protected |
Parameter to choose the second integration by parts in the variational formulation.
Definition at line 67 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), IsUsedSecondIntegration(), operator=(), and TPZMatMixedPoisson3D().
|
protected |
permeability tensor. Coeficient which multiplies the gradient operator
Definition at line 52 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), operator=(), SetPermeability(), and TPZMatMixedPoisson3D().
|
protected |
fluid viscosity
Definition at line 49 of file pzmatmixedpoisson3d.h.
Referenced by Contribute(), ContributeWithoutSecondIntegration(), and TPZMatMixedPoisson3D().