9 #ifndef __PZ__pzconvectionproblem__ 10 #define __PZ__pzconvectionproblem__ 76 virtual void Print(std::ostream & out)
override;
78 virtual std::string
Name()
override {
return "TPZMatConvectionProblem"; }
108 fTimeValue = TimeValue;
117 fRungeKuttaTwo =
true;
184 virtual int ClassId()
const override;
Defines the interface which material objects need to implement for discontinuous Galerkin formulation...
void SetTimeStep(REAL delt)
virtual int NSolutionVariables(int var) override
Returns the number of variables associated with the variable indexed by var.
REAL fXf
Forcing function value.
virtual ~TPZMatConvectionProblem()
clarg::argBool bc("-bc", "binary checkpoints", false)
TPZMatConvectionProblem & operator=(const TPZMatConvectionProblem ©)
void SetInternalFlux(REAL flux)
void GetParameters(REAL &rho, TPZVec< REAL > &convdir)
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...
void Errors(TPZVec< REAL > &x, TPZVec< STATE > &u, TPZFMatrix< STATE > &dudx, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux, 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 o...
virtual int NStateVariables() const override
Returns the number of state variables associated with the material.
int Dimension() const override
Returns the integrable dimension of the material.
virtual int ClassId() const override
Unique identifier for serialization purposes.
void GetTimeValue(REAL &TimeValue)
void SetParameters(REAL rho, TPZVec< REAL > &convdir)
REAL fRho
Coeficient which multiplies the temporal derivative.
virtual std::string Name() override
Returns the name of the material.
void SetTimeValue(REAL TimeValue)
TPZMatConvectionProblem()
virtual int VariableIndex(const std::string &name) override
Contains the TPZMatConvectionProblem class which implements a convection problem 2D with time depende...
This class defines the boundary condition for TPZMaterial objects.
virtual void Solution(TPZMaterialData &data, int var, TPZVec< STATE > &Solout) override
It return a solution to multiphysics simulation.
virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
It computes a contribution to stiffness matrix and load vector at one BC integration point...
virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override
It computes a contribution to stiffness matrix and load vector at one BC integration point...
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...
virtual void Print(std::ostream &out) override
Prints out the data associated with the material.
Contains the TPZDiscontinuousGalerkin class which implements the interface for discontinuous Galerkin...
void SetTrueRungeKuttaTwo()
int fDim
Problem dimension.
TPZVec< REAL > fConvDir
convection term (direction velocity)