OpenCAEPoro  v0.5.0
A simulator for multicomponent porous media flow
OCPFluidMethod.hpp
Go to the documentation of this file.
1 
12 #ifndef __OCPFLUIDMETHOD_HEADER__
13 #define __OCPFLUIDMETHOD_HEADER__
14 
15 // OpenCAEPoro header files
16 #include "LinearSystem.hpp"
17 #include "OCPControl.hpp"
18 #include "Reservoir.hpp"
19 #include "UtilOutput.hpp"
20 #include "UtilTiming.hpp"
21 
23 {
24 public:
25  void InitRock(Bulk& bk) const;
26  void CalRock(Bulk& bk) const;
27 };
28 
30 class IsoT_IMPEC : virtual public IsothermalMethod
31 {
32 public:
34  void Setup(Reservoir& rs, LinearSystem& ls, const OCPControl& ctrl);
36  void InitReservoir(Reservoir& rs) const;
38  void Prepare(Reservoir& rs, OCPControl& ctrl);
40  void AssembleMat(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
42  void SolveLinearSystem(LinearSystem& ls, Reservoir& rs, OCPControl& ctrl);
46  OCP_BOOL FinishNR(const Reservoir& rs);
47  void FinishStep(Reservoir& rs, OCPControl& ctrl);
48 
49 protected:
51  void InitFlash(Bulk& bk) const;
53  void CalKrPc(Bulk& bk) const;
55  void PassFlashValue(Bulk& bk, const OCP_USI& n) const;
56 
57 private:
59  void AllocateReservoir(Reservoir& rs);
61  void
62  AllocateLinearSystem(LinearSystem& ls, const Reservoir& rs, const OCPControl& ctrl);
64  void CalFlash(Bulk& bk);
66  void CalFlux(Reservoir& rs) const;
68  void CalBulkFlux(Reservoir& rs) const;
70  void MassConserve(Reservoir& rs, const OCP_DBL& dt) const;
72  void
73  AssembleMatBulks(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
75  void
76  AssembleMatWells(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
78  void AssembleMatInjWells(LinearSystem& ls,
79  const Bulk& bk,
80  const Well& wl,
81  const OCP_DBL& dt) const;
83  void AssembleMatProdWells(LinearSystem& ls,
84  const Bulk& bk,
85  const Well& wl,
86  const OCP_DBL& dt) const;
88  void GetSolution(Reservoir& rs, const vector<OCP_DBL>& u);
90  void ResetToLastTimeStep01(Reservoir& rs, OCPControl& ctrl);
91  void ResetToLastTimeStep02(Reservoir& rs, OCPControl& ctrl);
92  void ResetToLastTimeStep03(Reservoir& rs, OCPControl& ctrl);
94  void UpdateLastTimeStep(Reservoir& rs) const;
95 };
96 
98 class IsoT_FIM : virtual public IsothermalMethod
99 {
100 public:
102  void Setup(Reservoir& rs, LinearSystem& ls, const OCPControl& ctrl);
104  void InitReservoir(Reservoir& rs) const;
106  void Prepare(Reservoir& rs, const OCP_DBL& dt);
108  void AssembleMat(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
110  void SolveLinearSystem(LinearSystem& ls, Reservoir& rs, OCPControl& ctrl) const;
116  void FinishStep(Reservoir& rs, OCPControl& ctrl);
117 
118 protected:
120  void AllocateReservoir(Reservoir& rs);
122  void
123  AllocateLinearSystem(LinearSystem& ls, const Reservoir& rs, const OCPControl& ctrl);
125  void PassFlashValue(Bulk& bk, const OCP_USI& n) const;
127  void CalKrPc(Bulk& bk) const;
129  void CalRes(Reservoir& rs, const OCP_DBL& dt, const OCP_BOOL& resetRes0) const;
131  void
132  AssembleMatWells(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
133  void
134  AssembleMatWellsNew(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
136  void ResetToLastTimeStep(Reservoir& rs, OCPControl& ctrl);
138  void UpdateLastTimeStep(Reservoir& rs) const;
139 
140 private:
142  void InitFlash(Bulk& bk) const;
144  void CalFlash(Bulk& bk);
146  void
147  AssembleMatBulks(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
149  void AssembleMatInjWells(LinearSystem& ls,
150  const Bulk& bk,
151  const Well& wl,
152  const OCP_DBL& dt) const;
154  void AssembleMatProdWells(LinearSystem& ls,
155  const Bulk& bk,
156  const Well& wl,
157  const OCP_DBL& dt) const;
159  void
160  AssembleMatBulksNew(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
161  void AssembleMatBulksNewS(LinearSystem& ls,
162  const Reservoir& rs,
163  const OCP_DBL& dt) const;
165  void AssembleMatInjWellsNew(LinearSystem& ls,
166  const Bulk& bk,
167  const Well& wl,
168  const OCP_DBL& dt) const;
170  void AssembleMatProdWellsNew(LinearSystem& ls,
171  const Bulk& bk,
172  const Well& wl,
173  const OCP_DBL& dt) const;
175  void
176  GetSolution(Reservoir& rs, const vector<OCP_DBL>& u, const OCPControl& ctrl) const;
177 };
178 
179 class IsoT_FIMn : protected IsoT_FIM
180 {
181 public:
183  void Setup(Reservoir& rs, LinearSystem& ls, const OCPControl& ctrl);
185  void InitReservoir(Reservoir& rs) const;
187  void Prepare(Reservoir& rs, const OCP_DBL& dt);
189  void AssembleMat(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
191  void SolveLinearSystem(LinearSystem& ls, Reservoir& rs, OCPControl& ctrl) const;
197  void FinishStep(Reservoir& rs, OCPControl& ctrl) const;
198 
199 protected:
201  void AllocateReservoir(Reservoir& rs);
203  void InitFlash(Bulk& bk) const;
205  void CalFlash(Bulk& bk);
207  void PassFlashValue(Bulk& bk, const OCP_USI& n) const;
209  void
210  AssembleMatBulksNew(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
212  void
213  AssembleMatWellsNew(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
216  const Bulk& bk,
217  const Well& wl,
218  const OCP_DBL& dt) const;
221  const Bulk& bk,
222  const Well& wl,
223  const OCP_DBL& dt) const;
225  void
226  GetSolution(Reservoir& rs, const vector<OCP_DBL>& u, const OCPControl& ctrl) const;
228  void ResetToLastTimeStep(Reservoir& rs, OCPControl& ctrl);
230  void UpdateLastTimeStep(Reservoir& rs) const;
231 };
232 
233 class IsoT_AIMc : protected IsoT_IMPEC, protected IsoT_FIM
234 {
235 public:
237  void Setup(Reservoir& rs, LinearSystem& ls, const OCPControl& ctrl);
239  void InitReservoir(Reservoir& rs) const;
241  void Prepare(Reservoir& rs, const OCP_DBL& dt);
243  void AssembleMat(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
245  void SolveLinearSystem(LinearSystem& ls, Reservoir& rs, OCPControl& ctrl);
251  void FinishStep(Reservoir& rs, OCPControl& ctrl) const;
252 
253 protected:
255  void AllocateReservoir(Reservoir& rs);
257  void SetFIMBulk(Reservoir& rs);
259  void CalFlashEp(Bulk& bk);
261  void CalFlashEa(Bulk& bk);
263  void CalFlashI(Bulk& bk);
265  void PassFlashValueEp(Bulk& bk, const OCP_USI& n);
267  void CalKrPcE(Bulk& bk);
269  void CalKrPcI(Bulk& bk);
271  void
272  AssembleMatBulks(LinearSystem& ls, const Reservoir& rs, const OCP_DBL& dt) const;
274  void
275  GetSolution(Reservoir& rs, const vector<OCP_DBL>& u, const OCPControl& ctrl) const;
277  void ResetToLastTimeStep(Reservoir& rs, OCPControl& ctrl);
279  void UpdateLastTimeStep(Reservoir& rs) const;
280 };
281 
282 #endif /* end if __OCPFLUIDMETHOD_HEADER__ */
283 
284 /*----------------------------------------------------------------------------*/
285 /* Brief Change History of This File */
286 /*----------------------------------------------------------------------------*/
287 /* Author Date Actions */
288 /*----------------------------------------------------------------------------*/
289 /* Shizhe Li Oct/01/2021 Create file */
290 /* Chensong Zhang Oct/15/2021 Format file */
291 /*----------------------------------------------------------------------------*/
Linear solver class declaration.
double OCP_DBL
Double precision.
Definition: OCPConst.hpp:27
unsigned int OCP_USI
Long unsigned integer.
Definition: OCPConst.hpp:25
unsigned int OCP_BOOL
OCP_BOOL in OCP.
Definition: OCPConst.hpp:29
OCPControl class declaration.
Reservoir class declaration.
Supply basic tools used to output files.
Elapsed wall-time and CPU-cycles declaration.
Physical information of each active reservoir bulk.
Definition: Bulk.hpp:104
void GetSolution(Reservoir &rs, const vector< OCP_DBL > &u, const OCPControl &ctrl) const
Update P, Ni, BHP after linear system is solved.
void AssembleMatBulks(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble linear system for bulks.
void CalFlashEa(Bulk &bk)
Perform flash calculation with Ni for Explicit bulk – Update all properties.
void SolveLinearSystem(LinearSystem &ls, Reservoir &rs, OCPControl &ctrl)
Solve the linear system.
void CalFlashEp(Bulk &bk)
Perform flash calculation with Ni for Explicit bulk – Update partial properties.
void CalKrPcI(Bulk &bk)
Calculate relative permeability and capillary pressure for Implicit bulk.
void Setup(Reservoir &rs, LinearSystem &ls, const OCPControl &ctrl)
Setup AIMc.
void AssembleMat(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble Matrix.
OCP_BOOL FinishNR(Reservoir &rs, OCPControl &ctrl)
Finish a Newton-Raphson iteration.
void CalFlashI(Bulk &bk)
Perform flash calculation with Ni for Implicit bulk.
void SetFIMBulk(Reservoir &rs)
Determine which bulk are treated Implicit.
OCP_BOOL UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
void UpdateLastTimeStep(Reservoir &rs) const
Update values of last step for AIMc.
void PassFlashValueEp(Bulk &bk, const OCP_USI &n)
Pass flash value needed for Explicit bulk – Update partial properties.
void CalKrPcE(Bulk &bk)
Calculate relative permeability and capillary pressure for Explicit bulk.
void ResetToLastTimeStep(Reservoir &rs, OCPControl &ctrl)
Reset variables to last time step.
void Prepare(Reservoir &rs, const OCP_DBL &dt)
Prepare for Assembling matrix.
void InitReservoir(Reservoir &rs) const
Init.
void FinishStep(Reservoir &rs, OCPControl &ctrl) const
Finish a time step.
void AllocateReservoir(Reservoir &rs)
Allocate memory for reservoir.
IsoT_FIM is FIM (Fully Implicit Method).
void Setup(Reservoir &rs, LinearSystem &ls, const OCPControl &ctrl)
Setup FIM.
void AssembleMat(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble Matrix.
void UpdateLastTimeStep(Reservoir &rs) const
Update values of last step for FIM.
void InitReservoir(Reservoir &rs) const
Init.
OCP_BOOL FinishNR(Reservoir &rs, OCPControl &ctrl)
Finish a Newton-Raphson iteration.
void AllocateLinearSystem(LinearSystem &ls, const Reservoir &rs, const OCPControl &ctrl)
Allocate memory for linear system.
void Prepare(Reservoir &rs, const OCP_DBL &dt)
Prepare for Assembling matrix.
void FinishStep(Reservoir &rs, OCPControl &ctrl)
Finish a time step.
void PassFlashValue(Bulk &bk, const OCP_USI &n) const
Pass value needed for FIM from flash to bulk.
void SolveLinearSystem(LinearSystem &ls, Reservoir &rs, OCPControl &ctrl) const
Solve the linear system.
OCP_BOOL UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
void CalRes(Reservoir &rs, const OCP_DBL &dt, const OCP_BOOL &resetRes0) const
Calculate residual.
void CalKrPc(Bulk &bk) const
Calculate relative permeability and capillary pressure needed for FIM.
void ResetToLastTimeStep(Reservoir &rs, OCPControl &ctrl)
Reset variables to last time step.
void AllocateReservoir(Reservoir &rs)
Allocate memory for reservoir.
void AssembleMatWells(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble linear system for wells.
void AssembleMatBulksNew(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble linear system for bulks.
void FinishStep(Reservoir &rs, OCPControl &ctrl) const
Finish a time step.
void Prepare(Reservoir &rs, const OCP_DBL &dt)
Prepare for Assembling matrix.
void AssembleMatProdWellsNew(LinearSystem &ls, const Bulk &bk, const Well &wl, const OCP_DBL &dt) const
Assemble linear system for production wells.
void GetSolution(Reservoir &rs, const vector< OCP_DBL > &u, const OCPControl &ctrl) const
Update P, Ni, BHP after linear system is solved.
void AllocateReservoir(Reservoir &rs)
Allocate memory for reservoir.
void ResetToLastTimeStep(Reservoir &rs, OCPControl &ctrl)
Reset variables to last time step.
OCP_BOOL UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
void AssembleMatInjWellsNew(LinearSystem &ls, const Bulk &bk, const Well &wl, const OCP_DBL &dt) const
Assemble linear system for injection wells.
void AssembleMatWellsNew(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble linear system for wells.
void InitFlash(Bulk &bk) const
Perform Flash with Sj and calculate values needed for FIMn.
void PassFlashValue(Bulk &bk, const OCP_USI &n) const
Pass value needed for FIMn from flash to bulk.
void AssembleMat(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble Matrix.
void Setup(Reservoir &rs, LinearSystem &ls, const OCPControl &ctrl)
Setup FIM.
void InitReservoir(Reservoir &rs) const
Init.
void CalFlash(Bulk &bk)
Perform Flash with Ni and calculate values needed for FIMn.
void SolveLinearSystem(LinearSystem &ls, Reservoir &rs, OCPControl &ctrl) const
Solve the linear system.
OCP_BOOL FinishNR(Reservoir &rs, OCPControl &ctrl)
Finish a Newton-Raphson iteration.
void UpdateLastTimeStep(Reservoir &rs) const
Update values of last step for FIMn.
IsoT_IMPEC is IMPEC (implicit pressure explict saturation) method.
void Prepare(Reservoir &rs, OCPControl &ctrl)
Prepare for Assembling matrix.
void SolveLinearSystem(LinearSystem &ls, Reservoir &rs, OCPControl &ctrl)
Solve the linear system.
void AssembleMat(LinearSystem &ls, const Reservoir &rs, const OCP_DBL &dt) const
Assemble Matrix.
void InitFlash(Bulk &bk) const
Perform Flash with Sj and calculate values needed for FIM.
void Setup(Reservoir &rs, LinearSystem &ls, const OCPControl &ctrl)
Setup IMPEC.
OCP_BOOL UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
void PassFlashValue(Bulk &bk, const OCP_USI &n) const
Pass value needed for FIM from flash to bulk.
void CalKrPc(Bulk &bk) const
Calculate relative permeability and capillary pressure needed for FIM.
OCP_BOOL FinishNR(const Reservoir &rs)
Determine if NR iteration finishes.
void InitReservoir(Reservoir &rs) const
Init.
Linear solvers for discrete systems.
All control parameters except for well controllers.
Definition: OCPControl.hpp:94
Definition: Well.hpp:45