blackoilnewtonmethod.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_BLACK_OIL_NEWTON_METHOD_HH
29#define EWOMS_BLACK_OIL_NEWTON_METHOD_HH
30
31#include "blackoilproperties.hh"
32
33#include <opm/common/Exceptions.hpp>
34
38
39namespace Opm::Properties {
40
41template <class TypeTag, class MyTypeTag>
42struct DiscNewtonMethod;
43
44template<class TypeTag, class MyTypeTag>
45struct DpMaxRel { using type = UndefinedProperty; };
46template<class TypeTag, class MyTypeTag>
47struct DsMax { using type = UndefinedProperty; };
48template<class TypeTag, class MyTypeTag>
50template<class TypeTag, class MyTypeTag>
52template<class TypeTag, class MyTypeTag>
54template<class TypeTag, class MyTypeTag>
56template<class TypeTag, class MyTypeTag>
58template<class TypeTag, class MyTypeTag>
60template<class TypeTag, class MyTypeTag>
62template<class TypeTag, class MyTypeTag>
64template<class TypeTag, class MyTypeTag>
66template<class TypeTag>
67struct DpMaxRel<TypeTag, TTag::NewtonMethod>
68{
70 static constexpr type value = 0.3;
71};
72template<class TypeTag>
73struct DsMax<TypeTag, TTag::NewtonMethod>
74{
76 static constexpr type value = 0.2;
77};
78template<class TypeTag>
80{
82 static constexpr type value = 1e-5;
83};
84template<class TypeTag>
85struct ProjectSaturations<TypeTag, TTag::NewtonMethod> { static constexpr bool value = false; };
86template<class TypeTag>
87struct MaxTemperatureChange<TypeTag, TTag::NewtonMethod>
88{
90 static constexpr type value = 5; //Kelvin
91};
92template<class TypeTag>
93struct TemperatureMax<TypeTag, TTag::NewtonMethod>
94{
96 static constexpr type value = 1e9; //Kelvin
97};
98template<class TypeTag>
99struct TemperatureMin<TypeTag, TTag::NewtonMethod>
100{
102 static constexpr type value = 0.0; //Kelvin
103};
104template<class TypeTag>
105struct PressureMax<TypeTag, TTag::NewtonMethod>
106{
108 static constexpr type value = 1e99; //Kelvin
109};
110template<class TypeTag>
111struct PressureMin<TypeTag, TTag::NewtonMethod>
112{
114 static constexpr type value = -1e99; //Kelvin
115};
116template<class TypeTag>
118{
120 static constexpr type value = 1.0;
121};
122template<class TypeTag>
123struct WaterOnlyThreshold<TypeTag, TTag::NewtonMethod>
124{
126 static constexpr type value = 1.0;
127};
128} // namespace Opm::Properties
129
130namespace Opm {
131
137template <class TypeTag>
138class BlackOilNewtonMethod : public GetPropType<TypeTag, Properties::DiscNewtonMethod>
139{
151
152 static const unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
153 static constexpr bool enableSaltPrecipitation = getPropValue<TypeTag, Properties::EnableSaltPrecipitation>();
154
155public:
156 BlackOilNewtonMethod(Simulator& simulator) : ParentType(simulator)
157 {
158 priVarOscilationThreshold_ = Parameters::get<TypeTag, Properties::PriVarOscilationThreshold>();
159 dpMaxRel_ = Parameters::get<TypeTag, Properties::DpMaxRel>();
160 dsMax_ = Parameters::get<TypeTag, Properties::DsMax>();
161 projectSaturations_ = Parameters::get<TypeTag, Properties::ProjectSaturations>();
162 maxTempChange_ = Parameters::get<TypeTag, Properties::MaxTemperatureChange>();
163 tempMax_ = Parameters::get<TypeTag, Properties::TemperatureMax>();
164 tempMin_ = Parameters::get<TypeTag, Properties::TemperatureMin>();
165 pressMax_ = Parameters::get<TypeTag, Properties::PressureMax>();
166 pressMin_ = Parameters::get<TypeTag, Properties::PressureMin>();
167 waterSaturationMax_ = Parameters::get<TypeTag, Properties::MaximumWaterSaturation>();
168 waterOnlyThreshold_ = Parameters::get<TypeTag, Properties::WaterOnlyThreshold>();
169 }
170
175 {
176 ParentType::finishInit();
177
178 wasSwitched_.resize(this->model().numTotalDof());
179 std::fill(wasSwitched_.begin(), wasSwitched_.end(), false);
180 }
181
185 static void registerParameters()
186 {
187 ParentType::registerParameters();
188
189 Parameters::registerParam<TypeTag, Properties::DpMaxRel>
190 ("Maximum relative change of pressure in a single iteration");
191 Parameters::registerParam<TypeTag, Properties::DsMax>
192 ("Maximum absolute change of any saturation in a single iteration");
193 Parameters::registerParam<TypeTag, Properties::PriVarOscilationThreshold>
194 ("The threshold value for the primary variable switching conditions "
195 "after its meaning has switched to hinder oscilations");
196 Parameters::registerParam<TypeTag, Properties::ProjectSaturations>
197 ("Option for doing saturation projection");
198 Parameters::registerParam<TypeTag, Properties::MaxTemperatureChange>
199 ("Maximum absolute change of temperature in a single iteration");
200 Parameters::registerParam<TypeTag, Properties::TemperatureMax>
201 ("Maximum absolute temperature");
202 Parameters::registerParam<TypeTag, Properties::TemperatureMin>
203 ("Minimum absolute temperature");
204 Parameters::registerParam<TypeTag, Properties::PressureMax>
205 ("Maximum absolute pressure");
206 Parameters::registerParam<TypeTag, Properties::PressureMin>
207 ("Minimum absolute pressure");
208 Parameters::registerParam<TypeTag, Properties::MaximumWaterSaturation>
209 ("Maximum water saturation");
210 Parameters::registerParam<TypeTag, Properties::WaterOnlyThreshold>
211 ("Cells with water saturation above or equal is considered one-phase water only");
212 }
213
218 unsigned numPriVarsSwitched() const
219 { return numPriVarsSwitched_; }
220
221protected:
224
229 {
230 numPriVarsSwitched_ = 0;
231 ParentType::beginIteration_();
232 }
233
237 void endIteration_(SolutionVector& uCurrentIter,
238 const SolutionVector& uLastIter)
239 {
240#if HAVE_MPI
241 // in the MPI enabled case we need to add up the number of DOF
242 // for which the interpretation changed over all processes.
243 int localSwitched = numPriVarsSwitched_;
244 MPI_Allreduce(&localSwitched,
245 &numPriVarsSwitched_,
246 /*num=*/1,
247 MPI_INT,
248 MPI_SUM,
249 MPI_COMM_WORLD);
250#endif // HAVE_MPI
251
252 this->simulator_.model().newtonMethod().endIterMsg()
253 << ", num switched=" << numPriVarsSwitched_;
254
255 ParentType::endIteration_(uCurrentIter, uLastIter);
256 }
257
258public:
259 void update_(SolutionVector& nextSolution,
260 const SolutionVector& currentSolution,
261 const GlobalEqVector& solutionUpdate,
262 const GlobalEqVector& currentResidual)
263 {
264 const auto& comm = this->simulator_.gridView().comm();
265
266 int succeeded;
267 try {
268 ParentType::update_(nextSolution,
269 currentSolution,
270 solutionUpdate,
271 currentResidual);
272 succeeded = 1;
273 }
274 catch (...) {
275 succeeded = 0;
276 }
277 succeeded = comm.min(succeeded);
278
279 if (!succeeded)
280 throw NumericalProblem("A process did not succeed in adapting the primary variables");
281
282 numPriVarsSwitched_ = comm.sum(numPriVarsSwitched_);
283 }
284
285 template <class DofIndices>
286 void update_(SolutionVector& nextSolution,
287 const SolutionVector& currentSolution,
288 const GlobalEqVector& solutionUpdate,
289 const GlobalEqVector& currentResidual,
290 const DofIndices& dofIndices)
291 {
292 const auto zero = 0.0 * solutionUpdate[0];
293 for (auto dofIdx : dofIndices) {
294 if (solutionUpdate[dofIdx] == zero) {
295 continue;
296 }
298 nextSolution[dofIdx],
299 currentSolution[dofIdx],
300 solutionUpdate[dofIdx],
301 currentResidual[dofIdx]);
302 }
303 }
304
305protected:
309 void updatePrimaryVariables_(unsigned globalDofIdx,
310 PrimaryVariables& nextValue,
311 const PrimaryVariables& currentValue,
312 const EqVector& update,
313 const EqVector& currentResidual)
314 {
315 static constexpr bool enableSolvent = Indices::solventSaturationIdx >= 0;
316 static constexpr bool enableExtbo = Indices::zFractionIdx >= 0;
317 static constexpr bool enablePolymer = Indices::polymerConcentrationIdx >= 0;
318 static constexpr bool enablePolymerWeight = Indices::polymerMoleWeightIdx >= 0;
319 static constexpr bool enableEnergy = Indices::temperatureIdx >= 0;
320 static constexpr bool enableFoam = Indices::foamConcentrationIdx >= 0;
321 static constexpr bool enableBrine = Indices::saltConcentrationIdx >= 0;
322 static constexpr bool enableMICP = Indices::microbialConcentrationIdx >= 0;
323
324 currentValue.checkDefined();
325 Valgrind::CheckDefined(update);
326 Valgrind::CheckDefined(currentResidual);
327
328 // saturation delta for each phase
329 Scalar deltaSw = 0.0;
330 Scalar deltaSo = 0.0;
331 Scalar deltaSg = 0.0;
332 Scalar deltaSs = 0.0;
333
334 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
335 {
336 deltaSw = update[Indices::waterSwitchIdx];
337 deltaSo -= deltaSw;
338 }
339 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
340 {
341 deltaSg = update[Indices::compositionSwitchIdx];
342 deltaSo -= deltaSg;
343 }
344 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
345 deltaSs = update[Indices::solventSaturationIdx];
346 deltaSo -= deltaSs;
347 }
348
349 // maximum saturation delta
350 Scalar maxSatDelta = std::max(std::abs(deltaSg), std::abs(deltaSo));
351 maxSatDelta = std::max(maxSatDelta, std::abs(deltaSw));
352 maxSatDelta = std::max(maxSatDelta, std::abs(deltaSs));
353
354 // scaling factor for saturation deltas to make sure that none of them exceeds
355 // the specified threshold value.
356 Scalar satAlpha = 1.0;
357 if (maxSatDelta > dsMax_)
358 satAlpha = dsMax_/maxSatDelta;
359
360 for (int pvIdx = 0; pvIdx < int(numEq); ++pvIdx) {
361 // calculate the update of the current primary variable. For the black-oil
362 // model we limit the pressure delta relative to the pressure's current
363 // absolute value (Default: 30%) and saturation deltas to an absolute change
364 // (Default: 20%). Further, we ensure that the R factors, solvent
365 // "saturation" and polymer concentration do not become negative after the
366 // update.
367 Scalar delta = update[pvIdx];
368
369 // limit pressure delta
370 if (pvIdx == Indices::pressureSwitchIdx) {
371 if (std::abs(delta) > dpMaxRel_*currentValue[pvIdx])
372 delta = signum(delta)*dpMaxRel_*currentValue[pvIdx];
373 }
374 // water saturation delta
375 else if (pvIdx == Indices::waterSwitchIdx)
376 if (currentValue.primaryVarsMeaningWater() == PrimaryVariables::WaterMeaning::Sw)
377 delta *= satAlpha;
378 else {
379 //Ensure Rvw and Rsw factor does not become negative
380 if (delta > currentValue[ Indices::waterSwitchIdx])
381 delta = currentValue[ Indices::waterSwitchIdx];
382 }
383 else if (pvIdx == Indices::compositionSwitchIdx) {
384 // the switching primary variable for composition is tricky because the
385 // "reasonable" value ranges it exhibits vary widely depending on its
386 // interpretation since it can represent Sg, Rs or Rv. For now, we only
387 // limit saturation deltas and ensure that the R factors do not become
388 // negative.
389 if (currentValue.primaryVarsMeaningGas() == PrimaryVariables::GasMeaning::Sg)
390 delta *= satAlpha;
391 else {
392 //Ensure Rv and Rs factor does not become negative
393 if (delta > currentValue[Indices::compositionSwitchIdx])
394 delta = currentValue[Indices::compositionSwitchIdx];
395 }
396 }
397 else if (enableSolvent && pvIdx == Indices::solventSaturationIdx) {
398 // solvent saturation updates are also subject to the Appleyard chop
399 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss) {
400 delta *= satAlpha;
401 } else {
402 //Ensure Rssolw factor does not become negative
403 if (delta > currentValue[Indices::solventSaturationIdx])
404 delta = currentValue[Indices::solventSaturationIdx];
405 }
406 }
407 else if (enableExtbo && pvIdx == Indices::zFractionIdx) {
408 // z fraction updates are also subject to the Appleyard chop
409 const auto& curr = currentValue[Indices::zFractionIdx]; // or currentValue[pvIdx] given the block condition
410 delta = std::clamp(delta, curr - 1.0, curr);
411 }
412 else if (enablePolymerWeight && pvIdx == Indices::polymerMoleWeightIdx) {
413 const double sign = delta >= 0. ? 1. : -1.;
414 // maximum change of polymer molecular weight, the unit is MDa.
415 // applying this limit to stabilize the simulation. The value itself is still experimental.
416 const double maxMolarWeightChange = 100.0;
417 delta = sign * std::min(std::abs(delta), maxMolarWeightChange);
418 delta *= satAlpha;
419 }
420 else if (enableEnergy && pvIdx == Indices::temperatureIdx) {
421 const double sign = delta >= 0. ? 1. : -1.;
422 delta = sign * std::min(std::abs(delta), maxTempChange_);
423 }
424 else if (enableBrine && pvIdx == Indices::saltConcentrationIdx &&
425 enableSaltPrecipitation &&
426 currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp) {
427 const double maxSaltSaturationChange = 0.1;
428 const double sign = delta >= 0. ? 1. : -1.;
429 delta = sign * std::min(std::abs(delta), maxSaltSaturationChange);
430 }
431
432 // do the actual update
433 nextValue[pvIdx] = currentValue[pvIdx] - delta;
434
435 // keep the solvent saturation between 0 and 1
436 if (enableSolvent && pvIdx == Indices::solventSaturationIdx) {
437 if (currentValue.primaryVarsMeaningSolvent() == PrimaryVariables::SolventMeaning::Ss )
438 nextValue[pvIdx] = std::min(std::max(nextValue[pvIdx], 0.0), 1.0);
439 }
440
441 // keep the z fraction between 0 and 1
442 if (enableExtbo && pvIdx == Indices::zFractionIdx)
443 nextValue[pvIdx] = std::min(std::max(nextValue[pvIdx], 0.0), 1.0);
444
445 // keep the polymer concentration above 0
446 if (enablePolymer && pvIdx == Indices::polymerConcentrationIdx)
447 nextValue[pvIdx] = std::max(nextValue[pvIdx], 0.0);
448
449 if (enablePolymerWeight && pvIdx == Indices::polymerMoleWeightIdx) {
450 nextValue[pvIdx] = std::max(nextValue[pvIdx], 0.0);
451 const double polymerConcentration = nextValue[Indices::polymerConcentrationIdx];
452 if (polymerConcentration < 1.e-10)
453 nextValue[pvIdx] = 0.0;
454 }
455
456 // keep the foam concentration above 0
457 if (enableFoam && pvIdx == Indices::foamConcentrationIdx)
458 nextValue[pvIdx] = std::max(nextValue[pvIdx], 0.0);
459
460 if (enableBrine && pvIdx == Indices::saltConcentrationIdx) {
461 // keep the salt concentration above 0
462 if (!enableSaltPrecipitation || (enableSaltPrecipitation && currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Cs))
463 nextValue[pvIdx] = std::max(nextValue[pvIdx], 0.0);
464 // keep the salt saturation below upperlimit
465 if ((enableSaltPrecipitation && currentValue.primaryVarsMeaningBrine() == PrimaryVariables::BrineMeaning::Sp))
466 nextValue[pvIdx] = std::min(nextValue[pvIdx], 1.0-1.e-8);
467 }
468
469 // keep the temperature within given values
470 if (enableEnergy && pvIdx == Indices::temperatureIdx)
471 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], tempMin_, tempMax_);
472
473 if (pvIdx == Indices::pressureSwitchIdx) {
474 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], pressMin_, pressMax_);
475 }
476
477
478 // Limit the variables to [0, cmax] values to improve the convergence.
479 // For the microorganisms we set this value equal to the biomass density value.
480 // For the oxygen and urea we set this value to the maximum injected
481 // concentration (the urea concentration has been scaled by 10). For
482 // the biofilm and calcite, we set this value equal to the porosity minus the clogging tolerance.
483 if (enableMICP && pvIdx == Indices::microbialConcentrationIdx)
484 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], 0.0, MICPModule::densityBiofilm());
485 if (enableMICP && pvIdx == Indices::oxygenConcentrationIdx)
486 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], 0.0, MICPModule::maximumOxygenConcentration());
487 if (enableMICP && pvIdx == Indices::ureaConcentrationIdx)
488 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], 0.0, MICPModule::maximumUreaConcentration());
489 if (enableMICP && pvIdx == Indices::biofilmConcentrationIdx)
490 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], 0.0, MICPModule::phi()[globalDofIdx] - MICPModule::toleranceBeforeClogging());
491 if (enableMICP && pvIdx == Indices::calciteConcentrationIdx)
492 nextValue[pvIdx] = std::clamp(nextValue[pvIdx], 0.0, MICPModule::phi()[globalDofIdx] - MICPModule::toleranceBeforeClogging());
493 }
494
495 // switch the new primary variables to something which is physically meaningful.
496 // use a threshold value after a switch to make it harder to switch back
497 // immediately.
498 if (wasSwitched_[globalDofIdx])
499 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(), globalDofIdx, waterSaturationMax_, waterOnlyThreshold_, priVarOscilationThreshold_);
500 else
501 wasSwitched_[globalDofIdx] = nextValue.adaptPrimaryVariables(this->problem(), globalDofIdx, waterSaturationMax_, waterOnlyThreshold_);
502
503 if (wasSwitched_[globalDofIdx])
504 ++ numPriVarsSwitched_;
505 if(projectSaturations_){
506 nextValue.chopAndNormalizeSaturations();
507 }
508
509 nextValue.checkDefined();
510 }
511
512private:
513 int numPriVarsSwitched_;
514
515 Scalar priVarOscilationThreshold_;
516 Scalar waterSaturationMax_;
517 Scalar waterOnlyThreshold_;
518
519 Scalar dpMaxRel_;
520 Scalar dsMax_;
521 bool projectSaturations_;
522 Scalar maxTempChange_;
523 Scalar tempMax_;
524 Scalar tempMin_;
525 Scalar pressMax_;
526 Scalar pressMin_;
527
528 // keep track of cells where the primary variable meaning has changed
529 // to detect and hinder oscillations
530 std::vector<bool> wasSwitched_;
531};
532} // namespace Opm
533
534#endif
Contains the classes required to extend the black-oil model by MICP.
Declares the properties required by the black oil model.
Contains the high level supplements required to extend the black oil model by MICP.
Definition: blackoilmicpmodules.hh:56
static const Scalar densityBiofilm()
Definition: blackoilmicpmodules.hh:351
static const Scalar toleranceBeforeClogging()
Definition: blackoilmicpmodules.hh:426
static const Scalar maximumUreaConcentration()
Definition: blackoilmicpmodules.hh:396
static const std::vector< Scalar > phi()
Definition: blackoilmicpmodules.hh:436
static const Scalar maximumOxygenConcentration()
Definition: blackoilmicpmodules.hh:391
A newton solver which is specific to the black oil model.
Definition: blackoilnewtonmethod.hh:139
void beginIteration_()
Indicates the beginning of a Newton iteration.
Definition: blackoilnewtonmethod.hh:228
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: blackoilnewtonmethod.hh:185
void update_(SolutionVector &nextSolution, const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual)
Definition: blackoilnewtonmethod.hh:259
void finishInit()
Finialize the construction of the object.
Definition: blackoilnewtonmethod.hh:174
void endIteration_(SolutionVector &uCurrentIter, const SolutionVector &uLastIter)
Indicates that one Newton iteration was finished.
Definition: blackoilnewtonmethod.hh:237
unsigned numPriVarsSwitched() const
Returns the number of degrees of freedom for which the interpretation has changed for the most recent...
Definition: blackoilnewtonmethod.hh:218
BlackOilNewtonMethod(Simulator &simulator)
Definition: blackoilnewtonmethod.hh:156
void updatePrimaryVariables_(unsigned globalDofIdx, PrimaryVariables &nextValue, const PrimaryVariables &currentValue, const EqVector &update, const EqVector &currentResidual)
Update a single primary variables object.
Definition: blackoilnewtonmethod.hh:309
friend ParentType
Definition: blackoilnewtonmethod.hh:223
void update_(SolutionVector &nextSolution, const SolutionVector &currentSolution, const GlobalEqVector &solutionUpdate, const GlobalEqVector &currentResidual, const DofIndices &dofIndices)
Definition: blackoilnewtonmethod.hh:286
The multi-dimensional Newton method.
Definition: newtonmethod.hh:113
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:242
int signum(Scalar val)
Template function which returns the sign of a floating point value.
Definition: signum.hh:40
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:69
Definition: blackoilnewtonmethod.hh:45
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:75
Definition: blackoilnewtonmethod.hh:47
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:89
Definition: blackoilnewtonmethod.hh:53
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:119
Definition: blackoilnewtonmethod.hh:63
Specifies the type of the actual Newton method.
Definition: newtonmethodproperties.hh:32
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:107
Definition: blackoilnewtonmethod.hh:59
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:113
Definition: blackoilnewtonmethod.hh:61
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:81
Definition: blackoilnewtonmethod.hh:49
Definition: blackoilnewtonmethod.hh:51
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:95
Definition: blackoilnewtonmethod.hh:55
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:101
Definition: blackoilnewtonmethod.hh:57
a tag to mark properties as undefined
Definition: propertysystem.hh:40
GetPropType< TypeTag, Scalar > type
Definition: blackoilnewtonmethod.hh:125
Definition: blackoilnewtonmethod.hh:65