BlackoilPolymerModel_impl.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
3  Copyright 2014, 2015 Dr. Blatt - HPC-Simulation-Software & Services
4  Copyright 2014, 2015 Statoil ASA.
5  Copyright 2015 NTNU
6  Copyright 2015 IRIS AS
7 
8  This file is part of the Open Porous Media project (OPM).
9 
10  OPM is free software: you can redistribute it and/or modify
11  it under the terms of the GNU General Public License as published by
12  the Free Software Foundation, either version 3 of the License, or
13  (at your option) any later version.
14 
15  OPM is distributed in the hope that it will be useful,
16  but WITHOUT ANY WARRANTY; without even the implied warranty of
17  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  GNU General Public License for more details.
19 
20  You should have received a copy of the GNU General Public License
21  along with OPM. If not, see <http://www.gnu.org/licenses/>.
22 */
23 
24 #ifndef OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
25 #define OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
26 
28 
29 #include <opm/autodiff/AutoDiffBlock.hpp>
30 #include <opm/autodiff/AutoDiffHelpers.hpp>
31 #include <opm/autodiff/GridHelpers.hpp>
32 #include <opm/autodiff/BlackoilPropsAdInterface.hpp>
33 #include <opm/autodiff/GeoProps.hpp>
34 #include <opm/autodiff/WellDensitySegmented.hpp>
35 
36 #include <opm/core/grid.h>
37 #include <opm/core/linalg/LinearSolverInterface.hpp>
38 #include <opm/core/linalg/ParallelIstlInformation.hpp>
39 #include <opm/core/props/rock/RockCompressibility.hpp>
40 #include <opm/common/ErrorMacros.hpp>
41 #include <opm/common/Exceptions.hpp>
42 #include <opm/core/utility/Units.hpp>
43 #include <opm/core/well_controls.h>
44 #include <opm/core/utility/parameters/ParameterGroup.hpp>
45 
46 #include <cassert>
47 #include <cmath>
48 #include <iostream>
49 #include <iomanip>
50 #include <limits>
51 
52 namespace Opm {
53 
54 
55 
56  namespace detail {
57 
58  template <class PU>
59  int polymerPos(const PU& pu)
60  {
61  const int maxnp = Opm::BlackoilPhases::MaxNumPhases;
62  int pos = 0;
63  for (int phase = 0; phase < maxnp; ++phase) {
64  if (pu.phase_used[phase]) {
65  pos++;
66  }
67  }
68 
69  return pos;
70  }
71 
72  } // namespace detail
73 
74 
75 
76  template <class Grid>
77  BlackoilPolymerModel<Grid>::BlackoilPolymerModel(const typename Base::ModelParameters& param,
78  const Grid& grid,
79  const BlackoilPropsAdInterface& fluid,
80  const DerivedGeology& geo,
81  const RockCompressibility* rock_comp_props,
82  const PolymerPropsAd& polymer_props_ad,
83  const Wells* wells,
84  const NewtonIterationBlackoilInterface& linsolver,
85  EclipseStateConstPtr eclipse_state,
86  const bool has_disgas,
87  const bool has_vapoil,
88  const bool has_polymer,
89  const bool has_plyshlog,
90  const bool has_shrate,
91  const std::vector<double>& wells_rep_radius,
92  const std::vector<double>& wells_perf_length,
93  const std::vector<double>& wells_bore_diameter,
94  const bool terminal_output)
95  : Base(param, grid, fluid, geo, rock_comp_props, wells, linsolver, eclipse_state,
96  has_disgas, has_vapoil, terminal_output),
97  polymer_props_ad_(polymer_props_ad),
98  has_polymer_(has_polymer),
99  has_plyshlog_(has_plyshlog),
100  has_shrate_(has_shrate),
101  poly_pos_(detail::polymerPos(fluid.phaseUsage())),
102  wells_rep_radius_(wells_rep_radius),
103  wells_perf_length_(wells_perf_length),
104  wells_bore_diameter_(wells_bore_diameter)
105  {
106  if (has_polymer_) {
107  if (!active_[Water]) {
108  OPM_THROW(std::logic_error, "Polymer must solved in water!\n");
109  }
110  // If deck has polymer, residual_ should contain polymer equation.
111  rq_.resize(fluid_.numPhases() + 1);
112  residual_.material_balance_eq.resize(fluid_.numPhases() + 1, ADB::null());
113  Base::material_name_.push_back("Polymer");
114  assert(poly_pos_ == fluid_.numPhases());
115  }
116  }
117 
118 
119 
120  template <class Grid>
121  void
123  prepareStep(const double dt,
124  ReservoirState& reservoir_state,
125  WellState& well_state)
126  {
127  Base::prepareStep(dt, reservoir_state, well_state);
128  // Initial max concentration of this time step from PolymerBlackoilState.
129  cmax_ = Eigen::Map<const V>(reservoir_state.maxconcentration().data(), Opm::AutoDiffGrid::numCells(grid_));
130  }
131 
132 
133 
134 
135  template <class Grid>
136  void
138  afterStep(const double /* dt */,
139  ReservoirState& reservoir_state,
140  WellState& /* well_state */)
141  {
142  computeCmax(reservoir_state);
143  }
144 
145 
146 
147 
148 
149  template <class Grid>
150  void
152  {
153  Base::makeConstantState(state);
154  state.concentration = ADB::constant(state.concentration.value());
155  }
156 
157 
158 
159 
160 
161  template <class Grid>
162  std::vector<V>
164  const WellState& xw) const
165  {
166  std::vector<V> vars0 = Base::variableStateInitials(x, xw);
167  assert(int(vars0.size()) == fluid_.numPhases() + 2);
168 
169  // Initial polymer concentration.
170  if (has_polymer_) {
171  assert (not x.concentration().empty());
172  const int nc = x.concentration().size();
173  const V c = Eigen::Map<const V>(&x.concentration()[0], nc);
174  // Concentration belongs after other reservoir vars but before well vars.
175  auto concentration_pos = vars0.begin() + fluid_.numPhases();
176  assert(concentration_pos == vars0.end() - 2);
177  vars0.insert(concentration_pos, c);
178  }
179  return vars0;
180  }
181 
182 
183 
184 
185 
186  template <class Grid>
187  std::vector<int>
189  {
190  std::vector<int> ind = Base::variableStateIndices();
191  assert(ind.size() == 5);
192  if (has_polymer_) {
193  ind.resize(6);
194  // Concentration belongs after other reservoir vars but before well vars.
195  ind[Concentration] = fluid_.numPhases();
196  // Concentration is pushing back the well vars.
197  ++ind[Qs];
198  ++ind[Bhp];
199  }
200  return ind;
201  }
202 
203 
204 
205 
206  template <class Grid>
209  const std::vector<int>& indices,
210  std::vector<ADB>& vars) const
211  {
212  SolutionState state = Base::variableStateExtractVars(x, indices, vars);
213  if (has_polymer_) {
214  state.concentration = std::move(vars[indices[Concentration]]);
215  }
216  return state;
217  }
218 
219 
220 
221 
222 
223  template <class Grid>
224  void
226  const int aix )
227  {
228  Base::computeAccum(state, aix);
229 
230  // Compute accumulation of polymer equation only if needed.
231  if (has_polymer_) {
232  const ADB& press = state.pressure;
233  const std::vector<ADB>& sat = state.saturation;
234  const ADB& c = state.concentration;
235  const ADB pv_mult = poroMult(press); // also computed in Base::computeAccum, could be optimized.
236  const Opm::PhaseUsage& pu = fluid_.phaseUsage();
237  // compute polymer properties.
238  const ADB cmax = ADB::constant(cmax_, state.concentration.blockPattern());
239  const ADB ads = polymer_props_ad_.adsorption(state.concentration, cmax);
240  const double rho_rock = polymer_props_ad_.rockDensity();
241  const V phi = Eigen::Map<const V>(&fluid_.porosity()[0], AutoDiffGrid::numCells(grid_));
242  const double dead_pore_vol = polymer_props_ad_.deadPoreVol();
243  // Compute polymer accumulation term.
244  rq_[poly_pos_].accum[aix] = pv_mult * rq_[pu.phase_pos[Water]].b * sat[pu.phase_pos[Water]] * c * (1. - dead_pore_vol)
245  + pv_mult * rho_rock * (1. - phi) / phi * ads;
246  }
247 
248  }
249 
250 
251 
252 
253 
254 
255  template <class Grid>
257  {
258  const int nc = AutoDiffGrid::numCells(grid_);
259  V tmp = V::Zero(nc);
260  for (int i = 0; i < nc; ++i) {
261  tmp[i] = std::max(state.maxconcentration()[i], state.concentration()[i]);
262  }
263  std::copy(&tmp[0], &tmp[0] + nc, state.maxconcentration().begin());
264  }
265 
266 
267 
268 
269 
270  template <class Grid>
271  void
274  {
275  // Base::assembleMassBalanceEq(state);
276 
277  // Compute b_p and the accumulation term b_p*s_p for each phase,
278  // except gas. For gas, we compute b_g*s_g + Rs*b_o*s_o.
279  // These quantities are stored in rq_[phase].accum[1].
280  // The corresponding accumulation terms from the start of
281  // the timestep (b^0_p*s^0_p etc.) were already computed
282  // on the initial call to assemble() and stored in rq_[phase].accum[0].
283  computeAccum(state, 1);
284 
285  // Set up the common parts of the mass balance equations
286  // for each active phase.
287  const V transi = subset(geo_.transmissibility(), ops_.internal_faces);
288  const std::vector<ADB> kr = computeRelPerm(state);
289 
290 
291  if (has_plyshlog_) {
292  std::vector<double> water_vel;
293  std::vector<double> visc_mult;
294 
295  computeWaterShearVelocityFaces(transi, kr, state.canonical_phase_pressures, state, water_vel, visc_mult);
296  if ( !polymer_props_ad_.computeShearMultLog(water_vel, visc_mult, shear_mult_faces_) ) {
297  // std::cerr << " failed in calculating the shear-multiplier " << std::endl;
298  OPM_THROW(std::runtime_error, " failed in calculating the shear-multiplier. ");
299  }
300  }
301 
302  for (int phaseIdx = 0; phaseIdx < fluid_.numPhases(); ++phaseIdx) {
303  computeMassFlux(phaseIdx, transi, kr[canph_[phaseIdx]], state.canonical_phase_pressures[canph_[phaseIdx]], state);
304 
305  residual_.material_balance_eq[ phaseIdx ] =
306  pvdt_ * (rq_[phaseIdx].accum[1] - rq_[phaseIdx].accum[0])
307  + ops_.div*rq_[phaseIdx].mflux;
308  }
309 
310  // -------- Extra (optional) rs and rv contributions to the mass balance equations --------
311 
312  // Add the extra (flux) terms to the mass balance equations
313  // From gas dissolved in the oil phase (rs) and oil vaporized in the gas phase (rv)
314  // The extra terms in the accumulation part of the equation are already handled.
315  if (active_[ Oil ] && active_[ Gas ]) {
316  const int po = fluid_.phaseUsage().phase_pos[ Oil ];
317  const int pg = fluid_.phaseUsage().phase_pos[ Gas ];
318 
319  const UpwindSelector<double> upwindOil(grid_, ops_,
320  rq_[po].dh.value());
321  const ADB rs_face = upwindOil.select(state.rs);
322 
323  const UpwindSelector<double> upwindGas(grid_, ops_,
324  rq_[pg].dh.value());
325  const ADB rv_face = upwindGas.select(state.rv);
326 
327  residual_.material_balance_eq[ pg ] += ops_.div * (rs_face * rq_[po].mflux);
328  residual_.material_balance_eq[ po ] += ops_.div * (rv_face * rq_[pg].mflux);
329 
330  // OPM_AD_DUMP(residual_.material_balance_eq[ Gas ]);
331 
332  }
333 
334  // Add polymer equation.
335  if (has_polymer_) {
336  residual_.material_balance_eq[poly_pos_] = pvdt_ * (rq_[poly_pos_].accum[1] - rq_[poly_pos_].accum[0])
337  + ops_.div*rq_[poly_pos_].mflux;
338  }
339  }
340 
341 
342 
343 
344 
345  template <class Grid>
347  const SolutionState& state,
348  WellState& xw)
349 
350  {
351  Base::addWellContributionToMassBalanceEq(cq_s, state, xw);
352 
353  // Add well contributions to polymer mass balance equation
354  if (has_polymer_) {
355  const ADB mc = computeMc(state);
356  const int nc = xw.polymerInflow().size();
357  const V polyin = Eigen::Map<const V>(xw.polymerInflow().data(), nc);
358  const int nperf = wells().well_connpos[wells().number_of_wells];
359  const std::vector<int> well_cells(wells().well_cells, wells().well_cells + nperf);
360  const V poly_in_perf = subset(polyin, well_cells);
361  const V poly_mc_perf = subset(mc.value(), well_cells);
362  const ADB& cq_s_water = cq_s[fluid_.phaseUsage().phase_pos[Water]];
363  Selector<double> injector_selector(cq_s_water.value());
364  const V poly_perf = injector_selector.select(poly_in_perf, poly_mc_perf);
365  const ADB cq_s_poly = cq_s_water * poly_perf;
366  residual_.material_balance_eq[poly_pos_] -= superset(cq_s_poly, well_cells, nc);
367  }
368  }
369 
370 
371 
372 
373 
374 
375  template <class Grid>
377  ReservoirState& reservoir_state,
378  WellState& well_state)
379  {
380  if (has_polymer_) {
381  // Extract concentration change.
382  const int np = fluid_.numPhases();
383  const int nc = Opm::AutoDiffGrid::numCells(grid_);
384  const V zero = V::Zero(nc);
385  const int concentration_start = nc * np;
386  const V dc = subset(dx, Span(nc, 1, concentration_start));
387 
388  // Create new dx with the dc part deleted.
389  V modified_dx = V::Zero(dx.size() - nc);
390  modified_dx.head(concentration_start) = dx.head(concentration_start);
391  const int tail_len = dx.size() - concentration_start - nc;
392  modified_dx.tail(tail_len) = dx.tail(tail_len);
393 
394  // Call base version.
395  Base::updateState(modified_dx, reservoir_state, well_state);
396 
397  // Update concentration.
398  const V c_old = Eigen::Map<const V>(&reservoir_state.concentration()[0], nc, 1);
399  const V c = (c_old - dc).max(zero);
400  std::copy(&c[0], &c[0] + nc, reservoir_state.concentration().begin());
401  } else {
402  // Just forward call to base version.
403  Base::updateState(dx, reservoir_state, well_state);
404  }
405  }
406 
407 
408 
409 
410 
411  template <class Grid>
412  void
414  const V& transi,
415  const ADB& kr ,
416  const ADB& phasePressure,
417  const SolutionState& state)
418  {
419  Base::computeMassFlux(actph, transi, kr, phasePressure, state);
420 
421  // Polymer treatment.
422  const int canonicalPhaseIdx = canph_[ actph ];
423  if (canonicalPhaseIdx == Water) {
424  if (has_polymer_) {
425  const std::vector<PhasePresence>& cond = phaseCondition();
426  const ADB tr_mult = transMult(state.pressure);
427  const ADB mu = fluidViscosity(canonicalPhaseIdx, phasePressure, state.temperature, state.rs, state.rv, cond);
428  const ADB cmax = ADB::constant(cmax_, state.concentration.blockPattern());
429  const ADB mc = computeMc(state);
430  const ADB krw_eff = polymer_props_ad_.effectiveRelPerm(state.concentration, cmax, kr);
431  const ADB inv_wat_eff_visc = polymer_props_ad_.effectiveInvWaterVisc(state.concentration, mu.value().data());
432  // Reduce mobility of water phase by relperm reduction and effective viscosity increase.
433  rq_[actph].mob = tr_mult * krw_eff * inv_wat_eff_visc;
434  // Compute polymer mobility.
435  rq_[poly_pos_].mob = tr_mult * mc * krw_eff * inv_wat_eff_visc;
436  rq_[poly_pos_].b = rq_[actph].b;
437  rq_[poly_pos_].dh = rq_[actph].dh;
438  UpwindSelector<double> upwind(grid_, ops_, rq_[poly_pos_].dh.value());
439  // Compute polymer flux.
440  rq_[poly_pos_].mflux = upwind.select(rq_[poly_pos_].b * rq_[poly_pos_].mob) * (transi * rq_[poly_pos_].dh);
441  // Must recompute water flux since we have to use modified mobilities.
442  rq_[ actph ].mflux = upwind.select(rq_[actph].b * rq_[actph].mob) * (transi * rq_[actph].dh);
443 
444  // applying the shear-thinning factors
445  if (has_plyshlog_) {
446  V shear_mult_faces_v = Eigen::Map<V>(shear_mult_faces_.data(), shear_mult_faces_.size());
447  ADB shear_mult_faces_adb = ADB::constant(shear_mult_faces_v);
448  rq_[poly_pos_].mflux = rq_[poly_pos_].mflux / shear_mult_faces_adb;
449  rq_[actph].mflux = rq_[actph].mflux / shear_mult_faces_adb;
450  }
451  }
452  }
453  }
454 
455 
456 
457 
458  template <class Grid>
459  void
461  WellState& well_state,
462  const bool initial_assembly)
463  {
464  using namespace Opm::AutoDiffGrid;
465 
466  // Possibly switch well controls and updating well state to
467  // get reasonable initial conditions for the wells
468  updateWellControls(well_state);
469 
470  // Create the primary variables.
471  SolutionState state = variableState(reservoir_state, well_state);
472 
473  if (initial_assembly) {
474  // Create the (constant, derivativeless) initial state.
475  SolutionState state0 = state;
476  makeConstantState(state0);
477  // Compute initial accumulation contributions
478  // and well connection pressures.
479  computeAccum(state0, 0);
480  computeWellConnectionPressures(state0, well_state);
481  }
482 
483  // OPM_AD_DISKVAL(state.pressure);
484  // OPM_AD_DISKVAL(state.saturation[0]);
485  // OPM_AD_DISKVAL(state.saturation[1]);
486  // OPM_AD_DISKVAL(state.saturation[2]);
487  // OPM_AD_DISKVAL(state.rs);
488  // OPM_AD_DISKVAL(state.rv);
489  // OPM_AD_DISKVAL(state.qs);
490  // OPM_AD_DISKVAL(state.bhp);
491 
492  // -------- Mass balance equations --------
493  assembleMassBalanceEq(state);
494 
495  // -------- Well equations ----------
496  if ( ! wellsActive() ) {
497  return;
498  }
499 
500  V aliveWells;
501 
502  const int np = wells().number_of_phases;
503  std::vector<ADB> cq_s(np, ADB::null());
504 
505  const int nw = wells().number_of_wells;
506  const int nperf = wells().well_connpos[nw];
507  const std::vector<int> well_cells(wells().well_cells, wells().well_cells + nperf);
508 
509  std::vector<ADB> mob_perfcells(np, ADB::null());
510  std::vector<ADB> b_perfcells(np, ADB::null());
511  for (int phase = 0; phase < np; ++phase) {
512  mob_perfcells[phase] = subset(rq_[phase].mob, well_cells);
513  b_perfcells[phase] = subset(rq_[phase].b, well_cells);
514  }
515  if (param_.solve_welleq_initially_ && initial_assembly) {
516  // solve the well equations as a pre-processing step
517  Base::solveWellEq(mob_perfcells, b_perfcells, state, well_state);
518  }
519 
520  Base::computeWellFlux(state, mob_perfcells, b_perfcells, aliveWells, cq_s);
521 
522  if (has_plyshlog_) {
523  std::vector<double> water_vel_wells;
524  std::vector<double> visc_mult_wells;
525 
526  const int water_pos = fluid_.phaseUsage().phase_pos[Water];
527  computeWaterShearVelocityWells(state, well_state, cq_s[water_pos], water_vel_wells, visc_mult_wells);
528 
529  if ( !polymer_props_ad_.computeShearMultLog(water_vel_wells, visc_mult_wells, shear_mult_wells_) ) {
530  OPM_THROW(std::runtime_error, " failed in calculating the shear factors for wells ");
531  }
532 
533  // applying the shear-thinning to the water phase
534  V shear_mult_wells_v = Eigen::Map<V>(shear_mult_wells_.data(), shear_mult_wells_.size());
535  ADB shear_mult_wells_adb = ADB::constant(shear_mult_wells_v);
536  mob_perfcells[water_pos] = mob_perfcells[water_pos] / shear_mult_wells_adb;
537  }
538 
539  Base::computeWellFlux(state, mob_perfcells, b_perfcells, aliveWells, cq_s);
540  Base::updatePerfPhaseRatesAndPressures(cq_s, state, well_state);
541  Base::addWellFluxEq(cq_s, state);
542  addWellContributionToMassBalanceEq(cq_s, state, well_state);
543  addWellControlEq(state, well_state, aliveWells);
544  }
545 
546 
547 
548 
549  template <class Grid>
550  ADB
552  {
553  return polymer_props_ad_.polymerWaterVelocityRatio(state.concentration);
554  }
555 
556 
557 
558 
559  template<class Grid>
560  void
561  BlackoilPolymerModel<Grid>::computeWaterShearVelocityFaces(const V& transi, const std::vector<ADB>& kr,
562  const std::vector<ADB>& phasePressure, const SolutionState& state,
563  std::vector<double>& water_vel, std::vector<double>& visc_mult)
564  {
565 
566  const int phase = fluid_.phaseUsage().phase_pos[Water]; // water position
567 
568  const int canonicalPhaseIdx = canph_[phase];
569 
570  const std::vector<PhasePresence> cond = phaseCondition();
571 
572  const ADB tr_mult = transMult(state.pressure);
573  const ADB mu = fluidViscosity(canonicalPhaseIdx, phasePressure[canonicalPhaseIdx], state.temperature, state.rs, state.rv, cond);
574  rq_[phase].mob = tr_mult * kr[canonicalPhaseIdx] / mu;
575 
576  // compute gravity potensial using the face average as in eclipse and MRST
577  const ADB rho = fluidDensity(canonicalPhaseIdx, rq_[phase].b, state.rs, state.rv);
578  const ADB rhoavg = ops_.caver * rho;
579  rq_[ phase ].dh = ops_.ngrad * phasePressure[ canonicalPhaseIdx ] - geo_.gravity()[2] * (rhoavg * (ops_.ngrad * geo_.z().matrix()));
580  if (use_threshold_pressure_) {
581  applyThresholdPressures(rq_[ phase ].dh);
582  }
583 
584  const ADB& b = rq_[ phase ].b;
585  const ADB& mob = rq_[ phase ].mob;
586  const ADB& dh = rq_[ phase ].dh;
587  UpwindSelector<double> upwind(grid_, ops_, dh.value());
588 
589  const ADB cmax = ADB::constant(cmax_, state.concentration.blockPattern());
590  const ADB mc = computeMc(state);
591  ADB krw_eff = polymer_props_ad_.effectiveRelPerm(state.concentration,
592  cmax,
593  kr[canonicalPhaseIdx]);
594  ADB inv_wat_eff_visc = polymer_props_ad_.effectiveInvWaterVisc(state.concentration, mu.value().data());
595  rq_[ phase ].mob = tr_mult * krw_eff * inv_wat_eff_visc;
596 
597  const V& polymer_conc = state.concentration.value();
598  V visc_mult_cells = polymer_props_ad_.viscMult(polymer_conc);
599  V visc_mult_faces = upwind.select(visc_mult_cells);
600 
601  size_t nface = visc_mult_faces.size();
602  visc_mult.resize(nface);
603  std::copy(visc_mult_faces.data(), visc_mult_faces.data() + nface, visc_mult.begin());
604 
605  rq_[ phase ].mflux = (transi * upwind.select(b * mob)) * dh;
606 
607 
608  const auto& b_faces_adb = upwind.select(b);
609  std::vector<double> b_faces(b_faces_adb.value().data(), b_faces_adb.value().data() + b_faces_adb.size());
610 
611  const auto& internal_faces = ops_.internal_faces;
612 
613  std::vector<double> internal_face_areas;
614  internal_face_areas.resize(internal_faces.size());
615 
616  for (int i = 0; i < internal_faces.size(); ++i) {
617  internal_face_areas[i] = grid_.face_areas[internal_faces[i]];
618  }
619 
620  const ADB phi = Opm::AutoDiffBlock<double>::constant(Eigen::Map<const V>(& fluid_.porosity()[0], AutoDiffGrid::numCells(grid_), 1));
621  const ADB phiavg_adb = ops_.caver * phi;
622 
623  std::vector<double> phiavg(phiavg_adb.value().data(), phiavg_adb.value().data() + phiavg_adb.size());
624 
625  water_vel.resize(nface);
626  std::copy(rq_[0].mflux.value().data(), rq_[0].mflux.value().data() + nface, water_vel.begin());
627 
628  for (size_t i = 0; i < nface; ++i) {
629  water_vel[i] = water_vel[i] / (b_faces[i] * phiavg[i] * internal_face_areas[i]);
630  }
631 
632  // for SHRATE keyword treatment
633  if (has_shrate_) {
634 
635  // get the upwind water saturation
636  const Opm::PhaseUsage pu = fluid_.phaseUsage();
637  const ADB& sw = state.saturation[pu.phase_pos[ Water ]];
638  const ADB& sw_upwind_adb = upwind.select(sw);
639  std::vector<double> sw_upwind(sw_upwind_adb.value().data(), sw_upwind_adb.value().data() + sw_upwind_adb.size());
640 
641  // get the absolute permeability for the faces
642  std::vector<double> perm;
643  perm.resize(transi.size());
644 
645  for (int i = 0; i < transi.size(); ++i) {
646  perm[i] = transi[i] / internal_faces[i];
647  }
648 
649  // get the upwind krw_eff
650  const ADB& krw_adb = upwind.select(krw_eff);
651  std::vector<double> krw_upwind(krw_adb.value().data(), krw_adb.value().data() + krw_adb.size());
652 
653  const double& shrate_const = polymer_props_ad_.shrate();
654 
655  const double epsilon = std::numeric_limits<double>::epsilon();
656  // std::cout << "espilon is " << epsilon << std::endl;
657  // std::cin.ignore();
658 
659  for (size_t i = 0; i < water_vel.size(); ++i) {
660  // assuming only when upwinding water saturation is not zero
661  // there will be non-zero water velocity
662  if (std::abs(water_vel[i]) < epsilon) {
663  continue;
664  }
665 
666  water_vel[i] *= shrate_const * std::sqrt(phiavg[i] / (perm[i] * sw_upwind[i] * krw_upwind[i]));
667 
668  }
669  }
670 
671  }
672 
673 
674 
675 
676  template<class Grid>
677  void
679  std::vector<double>& water_vel_wells, std::vector<double>& visc_mult_wells)
680  {
681  if( ! wellsActive() ) return ;
682 
683  const int nw = wells().number_of_wells;
684  const int nperf = wells().well_connpos[nw];
685  const std::vector<int> well_cells(wells().well_cells, wells().well_cells + nperf);
686 
687  water_vel_wells.resize(cq_sw.size());
688  std::copy(cq_sw.value().data(), cq_sw.value().data() + cq_sw.size(), water_vel_wells.begin());
689 
690  const V& polymer_conc = state.concentration.value();
691 
692  V visc_mult_cells = polymer_props_ad_.viscMult(polymer_conc);
693  V visc_mult_wells_v = subset(visc_mult_cells, well_cells);
694 
695  visc_mult_wells.resize(visc_mult_wells_v.size());
696  std::copy(visc_mult_wells_v.data(), visc_mult_wells_v.data() + visc_mult_wells_v.size(), visc_mult_wells.begin());
697 
698  const int water_pos = fluid_.phaseUsage().phase_pos[Water];
699  ADB b_perfcells = subset(rq_[water_pos].b, well_cells);
700 
701  const ADB& p_perfcells = subset(state.pressure, well_cells);
702  const V& cdp = well_perforation_pressure_diffs_;
703  const ADB perfpressure = (wops_.w2p * state.bhp) + cdp;
704  // Pressure drawdown (also used to determine direction of flow)
705  const ADB drawdown = p_perfcells - perfpressure;
706 
707  // selects injection perforations
708  V selectInjectingPerforations = V::Zero(nperf);
709  for (int c = 0; c < nperf; ++c) {
710  if (drawdown.value()[c] < 0) {
711  selectInjectingPerforations[c] = 1;
712  }
713  }
714 
715  // for the injection wells
716  for (size_t i = 0; i < well_cells.size(); ++i) {
717  if (xw.polymerInflow()[well_cells[i]] == 0. && selectInjectingPerforations[i] == 1) { // maybe comparison with epsilon threshold
718  visc_mult_wells[i] = 1.;
719  }
720  }
721 
722  const ADB phi = Opm::AutoDiffBlock<double>::constant(Eigen::Map<const V>(& fluid_.porosity()[0], AutoDiffGrid::numCells(grid_), 1));
723  const ADB phi_wells_adb = subset(phi, well_cells);
724 
725  std::vector<double> phi_wells(phi_wells_adb.value().data(), phi_wells_adb.value().data() + phi_wells_adb.size());
726 
727  std::vector<double> b_wells(b_perfcells.value().data(), b_perfcells.value().data() + b_perfcells.size());
728 
729  for (size_t i = 0; i < water_vel_wells.size(); ++i) {
730  water_vel_wells[i] = b_wells[i] * water_vel_wells[i] / (phi_wells[i] * 2. * M_PI * wells_rep_radius_[i] * wells_perf_length_[i]);
731  // TODO: CHECK to make sure this formulation is corectly used. Why muliplied by bW.
732  // Although this formulation works perfectly with the tests compared with other formulations
733  }
734 
735  // for SHRATE treatment
736  if (has_shrate_) {
737  const double& shrate_const = polymer_props_ad_.shrate();
738  for (size_t i = 0; i < water_vel_wells.size(); ++i) {
739  water_vel_wells[i] = shrate_const * water_vel_wells[i] / wells_bore_diameter_[i];
740  }
741  }
742 
743  return;
744  }
745 
746 } // namespace Opm
747 
748 #endif // OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
void assembleMassBalanceEq(const SolutionState &state)
Definition: BlackoilPolymerModel_impl.hpp:273
void computeWaterShearVelocityWells(const SolutionState &state, WellState &xw, const ADB &cq_sw, std::vector< double > &water_vel_wells, std::vector< double > &visc_mult_wells)
Definition: BlackoilPolymerModel_impl.hpp:678
Definition: CompressibleTpfaPolymer.hpp:32
std::vector< double > & sat
Definition: GravityColumnSolverPolymer_impl.hpp:76
std::vector< double > & cmax
Definition: GravityColumnSolverPolymer_impl.hpp:78
void makeConstantState(SolutionState &state) const
Definition: BlackoilPolymerModel_impl.hpp:151
void computeCmax(ReservoirState &state)
Definition: BlackoilPolymerModel_impl.hpp:256
void prepareStep(const double dt, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilPolymerModel_impl.hpp:123
BlackoilModelBase< Grid, BlackoilPolymerModel< Grid > > Base
Definition: BlackoilPolymerModel.hpp:49
void computeAccum(const SolutionState &state, const int aix)
Definition: BlackoilPolymerModel_impl.hpp:225
const int poly_pos_
Definition: BlackoilPolymerModel.hpp:141
ADB computeMc(const SolutionState &state) const
Definition: BlackoilPolymerModel_impl.hpp:551
BlackoilPolymerModel(const typename Base::ModelParameters &param, const Grid &grid, const BlackoilPropsAdInterface &fluid, const DerivedGeology &geo, const RockCompressibility *rock_comp_props, const PolymerPropsAd &polymer_props_ad, const Wells *wells, const NewtonIterationBlackoilInterface &linsolver, EclipseStateConstPtr eclipse_state, const bool has_disgas, const bool has_vapoil, const bool has_polymer, const bool has_plyshlog, const bool has_shrate, const std::vector< double > &wells_rep_radius, const std::vector< double > &wells_perf_length, const std::vector< double > &wells_bore_diameter, const bool terminal_output)
Definition: BlackoilPolymerModel_impl.hpp:77
Base::SolutionState SolutionState
Definition: BlackoilPolymerModel.hpp:131
void assemble(const ReservoirState &reservoir_state, WellState &well_state, const bool initial_assembly)
Definition: BlackoilPolymerModel_impl.hpp:460
std::vector< V > variableStateInitials(const ReservoirState &x, const WellState &xw) const
Definition: BlackoilPolymerModel_impl.hpp:163
void computeWaterShearVelocityFaces(const V &transi, const std::vector< ADB > &kr, const std::vector< ADB > &phasePressure, const SolutionState &state, std::vector< double > &water_vel, std::vector< double > &visc_mult)
Definition: BlackoilPolymerModel_impl.hpp:561
Definition: PolymerPropsAd.hpp:32
void addWellContributionToMassBalanceEq(const std::vector< ADB > &cq_s, const SolutionState &state, WellState &xw)
Definition: BlackoilPolymerModel_impl.hpp:346
void updateState(const V &dx, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilPolymerModel_impl.hpp:376
std::vector< int > variableStateIndices() const
Definition: BlackoilPolymerModel_impl.hpp:188
Base::ReservoirState ReservoirState
Definition: BlackoilPolymerModel.hpp:50
int polymerPos(const PU &pu)
Definition: BlackoilPolymerModel_impl.hpp:59
SolutionState variableStateExtractVars(const ReservoirState &x, const std::vector< int > &indices, std::vector< ADB > &vars) const
Definition: BlackoilPolymerModel_impl.hpp:208
const bool has_polymer_
Definition: BlackoilPolymerModel.hpp:138
Base::WellState WellState
Definition: BlackoilPolymerModel.hpp:51
void computeMassFlux(const int actph, const V &transi, const ADB &kr, const ADB &p, const SolutionState &state)
Definition: BlackoilPolymerModel_impl.hpp:413
void afterStep(const double dt, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilPolymerModel_impl.hpp:138