OpenCAEPoro  0.2.0 Sep/22/2022
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 "UtilTiming.hpp"
20 #include "UtilOutput.hpp"
21 
23 class OCP_IMPEC
24 {
25 public:
27  void Setup(Reservoir& rs, LinearSystem& myLS, const OCPControl& ctrl);
28 
30  void InitReservoir(Reservoir& rs) const;
31 
33  void Prepare(Reservoir& rs, OCP_DBL& dt);
34 
36  void SolveLinearSystem(LinearSystem& myLS, Reservoir& rs, OCPControl& ctrl);
37 
39  bool UpdateProperty(Reservoir& rs, OCPControl& ctrl);
40  bool UpdateProperty01(Reservoir& rs, OCPControl& ctrl);
41 
43  bool FinishNR(const Reservoir& rs);
44  bool FinishNR01(Reservoir& rs, OCPControl& ctrl);
45 
46  void FinishStep(Reservoir& rs, OCPControl& ctrl);
47 
48 };
49 
51 class OCP_FIM
52 {
53 public:
55  void Setup(Reservoir& rs, LinearSystem& myLS, const OCPControl& ctrl);
56 
58  void InitReservoir(Reservoir& rs) const;
59 
61  void Prepare(Reservoir& rs, OCP_DBL& dt);
62 
64  void AssembleMat(LinearSystem& myLS, const Reservoir& rs, const OCP_DBL& dt) const;
65 
67  void SolveLinearSystem(LinearSystem& myLS, Reservoir& rs, OCPControl& ctrl) const;
68 
70  bool UpdateProperty(Reservoir& rs, OCPControl& ctrl);
71 
73  bool FinishNR(Reservoir& rs, OCPControl& ctrl);
74 
76  void FinishStep(Reservoir& rs, OCPControl& ctrl) const;
77 
78 
79 protected:
82 };
83 
84 
85 class OCP_FIMn : public OCP_FIM
86 {
87 public:
88 
90  void InitReservoir(Reservoir& rs) const;
91 
93  void AssembleMat(LinearSystem& myLS, const Reservoir& rs, const OCP_DBL& dt) const;
94 
96  void SolveLinearSystem(LinearSystem& myLS, Reservoir& rs, OCPControl& ctrl) const;
97 
99  bool UpdateProperty(Reservoir& rs, OCPControl& ctrl);
100 
101 };
102 
103 class OCP_AIMc : public OCP_FIM
104 {
105 public:
107  void Setup(Reservoir& rs, LinearSystem& myLS, const OCPControl& ctrl);
108 
110  void InitReservoir(Reservoir& rs) const;
111 
113  void Prepare(Reservoir& rs, OCP_DBL& dt);
114 
116  void AssembleMat(LinearSystem& myLS, const Reservoir& rs, const OCP_DBL& dt) const;
117 
119  void SolveLinearSystem(LinearSystem& myLS, Reservoir& rs, OCPControl& ctrl);
120 
122  bool UpdateProperty(Reservoir& rs, OCPControl& ctrl);
123 
125  bool FinishNR(Reservoir& rs, OCPControl& ctrl);
126 };
127 
129 class OCP_AIMs
130 {
131 public:
133  void Setup(Reservoir& rs, LinearSystem& myLS, const OCPControl& ctrl);
135  void Prepare(Reservoir& rs, OCP_DBL& dt);
137  void AssembleMat(LinearSystem& myLS, const Reservoir& rs, const OCP_DBL& dt);
139  void SolveLinearSystem(LinearSystem& myLS, Reservoir& rs, OCPControl& ctrl);
141  bool UpdateProperty(Reservoir& rs, OCPControl& ctrl);
143  bool FinishNR(Reservoir& rs, OCPControl& ctrl);
145  void FinishStep(Reservoir& rs, OCPControl& ctrl);
146 
147 private:
149  ResFIM resFIM;
150 
151 };
152 
154 class OCP_AIMt
155 {
156 public:
158  void Setup(Reservoir& rs, LinearSystem& myLS, LinearSystem& myAuxLS, const OCPControl& ctrl);
160  void Prepare(Reservoir& rs, OCP_DBL& dt);
162  bool UpdateProperty(Reservoir& rs, OCPControl& ctrl, LinearSystem& myAuxLS);
163 
164 private:
166  ResFIM resFIM;
167 
168 };
169 
170 
171 #endif /* end if __OCPFLUIDMETHOD_HEADER__ */
172 
173 /*----------------------------------------------------------------------------*/
174 /* Brief Change History of This File */
175 /*----------------------------------------------------------------------------*/
176 /* Author Date Actions */
177 /*----------------------------------------------------------------------------*/
178 /* Shizhe Li Oct/01/2021 Create file */
179 /* Chensong Zhang Oct/15/2021 Format file */
180 /*----------------------------------------------------------------------------*/
Linear solver class declaration.
double OCP_DBL
Double precision.
Definition: OCPConst.hpp:26
OCPControl class declaration.
Reservoir class declaration.
Supply basic tools used to output files.
Elapsed wall-time and CPU-cycles declaration.
Linear solvers for discrete systems.
All control parameters except for well controlers.
Definition: OCPControl.hpp:94
bool UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
void InitReservoir(Reservoir &rs) const
Init.
void SolveLinearSystem(LinearSystem &myLS, Reservoir &rs, OCPControl &ctrl)
Solve the linear system.
void Prepare(Reservoir &rs, OCP_DBL &dt)
Prepare for Assembling matrix.
bool FinishNR(Reservoir &rs, OCPControl &ctrl)
Finish a Newton-Raphson iteration.
void AssembleMat(LinearSystem &myLS, const Reservoir &rs, const OCP_DBL &dt) const
Assemble Matrix.
void Setup(Reservoir &rs, LinearSystem &myLS, const OCPControl &ctrl)
Setup AIMc.
perform AIM in space, that is, some grids will be implicit, others will be explicit at the same time ...
void AssembleMat(LinearSystem &myLS, const Reservoir &rs, const OCP_DBL &dt)
Assemble Matrix.
void Setup(Reservoir &rs, LinearSystem &myLS, const OCPControl &ctrl)
Setup AIMs.
void SolveLinearSystem(LinearSystem &myLS, Reservoir &rs, OCPControl &ctrl)
Solve the linear system.
void Prepare(Reservoir &rs, OCP_DBL &dt)
Prepare for Assembling matrix.
void FinishStep(Reservoir &rs, OCPControl &ctrl)
Finish a time step.
bool UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
bool FinishNR(Reservoir &rs, OCPControl &ctrl)
Finish a Newton-Raphson iteration.
perform AIM in time, that is, local FIM will be performed after global IMPEC performs
void Prepare(Reservoir &rs, OCP_DBL &dt)
Prepare for Assembling matrix.
bool UpdateProperty(Reservoir &rs, OCPControl &ctrl, LinearSystem &myAuxLS)
Update properties of fluids.
void Setup(Reservoir &rs, LinearSystem &myLS, LinearSystem &myAuxLS, const OCPControl &ctrl)
Setup AIMt.
OCP_FIM is FIM (Fully Implicit Method).
void AssembleMat(LinearSystem &myLS, const Reservoir &rs, const OCP_DBL &dt) const
Assemble Matrix.
void Prepare(Reservoir &rs, OCP_DBL &dt)
Prepare for Assembling matrix.
void SolveLinearSystem(LinearSystem &myLS, Reservoir &rs, OCPControl &ctrl) const
Solve the linear system.
void Setup(Reservoir &rs, LinearSystem &myLS, const OCPControl &ctrl)
Setup FIM.
void FinishStep(Reservoir &rs, OCPControl &ctrl) const
Finish a time step.
ResFIM resFIM
Resiual for FIM.
bool UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
void InitReservoir(Reservoir &rs) const
Init.
bool FinishNR(Reservoir &rs, OCPControl &ctrl)
Finish a Newton-Raphson iteration.
void SolveLinearSystem(LinearSystem &myLS, Reservoir &rs, OCPControl &ctrl) const
Solve the linear system.
bool UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
void InitReservoir(Reservoir &rs) const
Init.
void AssembleMat(LinearSystem &myLS, const Reservoir &rs, const OCP_DBL &dt) const
Assemble Matrix.
OCP_IMPEC is IMPEC (implict pressure explict saturation) method.
void InitReservoir(Reservoir &rs) const
Init.
void Setup(Reservoir &rs, LinearSystem &myLS, const OCPControl &ctrl)
Setup IMPEC.
void SolveLinearSystem(LinearSystem &myLS, Reservoir &rs, OCPControl &ctrl)
Solve the linear system.
void Prepare(Reservoir &rs, OCP_DBL &dt)
Prepare for Assembling matrix.
bool UpdateProperty(Reservoir &rs, OCPControl &ctrl)
Update properties of fluids.
bool FinishNR(const Reservoir &rs)
Determine if NR iteration finishes.