NeoPZ
pzblackoil2p3d.h
Go to the documentation of this file.
1 
6 #ifndef PZBLACKOIL2P3D_H
7 #define PZBLACKOIL2P3D_H
8 
9 #include "TPZMaterial.h"
10 #include "pzdiscgal.h"
11 
12 #ifdef _AUTODIFF
13 
14 #include "fad.h"
15 
16 
22 class TPZBlackOil2P3D : public TPZDiscontinuousGalerkin{
23 
24 public:
25 
26  typedef TFad<4, REAL> BFadREAL;
27 
28 protected:
29 
31  void Interpolate(std::map<REAL,REAL> &dados, double x, double &y, double &dy);
32  void Interpolate(std::map<REAL,REAL> &dados, BFadREAL x, BFadREAL &y);
33 
35  double fDeltaT;
36 
38  enum EState { ELastState = 0, ECurrentState = 1 };
39 
40  static EState gState;
41 
42  void testedoBo();
43 
44  void testeKrw();
45 
46 public:
47 
48  static void SetLastState(){ gState = ELastState; }
49  static void SetCurrentState(){ gState = ECurrentState; }
50 
56  TPZBlackOil2P3D(int id, double deltaT);
57 
59  ~TPZBlackOil2P3D();
60 
62  TPZBlackOil2P3D(const TPZBlackOil2P3D &cp);
63 
65  void SetTimeStep(double timestep){
66  this->fDeltaT = timestep;
67  }
68 
70  virtual TPZMaterial * NewMaterial() override;
71 
73  virtual int Dimension() const override { return 3; }
74 
76  virtual int NStateVariables() const override{ return 2; }
77 
79  virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
80 
89  virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override;
90 
92  virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
93 
95  virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc) override;
96 
97 
104  enum ESolutionVars { ENone = 0, EWaterPressure = 1, EOilPressure, EWaterSaturation, EOilSaturation, EDarcyVelocity };
105 
107  virtual int VariableIndex(const std::string &name) override;
108 
113  virtual int NSolutionVariables(int var) override;
114 
119  virtual void Solution(TPZVec<STATE> &Sol, TPZFMatrix<STATE> &DSol,
120  TPZFMatrix<REAL> &axes, int var, TPZVec<STATE> &Solout) override;
128  virtual void FillDataRequirements(TPZMaterialData &data) override {
129  data.SetAllRequirements(true);
130  data.fNeedsNeighborSol = false;
131  data.fNeedsNeighborCenter = false;
132  }
133 
139  virtual void FillDataRequirementsInterface(TPZMaterialData &data) override {
140  data.SetAllRequirements(true);
141  data.fNeedsSol = false;
142  }
143 
153  void Kro(double So, double &Kro, double &dKroSo);
154  void Kro(BFadREAL So, BFadREAL &Kro);
155 
160  void Krw(double So, double &Krw, double &dKrwSo);
161  void Krw(BFadREAL So, BFadREAL &Krw);
162 
167  void Bo(double po, double &Bo, double &dBoDpo);
168  void Bo(BFadREAL po, BFadREAL &Bo);
169 
171  void ViscOleo(double po, double &ViscOleo, double &dViscOleoDpo);
172  void ViscOleo(BFadREAL po, BFadREAL &ViscOleo);
173 
175  void PressaoCapilar(double So, double &pc, double &DpcDSo);
176  void PressaoCapilar(BFadREAL So, BFadREAL &pc);
177 
179  void Porosidade(double po, double &poros, double &dPorosDpo);
180  void Porosidade(BFadREAL po, BFadREAL &poros);
181 
183  double RhoOleoSC();
184 
186  double RhoAguaSC();
187 
189  double g();
190 
192  double Bw();
193 
195  double ViscAgua();
196 
198  void K(TPZFMatrix<REAL> &K);
199 
201 };
202 
203 #endif
204 
205 #endif
Defines the interface which material objects need to implement for discontinuous Galerkin formulation...
Definition: pzdiscgal.h:20
virtual void Solution(TPZMaterialData &data, int var, TPZVec< STATE > &Solout)
Returns the solution associated with the var index based on the finite element approximation.
clarg::argBool bc("-bc", "binary checkpoints", false)
virtual TPZMaterial * NewMaterial()
To create another material of the same type.
virtual int VariableIndex(const std::string &name)
Returns the variable index associated with the name.
void SetAllRequirements(bool set)
Set all flags at once.
virtual int Dimension() const =0
Returns the integrable dimension of the material.
virtual int NStateVariables() const =0
Returns the number of state variables associated with the material.
virtual void FillDataRequirementsInterface(TPZMaterialData &data) override
Fill material data parameter with necessary requirements for the ContributeInterface method...
Definition: pzdiscgal.cpp:22
This abstract class defines the behaviour which each derived class needs to implement.
Definition: TPZMaterial.h:39
virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override=0
It computes a contribution to the stiffness matrix and load vector at one BC integration point...
virtual void FillDataRequirements(TPZMaterialData &data)
Fill material data parameter with necessary requirements for the.
Definition: TPZMaterial.cpp:81
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
Definition: tfad.h:64
virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef)
It computes a contribution to stiffness matrix and load vector at one integration point...
Definition: pzdiscgal.cpp:30
virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)=0
It computes a contribution to stiffness matrix and load vector at one BC integration point...
virtual int NSolutionVariables(int var)
Returns the number of variables associated with the variable indexed by var.
Contains the TPZDiscontinuousGalerkin class which implements the interface for discontinuous Galerkin...
virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef) override=0
It computes a contribution to the stiffness matrix and load vector at one integration point...