FullyImplicitCompressiblePolymerSolver.hpp
Go to the documentation of this file.
1/*
2 Copyright 2014 SINTEF ICT, Applied Mathematics.
3 Copyright 2014 STATOIL ASA.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21#ifndef OPM_FULLYIMPLICITCOMPRESSIBLEPOLYMERSOLVER_HEADER_INCLUDED
22#define OPM_FULLYIMPLICITCOMPRESSIBLEPOLYMERSOLVER_HEADER_INCLUDED
23
24#include <opm/autodiff/AutoDiffBlock.hpp>
25#include <opm/autodiff/AutoDiffHelpers.hpp>
26#include <opm/autodiff/BlackoilPropsAdInterface.hpp>
27#include <opm/autodiff/NewtonIterationBlackoilInterface.hpp>
28#include <opm/autodiff/LinearisedBlackoilResidual.hpp>
32#include <opm/core/utility/parameters/ParameterGroup.hpp>
33
34struct UnstructuredGrid;
35struct Wells;
36
37namespace Opm {
38
39 class DerivedGeology;
40 class RockCompressibility;
41 class NewtonIterationBlackoilInterface;
42 class PolymerBlackoilState;
43 class WellStateFullyImplicitBlackoil;
44
54 {
55 public:
66 FullyImplicitCompressiblePolymerSolver(const UnstructuredGrid& grid ,
67 const BlackoilPropsAdInterface& fluid,
68 const DerivedGeology& geo ,
69 const RockCompressibility* rock_comp_props,
70 const PolymerPropsAd& polymer_props_ad,
71 const Wells& wells,
72 const NewtonIterationBlackoilInterface& linsolver);
73
84 int
85 step(const double dt,
88
90 int linearIterations() const;
91
93 typedef parameter::ParameterGroup SolverParameters;
94
97
99 double relativeChange(const PolymerBlackoilState& previous,
100 const PolymerBlackoilState& current ) const;
101
102 private:
103 typedef AutoDiffBlock<double> ADB;
104 typedef ADB::V V;
105 typedef ADB::M M;
106 typedef Eigen::Array<double,
107 Eigen::Dynamic,
108 Eigen::Dynamic,
109 Eigen::RowMajor> DataBlock;
110
111 struct ReservoirResidualQuant {
112 ReservoirResidualQuant();
113 std::vector<ADB> accum; // Accumulations
114 ADB mflux; // Mass flux (surface conditions)
115 ADB b; // Reciprocal FVF
116 ADB head; // Pressure drop across int. interfaces
117 ADB mob; // Phase mobility (per cell)
118 std::vector<ADB> ads; // Adsorption term.
119 };
120
121 struct SolutionState {
122 SolutionState(const int np);
123 ADB pressure;
124 ADB temperature;
125 std::vector<ADB> saturation;
126 ADB concentration;
127 ADB qs;
128 ADB bhp;
129 };
130
131 struct WellOps {
132 WellOps(const Wells& wells);
133 Eigen::SparseMatrix<double> w2p; // well -> perf (scatter)
134 Eigen::SparseMatrix<double> p2w; // perf -> well (gather)
135 };
136
137 enum { Water = BlackoilPropsAdInterface::Water,
138 Oil = BlackoilPropsAdInterface::Oil };
139
140 // Member data
141 const UnstructuredGrid& grid_;
142 const BlackoilPropsAdInterface& fluid_;
143 const DerivedGeology& geo_;
144 const RockCompressibility* rock_comp_props_;
145 const PolymerPropsAd& polymer_props_ad_;
146 const Wells& wells_;
147 const NewtonIterationBlackoilInterface& linsolver_;
148 const std::vector<int> cells_; // All grid cells
149 HelperOps ops_;
150 const WellOps wops_;
151 const M grav_;
152 V cmax_;
153 std::vector<PhasePresence> phaseCondition_;
154 std::vector<ReservoirResidualQuant> rq_;
155 // The mass_balance vector has one element for each active phase,
156 // each of which has size equal to the number of cells.
157 // The well_eq has size equal to the number of wells.
158 LinearisedBlackoilResidual residual_;
159
160 unsigned int newtonIterations_;
161 unsigned int linearIterations_;
162
163 // Private methods.
164 SolutionState
165 constantState(const PolymerBlackoilState& x,
166 const WellStateFullyImplicitBlackoil& xw);
167
168 SolutionState
169 variableState(const PolymerBlackoilState& x,
170 const WellStateFullyImplicitBlackoil& xw);
171
172 void
173 computeAccum(const SolutionState& state,
174 const int aix );
175
176 void
177 assemble(const double dt,
178 const PolymerBlackoilState& x,
179 const WellStateFullyImplicitBlackoil& xw,
180 const std::vector<double>& polymer_inflow);
181
182 V solveJacobianSystem() const;
183
184 void updateState(const V& dx,
185 PolymerBlackoilState& state,
186 WellStateFullyImplicitBlackoil& well_state) const;
187
188 std::vector<ADB>
189 computeRelPerm(const SolutionState& state) const;
190
191 std::vector<ADB>
192 computePressures(const SolutionState& state) const;
193
194 void
195 computeMassFlux(const int actph ,
196 const V& transi,
197 const std::vector<ADB>& kr ,
198 const SolutionState& state );
199
200 void
201 computeMassFlux(const V& trans,
202 const ADB& mc,
203 const ADB& kro,
204 const ADB& krw_eff,
205 const SolutionState& state);
206
207 std::vector<ADB>
208 computeFracFlow(const ADB& kro,
209 const ADB& krw_eff,
210 const ADB& c) const;
211
212 void
213 computeCmax(PolymerBlackoilState& state);
214
215 ADB
216 computeMc(const SolutionState& state) const;
217
218 ADB
219 rockPorosity(const ADB& p) const;
220
221 ADB
222 rockPermeability(const ADB& p) const;
223
224 double
225 residualNorm() const;
226
227 ADB
228 fluidViscosity(const int phase,
229 const ADB& p ,
230 const ADB& T ,
231 const std::vector<PhasePresence>& cond,
232 const std::vector<int>& cells) const;
233
234 ADB
235 fluidReciprocFVF(const int phase,
236 const ADB& p ,
237 const ADB& T ,
238 const std::vector<PhasePresence>& cond,
239 const std::vector<int>& cells) const;
240
241 ADB
242 fluidDensity(const int phase,
243 const ADB& p ,
244 const ADB& T ,
245 const std::vector<PhasePresence>& cond,
246 const std::vector<int>& cells) const;
247
248 ADB
249 poroMult(const ADB& p) const;
250
251 ADB
252 transMult(const ADB& p) const;
253
254 const std::vector<PhasePresence>
255 phaseCondition() const { return phaseCondition_; }
256
257 void
258 classifyCondition(const PolymerBlackoilState& state);
259 };
260} // namespace Opm
261
262
263#endif // OPM_FULLYIMPLICITCOMPRESSIBLEPOLYMERSOLVER_HEADER_INCLUDED
Definition: FullyImplicitCompressiblePolymerSolver.hpp:54
FullyImplicitCompressiblePolymerSolver(const UnstructuredGrid &grid, const BlackoilPropsAdInterface &fluid, const DerivedGeology &geo, const RockCompressibility *rock_comp_props, const PolymerPropsAd &polymer_props_ad, const Wells &wells, const NewtonIterationBlackoilInterface &linsolver)
double relativeChange(const PolymerBlackoilState &previous, const PolymerBlackoilState &current) const
Evaluate the relative changes in the physical variables.
int step(const double dt, PolymerBlackoilState &state, WellStateFullyImplicitBlackoilPolymer &wstate)
parameter::ParameterGroup SolverParameters
Not used by this class except to satisfy interface requirements.
Definition: FullyImplicitCompressiblePolymerSolver.hpp:93
const FullyImplicitCompressiblePolymerSolver & model() const
There is no separate model class for this solver, return itself.
Definition: PolymerBlackoilState.hpp:34
Definition: PolymerPropsAd.hpp:33
Definition: WellStateFullyImplicitBlackoilPolymer.hpp:29
Definition: CompressibleTpfaPolymer.hpp:33