NeoPZ
TPZEuler.h
Go to the documentation of this file.
1 
6 #ifndef TPZEULER_H
7 #define TPZEULER_H
8 
9 #include "TPZMaterial.h"
10 #include "eulerdif.h"
11 #include "pzfmatrix.h"
12 
13 class TPZBndCond;
14 class TPZEuler;
15 
16 
22 class TPZEuler : public TPZMaterial {
23 public:
24 
25 virtual int ClassId() const override;
26 
28  TPZEuler(TPZEuler & copy);
30  TPZEuler(int id, STATE deltat) ;
31 
36  void SetState(int state) {
37  fState = state;
38  }
39 
40 
42  virtual int Dimension() const override;
43 
45  virtual int NStateVariables() const override ;
46 
48  virtual int NFluxes() override {return 2;}
49 
56  virtual void Contribute(TPZMaterialData &data, REAL weight,
57  TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef) override;
59  virtual void Contribute(TPZMaterialData &data, REAL weight,
60  TPZFMatrix<STATE> &ef) override
61  {
62  TPZMaterial::Contribute(data,weight,ef);
63  }
64 
65  virtual void ContributeBC(TPZMaterialData &data,REAL weight,
67  virtual void ContributeBC(TPZMaterialData &data,REAL weight,
68  TPZFMatrix<STATE> &ef,TPZBndCond &bc) override
69  {
70  TPZMaterial::ContributeBC(data,weight,ef,bc);
71  }
72 
76  virtual void Print(std::ostream &out = std::cout) override;
77 
79  virtual int VariableIndex(const std::string &name) override;
80 
81  virtual int NSolutionVariables(int var) override;
82 protected:
83  virtual void Solution(TPZVec<STATE> &Sol,TPZFMatrix<STATE> &DSol,TPZFMatrix<REAL> &axes,int var,TPZVec<STATE> &Solout) override;
84 public:
85  virtual void Solution(TPZMaterialData &data,int var,TPZVec<STATE> &Solout) override
86  {
87  int numbersol = data.sol.size();
88  if (numbersol != 1) {
89  DebugStop();
90  }
91 
92  Solution(data.sol[0],data.dsol[0],data.axes,var,Solout);
93  }
94 
95 
96  virtual void Flux(TPZVec<REAL> &x, TPZVec<STATE> &Sol, TPZFMatrix<STATE> &DSol, TPZFMatrix<REAL> &axes, TPZVec<STATE> &flux) override {}
97 
99  virtual TPZMaterial * NewMaterial() override;
100 
102  virtual void SetData(std::istream &data) override;
103 
104 private:
105 
107  STATE fDeltaT;
108  int fState;
109 };
110 
111 #endif
static TEulerDiffusivity gEul
Definition: TPZEuler.h:106
virtual int NSolutionVariables(int var) override
Returns the number of variables associated with the variable indexed by var.
Definition: TPZEuler.cpp:36
clarg::argBool bc("-bc", "binary checkpoints", false)
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override
Computes contribution to the stiffness matrix and right hand side at an integration point...
Definition: TPZEuler.cpp:166
virtual int NFluxes() override
Return the number of components which form the flux function.
Definition: TPZEuler.h:48
TPZGradSolVec dsol
vector of the derivatives of the solution at the integration point
TPZEuler(TPZEuler &copy)
Copy constructor.
Definition: TPZEuler.cpp:280
virtual void Solution(TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, int var, TPZVec< STATE > &Solout) override
Definition: TPZEuler.cpp:22
virtual int ClassId() const override
Define the class id associated with the class.
Definition: TPZEuler.cpp:291
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)=0
It computes a contribution to the stiffness matrix and load vector at one integration point...
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
virtual int VariableIndex(const std::string &name) override
Returns the variable index associated with the name.
Definition: TPZEuler.cpp:42
int64_t size() const
Returns the number of elements of the vector.
Definition: pzvec.h:196
Contains TPZMatrixclass which implements full matrix (using column major representation).
#define DebugStop()
Returns a message to user put a breakpoint in.
Definition: pzerror.h:20
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
This class implements a linear scalar convection equation with modified SUPG difusion.
Definition: TPZEuler.h:22
Implements a numerical diffusivity coeficient for the SUPG method. Analysis: Solving process Analysis...
Definition: eulerdif.h:40
TPZFNMatrix< 9, REAL > axes
axes indicating the directions of the derivatives of the shapefunctions
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ef) override
Computes contribution to the right hand side at an integration point.
Definition: TPZEuler.h:59
void SetState(int state)
Set the state of the Euler material.
Definition: TPZEuler.h:36
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)=0
It computes a contribution to the stiffness matrix and load vector at one BC integration point...
virtual void SetData(std::istream &data) override
Reads data of the material from a istream (file data)
Definition: TPZEuler.cpp:13
virtual void Flux(TPZVec< REAL > &x, TPZVec< STATE > &Sol, TPZFMatrix< STATE > &DSol, TPZFMatrix< REAL > &axes, TPZVec< STATE > &flux) override
Computes the value of the flux function to be used by ZZ error estimator.
Definition: TPZEuler.h:96
virtual void Print(std::ostream &out=std::cout) override
Print out the data associated with the material.
Definition: TPZEuler.cpp:48
Contains the TEulerDiffusivity class which implements a numerical diffusivity coefficient for SUPG...
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.
Definition: TPZEuler.h:85
virtual TPZMaterial * NewMaterial() override
To create another material of the same type.
Definition: TPZEuler.cpp:17
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...
Definition: TPZEuler.cpp:51
int fState
Definition: TPZEuler.h:108
virtual int Dimension() const override
Returns the integrable dimension of the material.
Definition: TPZEuler.cpp:277
STATE fDeltaT
Definition: TPZEuler.h:107
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...
Definition: TPZEuler.h:67
TPZSolVec sol
vector of the solutions at the integration point
virtual int NStateVariables() const override
Returns the number of state variables associated with the material.
Definition: TPZEuler.cpp:274