TransportSolverTwophaseCompressiblePolymer.hpp
Go to the documentation of this file.
1/*
2 Copyright 2012 SINTEF ICT, Applied Mathematics.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_TRANSPORTSOLVERTWOPHASECOMPRESSIBLEPOLYMER_HEADER_INCLUDED
21#define OPM_TRANSPORTSOLVERTWOPHASECOMPRESSIBLEPOLYMER_HEADER_INCLUDED
22
24#include <opm/core/transport/reorder/ReorderSolverInterface.hpp>
25#include <opm/core/utility/linearInterpolation.hpp>
26#include <vector>
27#include <list>
28
29struct UnstructuredGrid;
30
31namespace {
32 class ResSOnCurve;
33 class ResCOnCurve;
34}
35
36namespace Opm
37{
38
39 class BlackoilPropertiesInterface;
40
44 class TransportSolverTwophaseCompressiblePolymer : public ReorderSolverInterface
45 {
46 public:
47
49 enum GradientMethod { Analytic, FinDif }; // Analytic is chosen (hard-coded)
50
62 TransportSolverTwophaseCompressiblePolymer(const UnstructuredGrid& grid,
63 const BlackoilPropertiesInterface& props,
64 const PolymerProperties& polyprops,
65 const SingleCellMethod method,
66 const double tol,
67 const int maxit);
68
71
92 void solve(const double* darcyflux,
93 const std::vector<double>& initial_pressure,
94 const std::vector<double>& pressure,
95 const std::vector<double>& temperature,
96 const double* porevolume0,
97 const double* porevolume,
98 const double* source,
99 const double* polymer_inflow_c,
100 const double dt,
101 std::vector<double>& saturation,
102 std::vector<double>& surfacevol,
103 std::vector<double>& concentration,
104 std::vector<double>& cmax);
105
108 void initGravity(const double* grav);
109
121 void solveGravity(const std::vector<std::vector<int> >& columns,
122 const double dt,
123 std::vector<double>& saturation,
124 std::vector<double>& surfacevol,
125 std::vector<double>& concentration,
126 std::vector<double>& cmax);
127
128
129
130
131
132 private:
133
134 const UnstructuredGrid& grid_;
135 const BlackoilPropertiesInterface& props_;
136 const PolymerProperties& polyprops_;
137 const double* darcyflux_; // one flux per grid face
138 const double* porevolume0_; // one volume per cell
139 const double* porevolume_; // one volume per cell
140 const double* source_; // one source per cell
141 const double* polymer_inflow_c_;
142 double dt_;
143 double tol_;
144 double maxit_;
145 SingleCellMethod method_;
146 double adhoc_safety_;
147
148 std::vector<double> saturation_; // one per cell, only water saturation!
149 std::vector<int> allcells_;
150 double* concentration_;
151 double* cmax_;
152 std::vector<double> fractionalflow_; // one per cell
153 std::vector<double> mc_; // one per cell
154 std::vector<double> visc_; // viscosity (without polymer, for given pressure)
155 std::vector<double> A_;
156 std::vector<double> A0_;
157 std::vector<double> smin_;
158 std::vector<double> smax_;
159
160 // For gravity segregation.
161 const double* gravity_;
162 std::vector<double> trans_;
163 std::vector<double> density_;
164 std::vector<double> gravflux_;
165 std::vector<double> mob_;
166 std::vector<double> cmax0_;
167
168 // For gravity segregation, column variables
169 std::vector<double> s0_;
170 std::vector<double> c0_;
171
172 // Storing the upwind and downwind graphs for experiments.
173 std::vector<int> ia_upw_;
174 std::vector<int> ja_upw_;
175 std::vector<int> ia_downw_;
176 std::vector<int> ja_downw_;
177
178 struct ResidualC;
179 struct ResidualS;
180
181 class ResidualCGrav;
182 class ResidualSGrav;
183
184 class ResidualEquation;
185 class ResSOnCurve;
186 class ResCOnCurve;
187
191
192
193 virtual void solveSingleCell(const int cell);
194 virtual void solveMultiCell(const int num_cells, const int* cells);
195 void solveSingleCellBracketing(int cell);
196 void solveSingleCellNewton(int cell, bool use_sc, bool use_explicit_step = false);
197 void solveSingleCellGradient(int cell);
198 void solveSingleCellGravity(const std::vector<int>& cells,
199 const int pos,
200 const double* gravflux);
201 int solveGravityColumn(const std::vector<int>& cells);
202
203 void initGravityDynamic();
204
205 void fracFlow(double s, double c, double cmax, int cell, double& ff) const;
206 void fracFlowWithDer(double s, double c, double cmax, int cell, double& ff,
207 double* dff_dsdc) const;
208 void fracFlowBoth(double s, double c, double cmax, int cell, double& ff,
209 double* dff_dsdc, bool if_with_der) const;
210 void computeMc(double c, double& mc) const;
211 void computeMcWithDer(double c, double& mc, double& dmc_dc) const;
212 void mobility(double s, double c, int cell, double* mob) const;
213 void scToc(const double* x, double* x_c) const;
214 #ifdef PROFILING
215 class Newton_Iter {
216 public:
217 bool res_s;
218 int cell;
219 double s;
220 double c;
221
222 Newton_Iter(bool res_s_val, int cell_val, double s_val, double c_val) {
223 res_s = res_s_val;
224 cell = cell_val;
225 s = s_val;
226 c = c_val;
227 }
228 };
229
230 std::list<Newton_Iter> res_counts;
231 #endif
232 };
233
234} // namespace Opm
235
236#endif // OPM_TRANSPORTSOLVERTWOPHASECOMPRESSIBLEPOLYMER_HEADER_INCLUDED
std::vector< double > & cmax
Definition: GravityColumnSolverPolymer_impl.hpp:78
Definition: PolymerProperties.hpp:35
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:45
friend class TransportSolverTwophaseCompressiblePolymer::ResidualEquation
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:188
void solveGravity(const std::vector< std::vector< int > > &columns, const double dt, std::vector< double > &saturation, std::vector< double > &surfacevol, std::vector< double > &concentration, std::vector< double > &cmax)
friend class TransportSolverTwophaseCompressiblePolymer::ResSOnCurve
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:189
void solve(const double *darcyflux, const std::vector< double > &initial_pressure, const std::vector< double > &pressure, const std::vector< double > &temperature, const double *porevolume0, const double *porevolume, const double *source, const double *polymer_inflow_c, const double dt, std::vector< double > &saturation, std::vector< double > &surfacevol, std::vector< double > &concentration, std::vector< double > &cmax)
SingleCellMethod
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:48
@ NewtonC
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:48
@ Bracketing
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:48
@ Newton
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:48
@ Gradient
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:48
TransportSolverTwophaseCompressiblePolymer(const UnstructuredGrid &grid, const BlackoilPropertiesInterface &props, const PolymerProperties &polyprops, const SingleCellMethod method, const double tol, const int maxit)
void setPreferredMethod(SingleCellMethod method)
Set the preferred method, Bracketing or Newton.
friend class TransportSolverTwophaseCompressiblePolymer::ResCOnCurve
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:190
GradientMethod
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:49
@ Analytic
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:49
@ FinDif
Definition: TransportSolverTwophaseCompressiblePolymer.hpp:49
Definition: CompressibleTpfaPolymer.hpp:33