OpenCAEPoro  0.2.0 Sep/22/2022
A simulator for multicomponent porous media flow
Reservoir.hpp
Go to the documentation of this file.
1 
12 #ifndef __RESERVOIR_HEADER__
13 #define __RESERVOIR_HEADER__
14 
15 // OpenCAEPoro header files
16 #include "Bulk.hpp"
17 #include "BulkConn.hpp"
18 #include "Grid.hpp"
19 #include "ParamRead.hpp"
20 #include "AllWells.hpp"
21 
30 class Reservoir
31 {
32  friend class OCPControl;
33  friend class Summary;
34  friend class CriticalInfo;
35  friend class DetailInfo;
36 
37  // temp
38  friend class OCP_IMPEC;
39  friend class OCP_FIM;
40  friend class OCP_AIMs;
41  friend class OCP_AIMt;
42  friend class OCP_AIMc;
43  friend class Solver;
44 
46  // General
48 
49 public:
52  void InputParam(ParamRead& param);
54  void Setup();
56  void ApplyControl(const USI& i);
58  void PrepareWell();
60  void CalWellFlux();
62  void CalWellTrans();
64  void CalVpore();
66  void CalKrPc();
68  void CalMaxChange();
70  void CalIPRT(const OCP_DBL& dt);
72  OCP_INT CheckP(const bool& bulkCheck = true, const bool& wellCheck = true);
74  bool CheckNi();
76  bool CheckVe(const OCP_DBL& Vlim) const;
78  OCP_USI GetBulkNum() const { return bulk.GetBulkNum(); }
80  OCP_USI GetMaxFIMBulk() const { return bulk.GetMaxFIMBulk(); }
82  USI GetWellNum() const { return allWells.GetWellNum(); }
84  USI GetComNum() const { return bulk.GetComNum(); }
85  void SetupWellBulk() { allWells.SetupWellBulk(bulk); }
86  void GetNTQT(const OCP_DBL& dt);
87 
88 private:
89  Grid grid;
90  Bulk bulk;
91  AllWells allWells;
92  BulkConn conn;
93 
95  // IMPEC
97 
98 public:
100  void AllocateAuxIMPEC();
102  void InitIMPEC();
104  OCP_DBL CalCFL(const OCP_DBL& dt);
106  void CalFLuxIMPEC();
108  void CalConnFluxIMPEC();
110  void MassConseveIMPEC(const OCP_DBL& dt);
112  void CalFlashIMPEC();
114  void UpdateLastStepIMPEC();
116  void AllocateMatIMPEC(LinearSystem& myLS) const;
118  void AssembleMatIMPEC(LinearSystem& myLS, const OCP_DBL& dt) const;
120  void GetSolutionIMPEC(const vector<OCP_DBL>& u);
122  void ResetWellIMPEC();
124  void ResetVal01IMPEC();
126  void ResetVal02IMPEC();
129  void ResetVal03IMPEC();
130 
131 private:
132 
133  OCP_DBL cfl;
134 
135 public:
137  // FIM
139 
141  void AllocateAuxFIM();
143  void InitFIM();
144  void InitFIM_n();
146  void CalFlashDerivFIM();
147  void CalFlashDerivFIM_n();
149  void CalKrPcDerivFIM();
151  void UpdateLastStepFIM();
153  void AllocateMatFIM(LinearSystem& myLS) const;
155  void AssembleMatFIM(LinearSystem& myLS, const OCP_DBL& dt) const;
156  void AssembleMatFIM_n(LinearSystem& myLS, const OCP_DBL& dt) const;
159  void GetSolutionFIM(const vector<OCP_DBL>& u, const OCP_DBL& dPmax,
160  const OCP_DBL& dSmax);
161  void GetSolutionFIM_n(const vector<OCP_DBL>& u, const OCP_DBL& dPmax,
162  const OCP_DBL& dSmax);
163  void GetSolution01FIM(const vector<OCP_DBL>& u);
165  void CalResFIM(ResFIM& resFIM, const OCP_DBL& dt);
167  void ResetFIM(const bool& flag);
169  OCP_DBL GetNRdPmax(){ return bulk.GetNRdPmax(); }
171  OCP_DBL GetNRdSmax(OCP_USI& index) { return bulk.CalNRdSmax(index); }
173  OCP_DBL GetNRdNmax() { return bulk.GetNRdNmax(); }
175  OCP_DBL GetNRdSmaxP() { return bulk.GetNRdSmaxP(); }
176  void PrintSolFIM(const string& outfile) const;
177  void ShowRes(const vector<OCP_DBL>& res) const;
178 
180  // AIMs, AIMt
182 
183 public:
185  void AllocateAuxAIMt();
187  void AllocateMatAIMt(LinearSystem& myLS) const;
189  void SetupFIMBulk(const bool& NRflag = false) { conn.SetupFIMBulk(bulk, NRflag); }
190  void AddFIMBulk() { conn.AddFIMBulk(bulk); }
191  void SetupFIMBulkBoundAIMs() { conn.SetupFIMBulkBoundAIMs(bulk); }
192 
194  void CalFlashDerivAIM(const bool& IfAIMs);
196  void CalKrPcDerivAIM(const bool& IfAIMs);
198  void CalResAIMt(ResFIM& resFIM, const OCP_DBL& dt);
200  void AssembleMatAIMt(LinearSystem& myLS, const OCP_DBL& dt) const;
203  void GetSolutionAIMt(const vector<OCP_DBL>& u, const OCP_DBL& dPmax,
204  const OCP_DBL& dSmax);
205 
207  void AllocateAuxAIMs();
209  void CalResAIMs(ResFIM& resFIM, const OCP_DBL& dt);
211  void AssembleMatAIMs(LinearSystem& myLS, vector<OCP_DBL>& res, const OCP_DBL& dt) const;
212  void GetSolutionAIMs(const vector<OCP_DBL>& u, const OCP_DBL& dPmax,
213  const OCP_DBL& dSmax);
214  void ResetValAIM();
215  OCP_DBL CalCFLAIM(const OCP_DBL& dt);
217  void UpdateLastStepAIM();
218 
220  void AllocateAuxAIMc();
222  void AssembleMatAIMc(LinearSystem& myLS, const OCP_DBL& dt) const;
224  void CalResAIMc(ResFIM& resFIM, const OCP_DBL& dt);
225  void CalFlashAIMc();
226  void CalFlashAIMc01();
227  void CalKrPcAIMc();
229  void CalFlashDerivAIMc();
231  void CalKrPcDerivAIMc();
232  void GetSolutionAIMc(const vector<OCP_DBL>& u, const OCP_DBL& dPmax,
233  const OCP_DBL& dSmax);
234  void InitAIMc();
235  void UpdatePj() { bulk.UpdatePj(); }
236 
237 };
238 
239 #endif /* end if __RESERVOIR_HEADER__ */
240 
241 /*----------------------------------------------------------------------------*/
242 /* Brief Change History of This File */
243 /*----------------------------------------------------------------------------*/
244 /* Author Date Actions */
245 /*----------------------------------------------------------------------------*/
246 /* Shizhe Li Oct/01/2021 Create file */
247 /* Chensong Zhang Oct/15/2021 Format file */
248 /*----------------------------------------------------------------------------*/
AllWells class declaration.
BulkConn class declaration.
Bulk class declaration.
Grid class declaration.
unsigned int USI
Generic unsigned integer.
Definition: OCPConst.hpp:22
double OCP_DBL
Double precision.
Definition: OCPConst.hpp:26
unsigned int OCP_USI
Long unsigned integer.
Definition: OCPConst.hpp:24
int OCP_INT
Long integer.
Definition: OCPConst.hpp:25
ParamRead class declaration.
void SetupWellBulk(Bulk &myBulk) const
Setup bulks which are penetrated by wells.
Definition: AllWells.cpp:139
USI GetWellNum() const
Return the num of wells.
Definition: AllWells.hpp:129
Properties and operations on connections between bulks (active grids).
Definition: BulkConn.hpp:58
Physical information of each active reservoir bulk.
Definition: Bulk.hpp:58
OCP_DBL GetNRdPmax()
Return NRdPmax.
Definition: Bulk.cpp:3061
OCP_DBL GetNRdSmaxP()
Return NRdSmaxP.
Definition: Bulk.cpp:3063
OCP_DBL CalNRdSmax(OCP_USI &index)
Calculate some auxiliary variable, for example, dSmax.
Definition: Bulk.cpp:3047
OCP_USI GetBulkNum() const
Return the number of bulks.
Definition: Bulk.hpp:145
USI GetComNum() const
Return the number of components.
Definition: Bulk.hpp:147
Collect important information of each time step for fast review.
Definition: OCPOutput.hpp:142
Collect more detailed information of each time step.
Definition: OCPOutput.hpp:165
Basic information of computational grid, including the rock properties.
Definition: Grid.hpp:76
Linear solvers for discrete systems.
All control parameters except for well controlers.
Definition: OCPControl.hpp:94
perform AIM in space, that is, some grids will be implicit, others will be explicit at the same time ...
perform AIM in time, that is, local FIM will be performed after global IMPEC performs
OCP_FIM is FIM (Fully Implicit Method).
OCP_IMPEC is IMPEC (implict pressure explict saturation) method.
Pre-processing unit for OpenCAEPoro for reading params from input files.
Definition: ParamRead.hpp:33
bool CheckVe(const OCP_DBL &Vlim) const
Check error between Fluids and Pores.
Definition: Reservoir.cpp:124
void AllocateMatIMPEC(LinearSystem &myLS) const
Allocate Maxmimum memory for internal Matirx for IMPEC.
Definition: Reservoir.cpp:222
void InitFIM()
Initialize the properties of Reservoir for FIM.
Definition: Reservoir.cpp:308
void ResetWellIMPEC()
Reset Well for IMPEC.
Definition: Reservoir.cpp:249
void AssembleMatAIMs(LinearSystem &myLS, vector< OCP_DBL > &res, const OCP_DBL &dt) const
Assemble Matrix for AIMs.
Definition: Reservoir.cpp:601
void AssembleMatFIM(LinearSystem &myLS, const OCP_DBL &dt) const
Assemble Matrix for FIM.
Definition: Reservoir.cpp:384
void AllocateAuxAIMt()
Allocate memory for auxiliary variables used for AIMt.
Definition: Reservoir.cpp:517
void UpdateLastStepFIM()
Update value of last step for FIM.
Definition: Reservoir.cpp:365
void AllocateAuxFIM()
Allocate memory for auxiliary variables used for FIM.
Definition: Reservoir.cpp:300
void ResetVal01IMPEC()
Reset Capillary Pressure, Flux for IMPEC.
Definition: Reservoir.cpp:259
void CalKrPc()
Calculate Relative Permeability and Capillary for each Bulk.
Definition: Reservoir.cpp:74
void Setup()
Setup static information for reservoir with input params.
Definition: Reservoir.cpp:28
void UpdateLastStepAIM()
Update value of last step for IMPEC.
Definition: Reservoir.cpp:640
void CalFlashIMPEC()
Calculate Flash For IMPEC.
Definition: Reservoir.cpp:205
void MassConseveIMPEC(const OCP_DBL &dt)
Calculate Ni according to Flux.
Definition: Reservoir.cpp:197
void GetSolutionFIM(const vector< OCP_DBL > &u, const OCP_DBL &dPmax, const OCP_DBL &dSmax)
Definition: Reservoir.cpp:409
void CalFlashDerivAIM(const bool &IfAIMs)
Calculate Flash for local FIM, some derivatives are needed.
Definition: Reservoir.cpp:537
OCP_DBL GetNRdSmax(OCP_USI &index)
Return NRdSmax.
Definition: Reservoir.hpp:171
void GetSolutionIMPEC(const vector< OCP_DBL > &u)
Return the Solution to Reservoir Pressure for IMPEC.
Definition: Reservoir.cpp:241
bool CheckNi()
Check if abnormal Pressure occurs.
Definition: Reservoir.cpp:117
OCP_DBL GetNRdSmaxP()
Return NRdSmaxP.
Definition: Reservoir.hpp:175
void InputParam(ParamRead &param)
Definition: Reservoir.cpp:18
void CalVpore()
Calculate pore of Bulks.
Definition: Reservoir.cpp:67
void PrepareWell()
Calculate Well Properties at the beginning of each time step.
Definition: Reservoir.cpp:46
void AllocateAuxIMPEC()
Allocate memory for auxiliary variables used for IMPEC.
Definition: Reservoir.cpp:142
void InitIMPEC()
Initialize the properties of Reservoir for IMPEC.
Definition: Reservoir.cpp:150
OCP_USI GetMaxFIMBulk() const
Return MaxNUMFIMBulk.
Definition: Reservoir.hpp:80
void ResetVal03IMPEC()
Definition: Reservoir.cpp:275
USI GetWellNum() const
Return the num of Well.
Definition: Reservoir.hpp:82
void CalKrPcDerivAIM(const bool &IfAIMs)
Calculate Relative Permeability and Capillary and some derivatives for each Bulk.
Definition: Reservoir.cpp:543
void GetSolutionAIMt(const vector< OCP_DBL > &u, const OCP_DBL &dPmax, const OCP_DBL &dSmax)
Definition: Reservoir.cpp:570
void UpdateLastStepIMPEC()
Update value of last step for IMPEC.
Definition: Reservoir.cpp:212
OCP_DBL GetNRdPmax()
Return NRdPmax.
Definition: Reservoir.hpp:169
void CalMaxChange()
Calculate Maximum Change of some reference variables for IMPEC.
Definition: Reservoir.cpp:81
void AllocateMatFIM(LinearSystem &myLS) const
Allocate Maxmimum memory for internal Matirx for FIM.
Definition: Reservoir.cpp:373
OCP_INT CheckP(const bool &bulkCheck=true, const bool &wellCheck=true)
Check if abnormal Pressure occurs.
Definition: Reservoir.cpp:98
void AssembleMatAIMc(LinearSystem &myLS, const OCP_DBL &dt) const
Assemble Matrix for AIMc.
Definition: Reservoir.cpp:659
void ResetFIM(const bool &flag)
Reset FIM.
Definition: Reservoir.cpp:471
USI GetComNum() const
Return the num of Components.
Definition: Reservoir.hpp:84
void CalFlashDerivAIMc()
Calculate Flash for local FIM, some derivatives are needed.
Definition: Reservoir.cpp:699
void CalConnFluxIMPEC()
Calculate flux between bulks.
Definition: Reservoir.cpp:190
void ApplyControl(const USI &i)
Apply the control of ith critical time point.
Definition: Reservoir.cpp:38
void CalWellTrans()
Calculate Trans of Wells.
Definition: Reservoir.cpp:60
void AllocateMatAIMt(LinearSystem &myLS) const
Allocate Maxmimum memory for internal Matirx for local FIM.
Definition: Reservoir.cpp:527
void AllocateAuxAIMc()
Allocate memory for auxiliary variables used for FIM.
Definition: Reservoir.cpp:650
void AssembleMatIMPEC(LinearSystem &myLS, const OCP_DBL &dt) const
Assemble Matrix for IMPEC.
Definition: Reservoir.cpp:232
void AssembleMatAIMt(LinearSystem &myLS, const OCP_DBL &dt) const
Assemble Matrix for AIMt -— local FIM here.
Definition: Reservoir.cpp:563
void SetupFIMBulk(const bool &NRflag=false)
Setup FIMBulk.
Definition: Reservoir.hpp:189
void AllocateAuxAIMs()
Allocate memory for auxiliary variables used for AIMs.
Definition: Reservoir.cpp:578
void CalIPRT(const OCP_DBL &dt)
Calculate num of Injection, Production.
Definition: Reservoir.cpp:89
OCP_DBL CalCFL(const OCP_DBL &dt)
Calcluate the CFL number, including bulks and wells for IMPEC.
Definition: Reservoir.cpp:170
OCP_DBL GetNRdNmax()
Return NRdNmax.
Definition: Reservoir.hpp:173
void ResetVal02IMPEC()
Reset Capillary Pressure, Moles of Componnets, Flux for IMPEC.
Definition: Reservoir.cpp:266
void CalKrPcDerivAIMc()
Calculate Relative Permeability and Capillary and some derivatives for each Bulk.
Definition: Reservoir.cpp:706
void CalKrPcDerivFIM()
Calculate Relative Permeability and Capillary and some derivatives for each Bulk.
Definition: Reservoir.cpp:358
void CalFlashDerivFIM()
Calculate Flash for FIM, some derivatives are needed.
Definition: Reservoir.cpp:344
void CalFLuxIMPEC()
Calculate flux between bulks, bulks and wells.
Definition: Reservoir.cpp:182
OCP_USI GetBulkNum() const
Return the num of Bulk.
Definition: Reservoir.hpp:78
void CalResAIMs(ResFIM &resFIM, const OCP_DBL &dt)
Calculate the Resiual for AIMs, it's also RHS of Linear System.
Definition: Reservoir.cpp:588
void CalResAIMc(ResFIM &resFIM, const OCP_DBL &dt)
Calculate the Resiual for FIM, it's also RHS of Linear System.
Definition: Reservoir.cpp:668
void CalWellFlux()
Calculate Flux between Bulk and Wells.
Definition: Reservoir.cpp:53
void CalResFIM(ResFIM &resFIM, const OCP_DBL &dt)
Calculate the Resiual for FIM, it's also RHS of Linear System.
Definition: Reservoir.cpp:434
void CalResAIMt(ResFIM &resFIM, const OCP_DBL &dt)
Calculate the Resiual for local FIM, it's also RHS of Linear System.
Definition: Reservoir.cpp:549
Solver class for overall solution methods.
Definition: Solver.hpp:21
The Summary class manages the output in the summary file.
Definition: OCPOutput.hpp:92