NeoPZ
pzmultiphase.h
Go to the documentation of this file.
1 //
2 // pzmultiphase.h
3 // PZ
4 //
5 // Created by Omar Duran on 19/08/2013.
6 // Copyright (c) 2012 LabMec-Unicamp. All rights reserved.
7 //
8 
9 #ifndef PZ_pzmultiphase_h
10 #define PZ_pzmultiphase_h
11 
12 #include "TPZMaterial.h"
13 #include "pzdiscgal.h"
14 #ifdef _AUTODIFF
15 #include "fad.h"
16 #endif
17 #include <iostream>
18 #include <fstream>
19 #include <string>
20 
31 
32 protected:
33 
35  int fDim;
36 
38  int fmatId;
39 
41  REAL ff;
42 
44  REAL fxi;
45 
47  enum EState { ELastState = 0, ECurrentState = 1 };
49 
50 #ifdef _AUTODIFF
51 
52  typedef TFad<2, REAL> BFadREAL;
53 
54 #endif
55 
56 public:
57 
58  bool fnewWS;
59 
60  TPZMultiphase();
61 
62  TPZMultiphase(int matid, int dim);
63 
64  virtual ~TPZMultiphase();
65 
67  TPZMultiphase(const TPZMultiphase &copy);
68 
70 
71  virtual void Print(std::ostream & out) override;
72 
73  virtual std::string Name() override { return "TPZMultiphase"; }
74 
75  virtual int Dimension() const override;
76 
77  virtual int NStateVariables() const override;
78 
79  virtual int MatId();
80 
81 
82  virtual void Contribute(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ef) override;
83 
84  virtual void Contribute(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
85 
86  virtual void Contribute(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
87 
88 
89  virtual void ContributeBC(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override;
90 
91  virtual void ContributeBC(TPZMaterialData &data, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override;
92 
93  virtual void ContributeBC(TPZVec<TPZMaterialData> &datavec, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc) override;
94 
95 
96  virtual void ContributeInterface(TPZMaterialData &data, TPZMaterialData &dataleft, TPZMaterialData &dataright, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef) override;
97 
98  virtual void ContributeInterface(TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, TPZVec<TPZMaterialData> &dataright, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef) override;
99 
100  virtual void ContributeInterface(TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, TPZVec<TPZMaterialData> &dataright, REAL weight, TPZFMatrix<STATE> &ef) override;
101 
102  virtual void ContributeInterface(TPZVec<TPZMaterialData> &datavec,TPZVec<TPZMaterialData> &dataleftvec,TPZVec<TPZMaterialData> &datarightvec,REAL weight,TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef) override;
103 
104 
105  virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix<STATE> &ef,TPZBndCond &bc) override;
106 
107  virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc) override;
108 
109  virtual void ContributeBCInterface(TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleftvec, REAL weight, TPZFMatrix<STATE> &ek, TPZFMatrix<STATE> &ef, TPZBndCond &bc) override;
110 
111  // Contribution for each state variable
112  virtual void ApplyUxD (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc);
113  virtual void ApplyUyD (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc);
114  virtual void ApplySigmaN (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc);
115  virtual void ApplyQnD (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc);
116  virtual void ApplyPN (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc);
117  virtual void ApplySin (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc);
118  virtual void ApplySout (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ek,TPZFMatrix<STATE> &ef,TPZBndCond &bc);
119  virtual void ApplySin (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ef,TPZBndCond &bc);
120  virtual void ApplySout (TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleft, REAL weight, TPZFMatrix<STATE> &ef,TPZBndCond &bc);
121 
122  virtual void FillDataRequirements(TPZVec<TPZMaterialData > &datavec) override;
123 
124  virtual void FillBoundaryConditionDataRequirement(int type,TPZVec<TPZMaterialData > &datavec) override;
125 
126 
127  virtual int VariableIndex(const std::string &name) override;
128 
129  virtual int NSolutionVariables(int var) override;
130 
131  virtual void Solution(TPZVec<TPZMaterialData> &datavec, int var, TPZVec<STATE> &Solout) override;
132 
133 
134  virtual void Solution(TPZMaterialData &data, TPZVec<TPZMaterialData> &dataleftvec, TPZVec<TPZMaterialData> &datarightvec, int var, TPZVec<STATE> &Solout, TPZCompEl * Left, TPZCompEl * Right) override
135  {
136  TPZDiscontinuousGalerkin::Solution(data,dataleftvec,datarightvec,var,Solout,Left,Right);
137  }
138 
139 
140  // Here is needed to spefify which models are used in this bi-phasic model see Advanced-Petroleum-Reservoir-Simulation M. Rafiqul Islam
141 
143  REAL fLref;
144 
146  REAL fKref;
147 
149  REAL fQref;
150 
152  REAL fPref;
153 
155  REAL fRhoref;
156 
158  REAL fEtaref;
159 
160 
163 
166 
168  bool fYorN;
169 
171  REAL fDeltaT;
172 
174  REAL fTime;
175 
177  REAL fTheta;
178 
180  REAL fGamma;
181 
183  void SetXiBalance(REAL xi){ this->fxi = xi;}
184 
186  void SetTimeStep(REAL timestep){ this->fDeltaT = timestep;}
187 
189  void SetTime(REAL time){ this->fTime = time;}
190 
192  void SetTScheme(REAL timegamma, REAL timetheta){ this->fTheta = timetheta; this->fGamma = timegamma;}
193 
194  void SetLastState(){ gState = ELastState;}
195 
196  void SetCurrentState(){ gState = ECurrentState;}
197 
199  void SetYorN(bool dummybool){ this->fYorN = dummybool;}
200 
205  STATE LameLambda();
206 
211  STATE LameLambdaU();
212 
217  STATE LameMu();
218 
223  STATE BiotAlpha();
224 
229  STATE Se();
230 
231 
237  void SetKMap(TPZStack< TPZFMatrix<REAL> > KabsoluteMap)
238  {
239  fKabsoluteMap = KabsoluteMap;
240  }
241 
243  void CapillaryPressure(REAL So, REAL &pc, REAL &DpcDSo);
244 
249  void Kro(REAL Sw, REAL &Kro, REAL &dKroDSw);
250 
255  void Krw(REAL Sw, REAL &Krw, REAL &dKrwSo);
256 
261  void SWaterstar(REAL &Swstar, REAL &Po, REAL &Sw);
262 
267  void fWaterstar(REAL &fWstar, REAL Pw, REAL Sw, REAL &dfWstarDPw, REAL &dfWstarDSw);
268 
273  void fOilstar(REAL &fOstar, REAL Pw, REAL Sw, REAL &dfOstarDPw, REAL &dfOstarDSw);
274 
279  void fstar(REAL &fStar, REAL Pw, REAL Sw, REAL Gdotn, REAL &dfstarDPw, REAL &dfstarDSw);
280 
285  void Porosity(REAL po, REAL &poros, REAL &dPorosDp);
286 
291  void RhoOil(REAL po, REAL &RhoOil, REAL &dRhoOilDpo);
292 
297  void RhoWater(REAL pw, REAL &RhoWater, REAL &dRhoWaterDpo);
298 
303  void OilViscosity(REAL po, REAL &OilViscosity, REAL &dOilViscosityDpo);
304 
309  void WaterViscosity(REAL po, REAL &WaterViscosity, REAL &dWaterViscosityDpo);
310 
315  void OilLabmda(REAL &OilLabmda, REAL Po, REAL Sw, REAL &dOilLabmdaDPo, REAL &dOilLabmdaDSw);
316 
321  void WaterLabmda(REAL &WaterLabmda, REAL Pw, REAL Sw, REAL &dWaterLabmdaDPw, REAL &dWaterLabmdaDSw);
322 
327  void Labmda(REAL &Labmda, REAL Pw, REAL Sw, REAL &dLabmdaDPw, REAL &dLabmdaDSw);
328 
333  void fOil(REAL &fOil, REAL Pw, REAL Sw, REAL &dfOilDPw, REAL &dfOilDSw);
334 
339  void fWater(REAL &fWater, REAL Pw, REAL Sw, REAL &dfWaterDPw, REAL &dfWaterDSw);
340 
341 #ifdef _AUTODIFF
342 
343  // Fad Methods ///////////////////////////////////////////////////////////////////////////////////////////////////////
344 
346  void CapillaryPressure(BFadREAL So, BFadREAL &pc);
347 
352  void Kro(BFadREAL Sw, BFadREAL &Kro);
353 
358  void Krw(BFadREAL Sw, BFadREAL &Krw);
359 
364  void Porosity(BFadREAL po, BFadREAL &poros);
365 
370  void RhoOil(BFadREAL po, BFadREAL &RhoOil);
371 
376  void RhoWater(BFadREAL pw, BFadREAL &RhoWater);
377 
382  void OilViscosity(BFadREAL po, BFadREAL &OilViscosity);
383 
388  void WaterViscosity(BFadREAL po, BFadREAL &WaterViscosity);
389 
394  void OilLabmda(BFadREAL OilLabmda, BFadREAL Po, BFadREAL &Sw);
395 
400  void WaterLabmda(BFadREAL WaterLabmda, BFadREAL Pw, BFadREAL &Sw);
401 
402 
407  void Labmda(BFadREAL Labmda, BFadREAL Pw, BFadREAL &Sw);
408 
413  void fOil(BFadREAL fOil, BFadREAL Pw, BFadREAL &Sw);
414 
415 
420  void fWater(BFadREAL fWater, BFadREAL Pw, BFadREAL &Sw);
421 
422 #endif
423 
424  // Fad Methods ///////////////////////////////////////////////////////////////////////////////////////////////////////
425 
426 
428  void SetLreference(REAL &Lref){ fLref = Lref;}
429 
431  void SetKreference(REAL &Kref){ fKref = Kref;}
432 
434  void SetPreference(REAL &Pref){ fPref = Pref;}
435 
437  void SetQreference(REAL &Qref){ fQref = Qref;}
438 
440  void SetRhoSCreference(REAL &Densityref){ fRhoref = Densityref;}
441 
443  void SetEtaSCreference(REAL &Viscosityref){ fEtaref = Viscosityref;}
444 
446  REAL RhoOilSC();
447 
449  REAL RhoWaterSC();
450 
453 
455  void K(TPZFMatrix<REAL> &Kab);
456 
459 
461  void LoadKMap(std::string MaptoRead);
462 
464  public:
465 virtual int ClassId() const override;
466 
467 };
468 
469 #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.
void OilLabmda(REAL &OilLabmda, REAL Po, REAL Sw, REAL &dOilLabmdaDPo, REAL &dOilLabmdaDSw)
Oil mobility. .
void Labmda(REAL &Labmda, REAL Pw, REAL Sw, REAL &dLabmdaDPw, REAL &dLabmdaDSw)
Bulk mobility. .
REAL fGamma
Parameter representing temporal scheme for conservation equation.
Definition: pzmultiphase.h:180
virtual void Contribute(TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ef) override
It computes a contribution to the stiffness matrix and load vector at one integration point to multip...
virtual int NStateVariables() const override
Returns the number of state variables associated with the material.
virtual int MatId()
virtual std::string Name() override
Returns the name of the material.
Definition: pzmultiphase.h:73
clarg::argBool bc("-bc", "binary checkpoints", false)
void OilViscosity(REAL po, REAL &OilViscosity, REAL &dOilViscosityDpo)
Oil viscosity. .
REAL ff
Definition of constants.
Definition: pzmultiphase.h:41
REAL RhoWaterSC()
Water density on standard conditions - kg/m3.
virtual void FillDataRequirements(TPZVec< TPZMaterialData > &datavec) override
Fill material data parameter with necessary requirements for the Contribute method. Here, in base class, all requirements are considered as necessary. Each derived class may optimize performance by selecting only the necessary data.
EState
State: Stiffness or Mass Matrix Calculations.
Definition: pzmultiphase.h:47
void CapillaryPressure(REAL So, REAL &pc, REAL &DpcDSo)
Capilar pressure. .
void WaterLabmda(REAL &WaterLabmda, REAL Pw, REAL Sw, REAL &dWaterLabmdaDPw, REAL &dWaterLabmdaDSw)
Water mobility. .
void Kro(REAL Sw, REAL &Kro, REAL &dKroDSw)
Oil relative permeability. .
virtual int Dimension() const override
Returns the integrable dimension of the material.
void RhoOil(REAL po, REAL &RhoOil, REAL &dRhoOilDpo)
void fOil(REAL &fOil, REAL Pw, REAL Sw, REAL &dfOilDPw, REAL &dfOilDSw)
Fractional oil flux. .
STATE BiotAlpha()
Biot parameter Parameter. .
REAL fRhoref
density reference - kg/m3
Definition: pzmultiphase.h:155
TPZStack< TPZFMatrix< REAL > > fKabsoluteMap
K map.
Definition: pzmultiphase.h:162
void WaterViscosity(REAL po, REAL &WaterViscosity, REAL &dWaterViscosityDpo)
Water viscosity. .
virtual void ApplyUyD(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
void fstar(REAL &fStar, REAL Pw, REAL Sw, REAL Gdotn, REAL &dfstarDPw, REAL &dfstarDSw)
Fractional product upwind function, Gdotn > 0 means water decrease (dSw/dt < 0), Gdotn < 0 means wate...
This class implements a simple vector storage scheme for a templated class T. Utility.
Definition: pzgeopoint.h:19
REAL fQref
Pressure reference - Pa.
Definition: pzmultiphase.h:149
virtual void ApplySin(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
virtual int ClassId() const override
Unique identifier for serialization purposes.
int fDim
Problem dimension.
Definition: pzmultiphase.h:35
REAL fLref
Characteristic length - m.
Definition: pzmultiphase.h:143
int fmatId
Material id.
Definition: pzmultiphase.h:38
virtual void ApplySout(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
void SetKreference(REAL &Kref)
Set permeability reference - m2.
Definition: pzmultiphase.h:431
void SetTimeStep(REAL timestep)
Defines simulation time step.
Definition: pzmultiphase.h:186
virtual void ApplySigmaN(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
STATE LameMu()
Lame Second Parameter. .
bool fYorN
Use or not K map.
Definition: pzmultiphase.h:168
void SetTime(REAL time)
Defines simulation time step.
Definition: pzmultiphase.h:189
virtual 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.
REAL fEtaref
viscosity reference - Pa s
Definition: pzmultiphase.h:158
virtual void Solution(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleftvec, TPZVec< TPZMaterialData > &datarightvec, int var, TPZVec< STATE > &Solout, TPZCompEl *Left, TPZCompEl *Right) override
Returns the solution associated with the var index based on the finite element approximation around o...
Definition: pzmultiphase.h:134
virtual void ApplyPN(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
virtual int VariableIndex(const std::string &name) override
This class defines the boundary condition for TPZMaterial objects.
Definition: pzbndcond.h:29
Definition: tfad.h:64
REAL fPref
Pressure reference - Pa.
Definition: pzmultiphase.h:152
void fWater(REAL &fWater, REAL Pw, REAL Sw, REAL &dfWaterDPw, REAL &dfWaterDSw)
Fractional water flux. .
virtual void ContributeBCInterface(TPZMaterialData &data, TPZMaterialData &dataleft, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
It computes a contribution to residual vector at one BC integration point.
void SetRhoSCreference(REAL &Densityref)
Set density reference - kg/m3.
Definition: pzmultiphase.h:440
Material to solve a 2d multiphase transport problem by multiphysics simulation.
Definition: pzmultiphase.h:30
REAL fDeltaT
Simulation time step.
Definition: pzmultiphase.h:171
void SetKMap(TPZStack< TPZFMatrix< REAL > > KabsoluteMap)
Definition: pzmultiphase.h:237
STATE Se()
//Se o 1/M coeficiente poroelastico de armazenamento a volume constante.
virtual void Print(std::ostream &out) override
Prints out the data associated with the material.
virtual void ApplyUxD(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
STATE LameLambdaU()
Undrained Lame First Parameter. .
void fOilstar(REAL &fOstar, REAL Pw, REAL Sw, REAL &dfOstarDPw, REAL &dfOstarDSw)
Fractional product function, oil decrease direction (dSw/dt > 0). .
void Porosity(REAL po, REAL &poros, REAL &dPorosDp)
Rock porosity. .
void fWaterstar(REAL &fWstar, REAL Pw, REAL Sw, REAL &dfWstarDPw, REAL &dfWstarDSw)
Fractional product function, water decrease direction (dSw/dt < 0). .
REAL fKref
Permeability reference - m2.
Definition: pzmultiphase.h:146
virtual ~TPZMultiphase()
void SetTScheme(REAL timegamma, REAL timetheta)
Defines stemporal scheme.
Definition: pzmultiphase.h:192
TPZMultiphase & operator=(const TPZMultiphase &copy)
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 integration point...
void SetQreference(REAL &Qref)
Set flux reference - Pa.
Definition: pzmultiphase.h:437
REAL RhoOilSC()
Oil density on standard conditions - kg/m3.
This class implements a stack object. Utility.
Definition: pzcheckmesh.h:14
void Krw(REAL Sw, REAL &Krw, REAL &dKrwSo)
Water relative permeability. .
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 t...
void SetLastState()
Definition: pzmultiphase.h:194
void SetLreference(REAL &Lref)
Set characteristic length - m.
Definition: pzmultiphase.h:428
virtual void ContributeBC(TPZVec< TPZMaterialData > &datavec, REAL weight, TPZFMatrix< STATE > &ef, TPZBndCond &bc) override
It computes a contribution to the stiffness matrix and load vector at one BC integration point to mul...
virtual void ApplyQnD(TPZMaterialData &data, TPZVec< TPZMaterialData > &dataleft, REAL weight, TPZFMatrix< STATE > &ek, TPZFMatrix< STATE > &ef, TPZBndCond &bc)
void LoadKMap(std::string MaptoRead)
Absolute permeability.
void K(TPZFMatrix< REAL > &Kab)
Absolute permeability.
Contains the TPZDiscontinuousGalerkin class which implements the interface for discontinuous Galerkin...
REAL fTheta
Parameter representing temporal scheme for transport equation.
Definition: pzmultiphase.h:177
REAL fTime
Simulation current time.
Definition: pzmultiphase.h:174
Defines the interface of a computational element. Computational Element.
Definition: pzcompel.h:59
virtual int NSolutionVariables(int var) override
Returns the number of variables associated with the variable indexed by var.
void SetEtaSCreference(REAL &Viscosityref)
Set viscosity reference - Pa s.
Definition: pzmultiphase.h:443
STATE LameLambda()
Lame First Parameter. .
void SWaterstar(REAL &Swstar, REAL &Po, REAL &Sw)
Water saturation maximum value of the fractional flow product function. .
TPZFMatrix< REAL > Gravity()
Gravity.
TPZFMatrix< REAL > Kinv(TPZFMatrix< REAL > &Kab)
Absolute permeability inverse.
void SetPreference(REAL &Pref)
Set pressure reference - Pa.
Definition: pzmultiphase.h:434
REAL fxi
Big number balance.
Definition: pzmultiphase.h:44
void SetXiBalance(REAL xi)
Defines simulation bc numeric balance.
Definition: pzmultiphase.h:183
void SetYorN(bool dummybool)
Defines simulation use a K map.
Definition: pzmultiphase.h:199
int fPlaneStress
plane stress condition
Definition: pzmultiphase.h:165
void RhoWater(REAL pw, REAL &RhoWater, REAL &dRhoWaterDpo)
void SetCurrentState()
Definition: pzmultiphase.h:196