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 
34 struct UnstructuredGrid;
35 struct Wells;
36 
37 namespace 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,
86  PolymerBlackoilState& state ,
88 
89  int newtonIterations() const;
90  int linearIterations() const;
91 
93  typedef parameter::ParameterGroup SolverParameters;
94 
95  private:
96  typedef AutoDiffBlock<double> ADB;
97  typedef ADB::V V;
98  typedef ADB::M M;
99  typedef Eigen::Array<double,
100  Eigen::Dynamic,
101  Eigen::Dynamic,
102  Eigen::RowMajor> DataBlock;
103 
104  struct ReservoirResidualQuant {
105  ReservoirResidualQuant();
106  std::vector<ADB> accum; // Accumulations
107  ADB mflux; // Mass flux (surface conditions)
108  ADB b; // Reciprocal FVF
109  ADB head; // Pressure drop across int. interfaces
110  ADB mob; // Phase mobility (per cell)
111  std::vector<ADB> ads; // Adsorption term.
112  };
113 
114  struct SolutionState {
115  SolutionState(const int np);
116  ADB pressure;
117  ADB temperature;
118  std::vector<ADB> saturation;
119  ADB concentration;
120  ADB qs;
121  ADB bhp;
122  };
123 
124  struct WellOps {
125  WellOps(const Wells& wells);
126  Eigen::SparseMatrix<double> w2p; // well -> perf (scatter)
127  Eigen::SparseMatrix<double> p2w; // perf -> well (gather)
128  };
129 
130  enum { Water = BlackoilPropsAdInterface::Water,
131  Oil = BlackoilPropsAdInterface::Oil };
132 
133  // Member data
134  const UnstructuredGrid& grid_;
135  const BlackoilPropsAdInterface& fluid_;
136  const DerivedGeology& geo_;
137  const RockCompressibility* rock_comp_props_;
138  const PolymerPropsAd& polymer_props_ad_;
139  const Wells& wells_;
140  const NewtonIterationBlackoilInterface& linsolver_;
141  const std::vector<int> cells_; // All grid cells
142  HelperOps ops_;
143  const WellOps wops_;
144  const M grav_;
145  V cmax_;
146  std::vector<PhasePresence> phaseCondition_;
147  std::vector<ReservoirResidualQuant> rq_;
148  // The mass_balance vector has one element for each active phase,
149  // each of which has size equal to the number of cells.
150  // The well_eq has size equal to the number of wells.
151  LinearisedBlackoilResidual residual_;
152 
153  unsigned int newtonIterations_;
154  unsigned int linearIterations_;
155 
156  // Private methods.
157  SolutionState
158  constantState(const PolymerBlackoilState& x,
159  const WellStateFullyImplicitBlackoil& xw);
160 
161  SolutionState
162  variableState(const PolymerBlackoilState& x,
163  const WellStateFullyImplicitBlackoil& xw);
164 
165  void
166  computeAccum(const SolutionState& state,
167  const int aix );
168 
169  void
170  assemble(const double dt,
171  const PolymerBlackoilState& x,
172  const WellStateFullyImplicitBlackoil& xw,
173  const std::vector<double>& polymer_inflow);
174 
175  V solveJacobianSystem() const;
176 
177  void updateState(const V& dx,
178  PolymerBlackoilState& state,
179  WellStateFullyImplicitBlackoil& well_state) const;
180 
181  std::vector<ADB>
182  computeRelPerm(const SolutionState& state) const;
183 
184  std::vector<ADB>
185  computePressures(const SolutionState& state) const;
186 
187  void
188  computeMassFlux(const int actph ,
189  const V& transi,
190  const std::vector<ADB>& kr ,
191  const SolutionState& state );
192 
193  void
194  computeMassFlux(const V& trans,
195  const ADB& mc,
196  const ADB& kro,
197  const ADB& krw_eff,
198  const SolutionState& state);
199 
200  std::vector<ADB>
201  computeFracFlow(const ADB& kro,
202  const ADB& krw_eff,
203  const ADB& c) const;
204 
205  void
206  computeCmax(PolymerBlackoilState& state);
207 
208  ADB
209  computeMc(const SolutionState& state) const;
210 
211  ADB
212  rockPorosity(const ADB& p) const;
213 
214  ADB
215  rockPermeability(const ADB& p) const;
216 
217  double
218  residualNorm() const;
219 
220  ADB
221  fluidViscosity(const int phase,
222  const ADB& p ,
223  const ADB& T ,
224  const std::vector<PhasePresence>& cond,
225  const std::vector<int>& cells) const;
226 
227  ADB
228  fluidReciprocFVF(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  fluidDensity(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  poroMult(const ADB& p) const;
243 
244  ADB
245  transMult(const ADB& p) const;
246 
247  const std::vector<PhasePresence>
248  phaseCondition() const { return phaseCondition_; }
249 
250  void
251  classifyCondition(const PolymerBlackoilState& state);
252  };
253 } // namespace Opm
254 
255 
256 #endif // OPM_FULLYIMPLICITCOMPRESSIBLEPOLYMERSOLVER_HEADER_INCLUDED
Definition: CompressibleTpfaPolymer.hpp:32
Definition: WellStateFullyImplicitBlackoilPolymer.hpp:28
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)
Definition: PolymerBlackoilState.hpp:33
Definition: PolymerPropsAd.hpp:32
Definition: FullyImplicitCompressiblePolymerSolver.hpp:53
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