NeoPZ
pzconvectionproblem.h
Go to the documentation of this file.
1 //
2 // pzconvectionproblem.h
3 // PZ
4 //
5 // Created by Agnaldo Farias on 4/2/13.
6 //
7 //
8 
9 #ifndef __PZ__pzconvectionproblem__
10 #define __PZ__pzconvectionproblem__
11 
12 #include "pzdiscgal.h"
13 //#include "pzfmatrix.h"
14 //#include "pzmaterial.h"
15 //#include "pzvec.h"
16 
17 #include <iostream>
18 
19 
26 /*
27  * \f$ phi*du/dt + div(ConvDir*u) = f (Eq. 1) \f$
28  *
29  */
30 
31 
33 
34 protected:
36  REAL fXf;
37 
39  int fDim;
40 
42  REAL fRho;
43 
46 
47  //time step
48  REAL fTimeStep;
49 
50  //time of simulation
51  REAL fTimeValue;
52 
53  int fmatId;
54 
55  //Second order Runge-Kutta
57 
59  enum EState { ELastState = 0, ECurrentState = 1 };
61 
62 
63 public:
64 
66 
67  TPZMatConvectionProblem(int matid, int dim);
68 
69  virtual ~TPZMatConvectionProblem();
70 
73 
75 
76  virtual void Print(std::ostream & out) override;
77 
78  virtual std::string Name() override { return "TPZMatConvectionProblem"; }
79 
80  int Dimension() const override {return fDim;}
81 
82  int MatId()
83  {
84  return fmatId;
85  }
86 
87  virtual int NStateVariables() const override;
88 
89  void SetInternalFlux(REAL flux);
90 
91  void SetParameters(REAL rho, TPZVec<REAL> &convdir);
92 
93  void GetParameters(REAL &rho, TPZVec<REAL> &convdir);
94 
95  void SetLastState(){
96  gState = ELastState;
97  }
98 
100  gState = ECurrentState;
101  }
102 
103  void SetTimeStep(REAL delt){
104  fTimeStep = delt;
105  }
106 
107  void SetTimeValue(REAL TimeValue){
108  fTimeValue = TimeValue;
109  }
110 
111  void GetTimeValue(REAL &TimeValue){
112  TimeValue = fTimeValue;
113  }
114 
115 
117  fRungeKuttaTwo = true;
118  }
119 
132  virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
133 
143  virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override;
144 
155  virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc) override;
166  virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef) override;
167 
168  virtual int VariableIndex(const std::string &name) override;
169 
170  virtual int NSolutionVariables(int var) override;
171 
178  virtual void Solution(TPZMaterialData &data, int var, TPZVec<STATE> &Solout) override;
179 
180  void Errors(TPZVec<REAL> &x,TPZVec<STATE> &u,
182  TPZVec<STATE> &u_exact,TPZFMatrix<STATE> &du_exact,TPZVec<REAL> &values) override;
183  public:
184 virtual int ClassId() const override;
185 
186 };
187 
188 #endif /* defined(__PZ__pzconvectionproblem__) */
Defines the interface which material objects need to implement for discontinuous Galerkin formulation...
Definition: pzdiscgal.h:20
virtual int NSolutionVariables(int var) override
Returns the number of variables associated with the variable indexed by var.
REAL fXf
Forcing function value.
clarg::argBool bc("-bc", "binary checkpoints", false)
TPZMatConvectionProblem & operator=(const TPZMatConvectionProblem &copy)
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)
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.
Definition: pzbndcond.h:29
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.
def values
Definition: rdt.py:119
Contains the TPZDiscontinuousGalerkin class which implements the interface for discontinuous Galerkin...
int fDim
Problem dimension.
TPZVec< REAL > fConvDir
convection term (direction velocity)