blackoilmodel.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 OPM_BLACK_OIL_MODEL_HPP
29#define OPM_BLACK_OIL_MODEL_HPP
30
31#include <opm/material/densead/Math.hpp>
32
33#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
34
48
50
54
56
57#include <cassert>
58#include <istream>
59#include <limits>
60#include <memory>
61#include <ostream>
62#include <sstream>
63#include <stdexcept>
64#include <string>
65#include <tuple>
66#include <vector>
67
68namespace Opm {
69
70template <class TypeTag>
71class BlackOilModel;
72
73}
74
75namespace Opm::Properties {
76
77namespace TTag {
78
81{ using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
82
83} // namespace TTag
84
86template<class TypeTag>
87struct LocalResidual<TypeTag, TTag::BlackOilModel>
89
91template<class TypeTag>
92struct NewtonMethod<TypeTag, TTag::BlackOilModel>
94
96template<class TypeTag>
97struct Model<TypeTag, TTag::BlackOilModel>
99
101template<class TypeTag>
102struct BaseProblem<TypeTag, TTag::BlackOilModel>
104
106template<class TypeTag>
107struct RateVector<TypeTag, TTag::BlackOilModel>
109
111template<class TypeTag>
112struct BoundaryRateVector<TypeTag, TTag::BlackOilModel>
114
116template<class TypeTag>
117struct PrimaryVariables<TypeTag, TTag::BlackOilModel>
119
121template<class TypeTag>
122struct IntensiveQuantities<TypeTag, TTag::BlackOilModel>
124
126template<class TypeTag>
127struct ExtensiveQuantities<TypeTag, TTag::BlackOilModel>
129
132template<class TypeTag>
133struct FluxModule<TypeTag, TTag::BlackOilModel>
135
137template<class TypeTag>
138struct Indices<TypeTag, TTag::BlackOilModel>
139{
141 getPropValue<TypeTag, Properties::EnableExtbo>(),
142 getPropValue<TypeTag, Properties::EnablePolymer>(),
143 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal,
144 getPropValue<TypeTag, Properties::EnableFoam>(),
145 getPropValue<TypeTag, Properties::EnableBrine>(),
146 /*PVOffset=*/0,
147 getPropValue<TypeTag, Properties::EnableBioeffects>()>;
148};
149
151template<class TypeTag>
152struct FluidSystem<TypeTag, TTag::BlackOilModel>
153{
154public:
157 using type = BlackOilFluidSystem<Scalar>;
158};
159
160// by default, all ECL extension modules are disabled
161template<class TypeTag>
162struct EnableSolvent<TypeTag, TTag::BlackOilModel>
163{ static constexpr bool value = false; };
164
165template<class TypeTag>
166struct EnableExtbo<TypeTag, TTag::BlackOilModel>
167{ static constexpr bool value = false; };
168
169template<class TypeTag>
170struct EnablePolymer<TypeTag, TTag::BlackOilModel>
171{ static constexpr bool value = false; };
172
173template<class TypeTag>
174struct EnablePolymerMW<TypeTag, TTag::BlackOilModel>
175{ static constexpr bool value = false; };
176
177template<class TypeTag>
178struct EnableFoam<TypeTag, TTag::BlackOilModel>
179{ static constexpr bool value = false; };
180
181template<class TypeTag>
182struct EnableBrine<TypeTag, TTag::BlackOilModel>
183{ static constexpr bool value = false; };
184
185template<class TypeTag>
186struct EnableVapwat<TypeTag, TTag::BlackOilModel>
187{ static constexpr bool value = false; };
188
189template<class TypeTag>
190struct EnableDisgasInWater<TypeTag, TTag::BlackOilModel>
191{ static constexpr bool value = false; };
192
193template<class TypeTag>
195{ static constexpr bool value = false; };
196
197template<class TypeTag>
198struct EnableBioeffects<TypeTag, TTag::BlackOilModel>
199{ static constexpr bool value = false; };
200
201template<class TypeTag>
202struct EnergyModuleType<TypeTag, TTag::BlackOilModel>
203{ static constexpr EnergyModules value = EnergyModules::NoTemperature; };
204
206template<class TypeTag>
207struct EnableDiffusion<TypeTag, TTag::BlackOilModel>
208{ static constexpr bool value = false; };
209
211template<class TypeTag>
212struct EnableDispersion<TypeTag, TTag::BlackOilModel>
213{ static constexpr bool value = false; };
214
215template<class TypeTag>
217{ static constexpr bool value = false; };
218
219template<class TypeTag>
220struct EnableMech<TypeTag, TTag::BlackOilModel>
221{ static constexpr bool value = false; };
222
223template<class TypeTag>
224struct RunAssemblyOnGpu<TypeTag, TTag::BlackOilModel>
225{ static constexpr bool value = false; };
226
233template<class TypeTag>
235{
236private:
238 static constexpr Scalar alpha =
239 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>() ? 1000.0 : 1.0;
240
241public:
242 using type = Scalar;
243 static constexpr Scalar value = 1.0/(30.0*4184.0*alpha);
244};
245
247template<class TypeTag>
249{
250private:
252 static constexpr Scalar alpha =
253 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>() ? 1000.0 : 1.0;
254
255public:
256 using type = Scalar;
257 static constexpr Scalar value = 1.0/(10.0*alpha);
258};
259
260// by default, ebos formulates the conservation equations in terms of mass not surface
261// volumes
262template<class TypeTag>
264{ static constexpr bool value = false; };
265
266} // namespace Opm::Properties
267
268namespace Opm {
269
333template<class TypeTag >
335 : public MultiPhaseBaseModel<TypeTag>
336{
337public:
341
342private:
343 using Implementation = GetPropType<TypeTag, Properties::Model>;
345
350
351 enum { numComponents = FluidSystem::numComponents };
352 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
353
354 static constexpr bool compositionSwitchEnabled =
355 Indices::compositionSwitchIdx != std::numeric_limits<unsigned>::max();
356 static constexpr bool enableBioeffects = getPropValue<TypeTag, Properties::EnableBioeffects>();
357 static constexpr bool enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>();
358 static constexpr bool enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>();
359 static constexpr bool enableExtbo = getPropValue<TypeTag, Properties::EnableExtbo>();
360 static constexpr bool enablePolymer = getPropValue<TypeTag, Properties::EnablePolymer>();
361 static constexpr bool enableSolvent = getPropValue<TypeTag, Properties::EnableSolvent>();
362 static constexpr EnergyModules energyModuleType = getPropValue<TypeTag, Properties::EnergyModuleType>();
363 static constexpr bool enableFullyImplicitThermal = energyModuleType == EnergyModules::FullyImplicitThermal;
364 static constexpr bool waterEnabled = Indices::waterEnabled;
365
366 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag, enableBioeffects>;
369 using EnergyModule = BlackOilEnergyModule<TypeTag, energyModuleType>;
370 using ExtboModule = BlackOilExtboModule<TypeTag, enableExtbo>;
371 using PolymerModule = BlackOilPolymerModule<TypeTag, enablePolymer>;
372 using SolventModule = BlackOilSolventModule<TypeTag, enableSolvent>;
373
374public:
376
377 explicit BlackOilModel(Simulator& simulator)
378 : ParentType(simulator)
379 {
380 eqWeights_.resize(numEq, 1.0);
381 }
382
386 static void registerParameters()
387 {
389
390 if constexpr (enableSolvent) {
391 SolventModule::registerParameters();
392 }
393 if constexpr (enableExtbo) {
394 ExtboModule::registerParameters();
395 }
396 if constexpr (enablePolymer) {
397 PolymerModule::registerParameters();
398 }
399 if constexpr (enableFullyImplicitThermal) {
400 EnergyModule::registerParameters();
401 }
402 if constexpr (enableDiffusion) {
403 DiffusionModule::registerParameters();
404 }
405 if constexpr (enableBioeffects) {
406 BioeffectsModule::registerParameters();
407 }
408
409 // register runtime parameters of the VTK output modules
412 if constexpr (enableDiffusion) {
414 }
415 }
416
420 static std::string name()
421 { return "blackoil"; }
422
426 std::string primaryVarName(unsigned pvIdx) const
427 {
428 if (pvIdx == Indices::waterSwitchIdx) {
429 return "water_switching";
430 }
431 else if (pvIdx == Indices::pressureSwitchIdx) {
432 return "pressure_switching";
433 }
434 else if (pvIdx == Indices::compositionSwitchIdx) {
435 return "composition_switching";
436 }
437
438 if constexpr (enableSolvent) {
439 if (SolventModule::primaryVarApplies(pvIdx)) {
440 return SolventModule::primaryVarName(pvIdx);
441 }
442 }
443
444 if constexpr (enableExtbo) {
445 if (ExtboModule::primaryVarApplies(pvIdx)) {
446 return ExtboModule::primaryVarName(pvIdx);
447 }
448 }
449
450 if constexpr (enablePolymer) {
451 if (PolymerModule::primaryVarApplies(pvIdx)) {
452 return PolymerModule::primaryVarName(pvIdx);
453 }
454 }
455
456 if constexpr (enableFullyImplicitThermal) {
457 if (EnergyModule::primaryVarApplies(pvIdx)) {
458 return EnergyModule::primaryVarName(pvIdx);
459 }
460 }
461
462 throw std::logic_error("Invalid primary variable index");
463 }
464
468 std::string eqName(int eqIdx) const
469 {
470 if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx + numComponents) {
471 std::ostringstream oss;
472 oss << "conti_" << FluidSystem::phaseName(eqIdx - Indices::conti0EqIdx);
473 return oss.str();
474 }
475
476 if constexpr (enableSolvent) {
477 if (SolventModule::eqApplies(eqIdx)) {
478 return SolventModule::eqName(eqIdx);
479 }
480 }
481
482 if constexpr (enableExtbo) {
483 if (ExtboModule::eqApplies(eqIdx)) {
484 return ExtboModule::eqName(eqIdx);
485 }
486 }
487
488 if constexpr (enablePolymer) {
489 if (PolymerModule::eqApplies(eqIdx)) {
490 return PolymerModule::eqName(eqIdx);
491 }
492 }
493
494 if constexpr (enableFullyImplicitThermal) {
495 if (EnergyModule::eqApplies(eqIdx)) {
496 return EnergyModule::eqName(eqIdx);
497 }
498 }
499
500 throw std::logic_error("Invalid equation index");
501 }
502
506 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
507 {
508 // do not care about the auxiliary equations as they are supposed to scale
509 // themselves
510 if (globalDofIdx >= this->numGridDof()) {
511 return 1.0;
512 }
513
514 // saturations are always in the range [0, 1]!
515 if (Indices::waterSwitchIdx == pvIdx) {
516 return 1.0;
517 }
518
519 // oil pressures usually are in the range of 100 to 500 bars for typical oil
520 // reservoirs (which is the only relevant application for the black-oil model).
521 else if (int(Indices::pressureSwitchIdx) == int(pvIdx)) {
522 return 1.0 / 300e5;
523 }
524
525 // deal with primary variables stemming from the solvent module
526 if constexpr (enableSolvent) {
527 if (SolventModule::primaryVarApplies(pvIdx)) {
528 return SolventModule::primaryVarWeight(pvIdx);
529 }
530 }
531
532 // deal with primary variables stemming from the extBO module
533 if constexpr (enableExtbo) {
534 if (ExtboModule::primaryVarApplies(pvIdx)) {
535 return ExtboModule::primaryVarWeight(pvIdx);
536 }
537 }
538
539 // deal with primary variables stemming from the polymer module
540 if constexpr (enablePolymer) {
541 if (PolymerModule::primaryVarApplies(pvIdx)) {
542 return PolymerModule::primaryVarWeight(pvIdx);
543 }
544 }
545
546 // deal with primary variables stemming from the energy module
547 if constexpr (enableFullyImplicitThermal) {
548 if (EnergyModule::primaryVarApplies(pvIdx)) {
549 return EnergyModule::primaryVarWeight(pvIdx);
550 }
551 }
552
553 // if the primary variable is either the gas saturation, Rs or Rv
554 assert(Indices::compositionSwitchIdx == pvIdx);
555
556 switch (this->solution(0)[globalDofIdx].primaryVarsMeaningGas()) {
557 case PrimaryVariables::GasMeaning::Sg: return 1.0; // gas saturation
558 case PrimaryVariables::GasMeaning::Rs: return 1.0 / 250.; // gas dissolution factor
559 case PrimaryVariables::GasMeaning::Rv: return 1.0 / 0.025; // oil vaporization factor
560 default: throw std::logic_error("Invalid primary variable meaning flag for gas");
561 }
562 }
563
570 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
571 {
572 // do not care about the auxiliary equations as they are supposed to scale
573 // themselves
574 if (globalDofIdx >= this->numGridDof()) {
575 return 1.0;
576 }
577
578 return eqWeights_[eqIdx];
579 }
580
581 void setEqWeight(unsigned eqIdx, Scalar value)
582 { eqWeights_[eqIdx] = value; }
583
592 template <class DofEntity>
593 void serializeEntity(std::ostream& outstream, const DofEntity& dof)
594 {
595 const unsigned dofIdx = static_cast<unsigned>(asImp_().dofMapper().index(dof));
596
597 // write phase state
598 if (!outstream.good()) {
599 throw std::runtime_error("Could not serialize degree of freedom " + std::to_string(dofIdx));
600 }
601
602 // write the primary variables
603 const auto& priVars = this->solution(/*timeIdx=*/0)[dofIdx];
604 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
605 outstream << priVars[eqIdx] << " ";
606 }
607
608 // write the pseudo primary variables
609 outstream << static_cast<int>(priVars.primaryVarsMeaningGas()) << " ";
610 outstream << static_cast<int>(priVars.primaryVarsMeaningWater()) << " ";
611 outstream << static_cast<int>(priVars.primaryVarsMeaningPressure()) << " ";
612
613 outstream << priVars.pvtRegionIndex() << " ";
614
615 if constexpr (enableSolvent) {
616 SolventModule::serializeEntity(asImp_(), outstream, dof);
617 }
618 if constexpr (enableExtbo) {
619 ExtboModule::serializeEntity(asImp_(), outstream, dof);
620 }
621 if constexpr (enablePolymer) {
622 PolymerModule::serializeEntity(asImp_(), outstream, dof);
623 }
624 if constexpr (enableFullyImplicitThermal) {
625 EnergyModule::serializeEntity(asImp_(), outstream, dof);
626 }
627 }
628
637 template <class DofEntity>
638 void deserializeEntity(std::istream& instream,
639 const DofEntity& dof)
640 {
641 const unsigned dofIdx = static_cast<unsigned>(asImp_().dofMapper().index(dof));
642
643 // read in the "real" primary variables of the DOF
644 auto& priVars = this->solution(/*timeIdx=*/0)[dofIdx];
645 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
646 if (!instream.good()) {
647 throw std::runtime_error("Could not deserialize degree of freedom " + std::to_string(dofIdx));
648 }
649 instream >> priVars[eqIdx];
650 }
651
652 // read the pseudo primary variables
653 unsigned primaryVarsMeaningGas;
654 instream >> primaryVarsMeaningGas;
655
656 unsigned primaryVarsMeaningWater;
657 instream >> primaryVarsMeaningWater;
658
659 unsigned primaryVarsMeaningPressure;
660 instream >> primaryVarsMeaningPressure;
661
662 unsigned pvtRegionIdx;
663 instream >> pvtRegionIdx;
664
665 if (!instream.good()) {
666 throw std::runtime_error("Could not deserialize degree of freedom " + std::to_string(dofIdx));
667 }
668
669 if constexpr (enableSolvent) {
670 SolventModule::deserializeEntity(asImp_(), instream, dof);
671 }
672 if constexpr (enableExtbo) {
673 ExtboModule::deserializeEntity(asImp_(), instream, dof);
674 }
675 if constexpr (enablePolymer) {
676 PolymerModule::deserializeEntity(asImp_(), instream, dof);
677 }
678 if constexpr (enableFullyImplicitThermal) {
679 EnergyModule::deserializeEntity(asImp_(), instream, dof);
680 }
681
682 using PVM_G = typename PrimaryVariables::GasMeaning;
683 using PVM_W = typename PrimaryVariables::WaterMeaning;
684 using PVM_P = typename PrimaryVariables::PressureMeaning;
685 priVars.setPrimaryVarsMeaningGas(static_cast<PVM_G>(primaryVarsMeaningGas));
686 priVars.setPrimaryVarsMeaningWater(static_cast<PVM_W>(primaryVarsMeaningWater));
687 priVars.setPrimaryVarsMeaningPressure(static_cast<PVM_P>(primaryVarsMeaningPressure));
688
689 priVars.setPvtRegionIndex(pvtRegionIdx);
690 }
691
699 template <class Restarter>
700 void deserialize(Restarter& res)
701 {
702 ParentType::deserialize(res);
703
704 // set the PVT indices of the primary variables. This is also done by writing
705 // them into the restart file and re-reading them, but it is better to calculate
706 // them from scratch because the input could have been changed in this regard...
707 ElementContext elemCtx(this->simulator_);
708 for (const auto& elem : elements(this->gridView())) {
709 elemCtx.updateStencil(elem);
710 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timIdx=*/0); ++dofIdx) {
711 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timIdx=*/0);
712 updatePvtRegionIndex_(this->solution(/*timeIdx=*/0)[globalDofIdx],
713 elemCtx,
714 dofIdx,
715 /*timeIdx=*/0);
716 }
717 }
718
719 this->solution(/*timeIdx=*/1) = this->solution(/*timeIdx=*/0);
720 }
721
722/*
723 // hack: this interferes with the static polymorphism trick
724protected:
725 friend ParentType;
726 friend Discretization;
727*/
728
729 template <class Context>
731 const Context& context,
732 unsigned dofIdx,
733 unsigned timeIdx)
734 { updatePvtRegionIndex_(priVars, context, dofIdx, timeIdx); }
735
737 {
739
740 // add the VTK output modules which make sense for the blackoil model
741 if constexpr (enableSolvent) {
742 SolventModule::registerOutputModules(asImp_(), this->simulator_);
743 }
744 if constexpr (enablePolymer) {
745 PolymerModule::registerOutputModules(asImp_(), this->simulator_);
746 }
747 if constexpr (enableFullyImplicitThermal) {
748 EnergyModule::registerOutputModules(asImp_(), this->simulator_);
749 }
750 if constexpr (enableBioeffects) {
751 BioeffectsModule::registerOutputModules(asImp_(), this->simulator_);
752 }
753
754 this->addOutputModule(std::make_unique<VtkBlackOilModule<TypeTag>>(this->simulator_));
755 this->addOutputModule(std::make_unique<VtkCompositionModule<TypeTag>>(this->simulator_));
756
757 if constexpr (enableDiffusion) {
758 this->addOutputModule(std::make_unique<VtkDiffusionModule<TypeTag>>(this->simulator_));
759 }
760 }
761
762private:
763 std::vector<Scalar> eqWeights_;
764
765 Implementation& asImp_()
766 { return *static_cast<Implementation*>(this); }
767
768 const Implementation& asImp_() const
769 { return *static_cast<const Implementation*>(this); }
770
771 template <class Context>
772 void updatePvtRegionIndex_(PrimaryVariables& priVars,
773 const Context& context,
774 unsigned dofIdx,
775 unsigned timeIdx)
776 {
777 const unsigned regionIdx = context.problem().pvtRegionIndex(context, dofIdx, timeIdx);
778 priVars.setPvtRegionIndex(regionIdx);
779 }
780};
781
782} // namespace Opm
783
784#endif // OPM_BLACK_OIL_MODEL_HPP
This file contains the default flux module of the blackoil model.
Contains classes extending the black-oil model. \detail This file holds dummy definitions,...
Declares the properties required by the black oil model.
Implements a boundary vector for the fully implicit black-oil model.
Definition: blackoilboundaryratevector.hh:55
Provides the auxiliary methods required for consideration of the diffusion equation.
Provides the auxiliary methods required for consideration of the dispersion equation.
This template class contains the data which is required to calculate the fluxes of the fluid phases o...
Definition: blackoilextensivequantities.hh:57
Contains the quantities which are are constant within a finite volume in the black-oil model.
Definition: blackoilintensivequantities.hh:80
Calculates the local residual of the black oil model.
Definition: blackoillocalresidual.hh:51
A fully-implicit black-oil flow model.
Definition: blackoilmodel.hh:336
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: blackoilmodel.hh:340
BlackOilModel(Simulator &simulator)
Definition: blackoilmodel.hh:377
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: blackoilmodel.hh:506
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: blackoilmodel.hh:570
void supplementInitialSolution_(PrimaryVariables &priVars, const Context &context, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilmodel.hh:730
void serializeEntity(std::ostream &outstream, const DofEntity &dof)
Write the current solution for a degree of freedom to a restart file.
Definition: blackoilmodel.hh:593
void registerOutputModules_()
Definition: blackoilmodel.hh:736
std::string eqName(int eqIdx) const
Given an equation index, return a human readable name.
Definition: blackoilmodel.hh:468
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: blackoilmodel.hh:339
static std::string name()
Definition: blackoilmodel.hh:420
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: blackoilmodel.hh:426
GetPropType< TypeTag, Properties::Indices > Indices
Definition: blackoilmodel.hh:338
void deserializeEntity(std::istream &instream, const DofEntity &dof)
Reads the current solution variables for a degree of freedom from a restart file.
Definition: blackoilmodel.hh:638
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: blackoilmodel.hh:375
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: blackoilmodel.hh:386
void setEqWeight(unsigned eqIdx, Scalar value)
Definition: blackoilmodel.hh:581
void deserialize(Restarter &res)
Deserializes the state of the model.
Definition: blackoilmodel.hh:700
A newton solver which is specific to the black oil model.
Definition: blackoilnewtonmethod.hpp:64
Represents the primary variables used by the black-oil model.
Definition: blackoilprimaryvariables.hh:70
Base class for all problems which use the black-oil model.
Definition: blackoilproblem.hh:43
Implements a vector representing mass, molar or volumetric rates for the black oil model.
Definition: blackoilratevector.hh:59
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:175
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:197
void registerOutputModules_()
Definition: multiphasebasemodel.hh:273
VTK output module for the black oil model's parameters.
Definition: vtkblackoilmodule.hpp:57
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkblackoilmodule.hpp:93
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hpp:57
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hpp:87
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Definition: vtkdiffusionmodule.hpp:58
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hpp:88
PressureMeaning
Definition: blackoilmeanings.hh:29
WaterMeaning
Definition: blackoilmeanings.hh:22
GasMeaning
Definition: blackoilmeanings.hh:35
Definition: blackoilmodel.hh:75
Definition: blackoilbioeffectsmodules.hh:45
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:233
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
The Opm property system, traits with inheritance.
Provides a Darcy flux module for the blackoil model.
Definition: blackoildarcyfluxmodule.hh:49
The primary variable and equation indices for the three-phase black-oil model.
Definition: blackoilvariableandequationindices.hh:49
The type of the base class for all problems which use this model.
Definition: fvbaseproperties.hh:90
Definition: blackoilproperties.hh:100
Similarly to the energy equation, a scaling is applied to the urea equation in MICP.
Definition: blackoilproperties.hh:104
Enable surface volume scaling.
Definition: blackoilproperties.hh:59
Type of object for specifying boundary conditions.
Definition: fvbaseproperties.hh:125
Enable the ECL-blackoil extension for bioeffects (biofilm/MICP)
Definition: blackoilproperties.hh:83
Enable the ECL-blackoil extension for salt.
Definition: blackoilproperties.hh:67
Enable convective mixing?
Definition: multiphasebaseproperties.hh:99
Enable diffusive fluxes?
Definition: multiphasebaseproperties.hh:91
Enable the ECL-blackoil extension for disolution of gas into water.
Definition: blackoilproperties.hh:79
Enable dispersive fluxes?
Definition: multiphasebaseproperties.hh:95
Enable the ECL-blackoil extension for extended BO. ("Second gas" - alternative approach)
Definition: blackoilproperties.hh:47
Enable the ECL-blackoil extension for foam.
Definition: blackoilproperties.hh:63
Definition: blackoilproperties.hh:86
Enable the tracking polymer molecular weight tracking and related functionalities.
Definition: blackoilproperties.hh:55
Enable the ECL-blackoil extension for polymer.
Definition: blackoilproperties.hh:51
Enable the ECL-blackoil extension for salt precipitation.
Definition: blackoilproperties.hh:71
Enable the ECL-blackoil extension for solvents. ("Second gas")
Definition: blackoilproperties.hh:43
Enable the ECL-blackoil extension for water evaporation.
Definition: blackoilproperties.hh:75
Specifies who temperature is modeled by the simulator.
Definition: blackoilproperties.hh:108
Data required to calculate a flux over a face.
Definition: fvbaseproperties.hh:164
GetPropType< TypeTag, Properties::Evaluation > Evaluation
Definition: blackoilmodel.hh:156
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: blackoilmodel.hh:155
BlackOilFluidSystem< Scalar > type
Definition: blackoilmodel.hh:157
The fluid systems including the information about the phases.
Definition: multiphasebaseproperties.hh:79
Specifies the relation used for velocity.
Definition: multiphasebaseproperties.hh:83
Enumerations used by the model.
Definition: multiphasebaseproperties.hh:51
The secondary variables within a sub-control volume.
Definition: fvbaseproperties.hh:139
The type of the local residual function.
Definition: fvbaseproperties.hh:100
The type of the model.
Definition: basicproperties.hh:92
Specifies the type of the actual Newton method.
Definition: newtonmethodproperties.hh:32
A vector of primary variables within a sub-control volume.
Definition: fvbaseproperties.hh:136
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:122
The discretization specific part of the intensive quantities.
Definition: fvbaseproperties.hh:148
The type tag for the black-oil problems.
Definition: blackoilmodel.hh:81
std::tuple< MultiPhaseBaseModel > InheritsFrom
Definition: blackoilmodel.hh:81