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
56
58
62
63#include <cassert>
64#include <istream>
65#include <memory>
66#include <ostream>
67#include <sstream>
68#include <stdexcept>
69#include <string>
70#include <tuple>
71#include <vector>
72
73namespace Opm {
74
75template <class TypeTag>
76class BlackOilModel;
77
78}
79
80namespace Opm::Properties {
81
82namespace TTag {
83
86{ using InheritsFrom = std::tuple<VtkBlackOilPolymer, MultiPhaseBaseModel>; };
87} // namespace TTag
88
90template<class TypeTag>
91struct LocalResidual<TypeTag, TTag::BlackOilModel>
93
95template<class TypeTag>
96struct NewtonMethod<TypeTag, TTag::BlackOilModel>
98
100template<class TypeTag>
101struct Model<TypeTag, TTag::BlackOilModel>
103
105template<class TypeTag>
106struct BaseProblem<TypeTag, TTag::BlackOilModel>
108
110template<class TypeTag>
111struct RateVector<TypeTag, TTag::BlackOilModel>
113
115template<class TypeTag>
116struct BoundaryRateVector<TypeTag, TTag::BlackOilModel>
118
120template<class TypeTag>
121struct PrimaryVariables<TypeTag, TTag::BlackOilModel>
123
125template<class TypeTag>
126struct IntensiveQuantities<TypeTag, TTag::BlackOilModel>
128
130template<class TypeTag>
131struct ExtensiveQuantities<TypeTag, TTag::BlackOilModel>
133
136template<class TypeTag>
137struct FluxModule<TypeTag, TTag::BlackOilModel>
139
141template<class TypeTag>
142struct Indices<TypeTag, TTag::BlackOilModel>
143{
145 getPropValue<TypeTag, Properties::EnableExtbo>(),
146 getPropValue<TypeTag, Properties::EnablePolymer>(),
147 getPropValue<TypeTag, Properties::EnergyModuleType>() == EnergyModules::FullyImplicitThermal,
148 getPropValue<TypeTag, Properties::EnableFoam>(),
149 getPropValue<TypeTag, Properties::EnableBrine>(),
150 /*PVOffset=*/0,
151 getPropValue<TypeTag, Properties::EnableBioeffects>()>;
152};
153
155template<class TypeTag>
156struct FluidSystem<TypeTag, TTag::BlackOilModel>
157{
158public:
161 using type = BlackOilFluidSystem<Scalar>;
162};
163
164// by default, all ECL extension modules are disabled
165template<class TypeTag>
166struct EnableSolvent<TypeTag, TTag::BlackOilModel>
167{ static constexpr bool value = false; };
168
169template<class TypeTag>
170struct EnableExtbo<TypeTag, TTag::BlackOilModel>
171{ static constexpr bool value = false; };
172
173template<class TypeTag>
174struct EnablePolymer<TypeTag, TTag::BlackOilModel>
175{ static constexpr bool value = false; };
176
177template<class TypeTag>
178struct EnablePolymerMW<TypeTag, TTag::BlackOilModel>
179{ static constexpr bool value = false; };
180
181template<class TypeTag>
182struct EnableFoam<TypeTag, TTag::BlackOilModel>
183{ static constexpr bool value = false; };
184
185template<class TypeTag>
186struct EnableBrine<TypeTag, TTag::BlackOilModel>
187{ static constexpr bool value = false; };
188
189template<class TypeTag>
190struct EnableVapwat<TypeTag, TTag::BlackOilModel>
191{ static constexpr bool value = false; };
192
193template<class TypeTag>
194struct EnableDisgasInWater<TypeTag, TTag::BlackOilModel>
195{ static constexpr bool value = false; };
196
197template<class TypeTag>
199{ static constexpr bool value = false; };
200
201template<class TypeTag>
202struct EnableBioeffects<TypeTag, TTag::BlackOilModel>
203{ static constexpr bool value = false; };
204
205template<class TypeTag>
206struct EnergyModuleType<TypeTag, TTag::BlackOilModel>
207{ static constexpr EnergyModules value = EnergyModules::NoTemperature; };
208
210template<class TypeTag>
211struct EnableDiffusion<TypeTag, TTag::BlackOilModel>
212{ static constexpr bool value = false; };
213
215template<class TypeTag>
216struct EnableDispersion<TypeTag, TTag::BlackOilModel>
217{ static constexpr bool value = false; };
218
219template<class TypeTag>
221{ static constexpr bool value = false; };
222
223template<class TypeTag>
224struct EnableMech<TypeTag, TTag::BlackOilModel>
225{ static constexpr bool value = false; };
226
227template<class TypeTag>
228struct RunAssemblyOnGpu<TypeTag, TTag::BlackOilModel>
229{ static constexpr bool value = false; };
230
237template<class TypeTag>
239{
240private:
242 static constexpr Scalar alpha =
243 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>() ? 1000.0 : 1.0;
244
245public:
246 using type = Scalar;
247 static constexpr Scalar value = 1.0/(30.0*4184.0*alpha);
248};
249
251template<class TypeTag>
253{
254private:
256 static constexpr Scalar alpha =
257 getPropValue<TypeTag, Properties::BlackoilConserveSurfaceVolume>() ? 1000.0 : 1.0;
258
259public:
260 using type = Scalar;
261 static constexpr Scalar value = 1.0/(10.0*alpha);
262};
263
264// by default, ebos formulates the conservation equations in terms of mass not surface
265// volumes
266template<class TypeTag>
268{ static constexpr bool value = false; };
269
270} // namespace Opm::Properties
271
272namespace Opm {
273
337template<class TypeTag >
339 : public MultiPhaseBaseModel<TypeTag>
340{
341public:
345
346private:
347 using Implementation = GetPropType<TypeTag, Properties::Model>;
349
354
355 enum { numComponents = FluidSystem::numComponents };
356 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
357 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
358 enum { enableDispersion = getPropValue<TypeTag, Properties::EnableDispersion>() };
359
360 static constexpr bool compositionSwitchEnabled = Indices::compositionSwitchIdx >= 0;
361 static constexpr bool waterEnabled = Indices::waterEnabled;
362
363 using SolventModule = BlackOilSolventModule<TypeTag>;
364 using ExtboModule = BlackOilExtboModule<TypeTag>;
365 using PolymerModule = BlackOilPolymerModule<TypeTag>;
366 using EnergyModule = BlackOilEnergyModule<TypeTag>;
367 using DiffusionModule = BlackOilDiffusionModule<TypeTag, enableDiffusion>;
368 using DispersionModule = BlackOilDispersionModule<TypeTag, enableDispersion>;
369 using BioeffectsModule = BlackOilBioeffectsModule<TypeTag>;
370
371public:
373
374 explicit BlackOilModel(Simulator& simulator)
375 : ParentType(simulator)
376 {
377 eqWeights_.resize(numEq, 1.0);
378 }
379
383 static void registerParameters()
384 {
386
391 DiffusionModule::registerParameters();
393
394 // register runtime parameters of the VTK output modules
398 }
399
403 static std::string name()
404 { return "blackoil"; }
405
409 std::string primaryVarName(int pvIdx) const
410 {
411 if (pvIdx == Indices::waterSwitchIdx) {
412 return "water_switching";
413 }
414 else if (pvIdx == Indices::pressureSwitchIdx) {
415 return "pressure_switching";
416 }
417 else if (static_cast<int>(pvIdx) == Indices::compositionSwitchIdx) {
418 return "composition_switching";
419 }
420 else if (SolventModule::primaryVarApplies(pvIdx)) {
421 return SolventModule::primaryVarName(pvIdx);
422 }
423 else if (ExtboModule::primaryVarApplies(pvIdx)) {
424 return ExtboModule::primaryVarName(pvIdx);
425 }
426 else if (PolymerModule::primaryVarApplies(pvIdx)) {
427 return PolymerModule::primaryVarName(pvIdx);
428 }
429 else if (EnergyModule::primaryVarApplies(pvIdx)) {
430 return EnergyModule::primaryVarName(pvIdx);
431 }
432 else {
433 throw std::logic_error("Invalid primary variable index");
434 }
435 }
436
440 std::string eqName(int eqIdx) const
441 {
442 if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx + numComponents) {
443 std::ostringstream oss;
444 oss << "conti_" << FluidSystem::phaseName(eqIdx - Indices::conti0EqIdx);
445 return oss.str();
446 }
447 else if (SolventModule::eqApplies(eqIdx)) {
448 return SolventModule::eqName(eqIdx);
449 }
450 else if (ExtboModule::eqApplies(eqIdx)) {
451 return ExtboModule::eqName(eqIdx);
452 }
453 else if (PolymerModule::eqApplies(eqIdx)) {
454 return PolymerModule::eqName(eqIdx);
455 }
456 else if (EnergyModule::eqApplies(eqIdx)) {
457 return EnergyModule::eqName(eqIdx);
458 }
459 else {
460 throw std::logic_error("Invalid equation index");
461 }
462 }
463
467 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
468 {
469 // do not care about the auxiliary equations as they are supposed to scale
470 // themselves
471 if (globalDofIdx >= this->numGridDof()) {
472 return 1.0;
473 }
474
475 // saturations are always in the range [0, 1]!
476 if (int(Indices::waterSwitchIdx) == int(pvIdx)) {
477 return 1.0;
478 }
479
480 // oil pressures usually are in the range of 100 to 500 bars for typical oil
481 // reservoirs (which is the only relevant application for the black-oil model).
482 else if (int(Indices::pressureSwitchIdx) == int(pvIdx)) {
483 return 1.0 / 300e5;
484 }
485
486 // deal with primary variables stemming from the solvent module
487 else if (SolventModule::primaryVarApplies(pvIdx)) {
489 }
490
491 // deal with primary variables stemming from the extBO module
492 else if (ExtboModule::primaryVarApplies(pvIdx)) {
493 return ExtboModule::primaryVarWeight(pvIdx);
494 }
495
496 // deal with primary variables stemming from the polymer module
497 else if (PolymerModule::primaryVarApplies(pvIdx)) {
499 }
500
501 // deal with primary variables stemming from the energy module
502 else if (EnergyModule::primaryVarApplies(pvIdx)) {
503 return EnergyModule::primaryVarWeight(pvIdx);
504 }
505
506 // if the primary variable is either the gas saturation, Rs or Rv
507 assert(int(Indices::compositionSwitchIdx) == int(pvIdx));
508
509 switch (this->solution(0)[globalDofIdx].primaryVarsMeaningGas()) {
510 case PrimaryVariables::GasMeaning::Sg: return 1.0; // gas saturation
511 case PrimaryVariables::GasMeaning::Rs: return 1.0 / 250.; // gas dissolution factor
512 case PrimaryVariables::GasMeaning::Rv: return 1.0 / 0.025; // oil vaporization factor
513 default: throw std::logic_error("Invalid primary variable meaning flag for gas");
514 }
515 }
516
523 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
524 {
525 // do not care about the auxiliary equations as they are supposed to scale
526 // themselves
527 if (globalDofIdx >= this->numGridDof()) {
528 return 1.0;
529 }
530
531 return eqWeights_[eqIdx];
532 }
533
534 void setEqWeight(unsigned eqIdx, Scalar value)
535 { eqWeights_[eqIdx] = value; }
536
545 template <class DofEntity>
546 void serializeEntity(std::ostream& outstream, const DofEntity& dof)
547 {
548 const unsigned dofIdx = static_cast<unsigned>(asImp_().dofMapper().index(dof));
549
550 // write phase state
551 if (!outstream.good()) {
552 throw std::runtime_error("Could not serialize degree of freedom " + std::to_string(dofIdx));
553 }
554
555 // write the primary variables
556 const auto& priVars = this->solution(/*timeIdx=*/0)[dofIdx];
557 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
558 outstream << priVars[eqIdx] << " ";
559 }
560
561 // write the pseudo primary variables
562 outstream << static_cast<int>(priVars.primaryVarsMeaningGas()) << " ";
563 outstream << static_cast<int>(priVars.primaryVarsMeaningWater()) << " ";
564 outstream << static_cast<int>(priVars.primaryVarsMeaningPressure()) << " ";
565
566 outstream << priVars.pvtRegionIndex() << " ";
567
568 SolventModule::serializeEntity(asImp_(), outstream, dof);
569 ExtboModule::serializeEntity(asImp_(), outstream, dof);
570 PolymerModule::serializeEntity(asImp_(), outstream, dof);
571 EnergyModule::serializeEntity(asImp_(), outstream, dof);
572 }
573
582 template <class DofEntity>
583 void deserializeEntity(std::istream& instream,
584 const DofEntity& dof)
585 {
586 const unsigned dofIdx = static_cast<unsigned>(asImp_().dofMapper().index(dof));
587
588 // read in the "real" primary variables of the DOF
589 auto& priVars = this->solution(/*timeIdx=*/0)[dofIdx];
590 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
591 if (!instream.good()) {
592 throw std::runtime_error("Could not deserialize degree of freedom " + std::to_string(dofIdx));
593 }
594 instream >> priVars[eqIdx];
595 }
596
597 // read the pseudo primary variables
598 unsigned primaryVarsMeaningGas;
599 instream >> primaryVarsMeaningGas;
600
601 unsigned primaryVarsMeaningWater;
602 instream >> primaryVarsMeaningWater;
603
604 unsigned primaryVarsMeaningPressure;
605 instream >> primaryVarsMeaningPressure;
606
607 unsigned pvtRegionIdx;
608 instream >> pvtRegionIdx;
609
610 if (!instream.good()) {
611 throw std::runtime_error("Could not deserialize degree of freedom " + std::to_string(dofIdx));
612 }
613
614 SolventModule::deserializeEntity(asImp_(), instream, dof);
615 ExtboModule::deserializeEntity(asImp_(), instream, dof);
616 PolymerModule::deserializeEntity(asImp_(), instream, dof);
617 EnergyModule::deserializeEntity(asImp_(), instream, dof);
618
619 using PVM_G = typename PrimaryVariables::GasMeaning;
620 using PVM_W = typename PrimaryVariables::WaterMeaning;
621 using PVM_P = typename PrimaryVariables::PressureMeaning;
622 priVars.setPrimaryVarsMeaningGas(static_cast<PVM_G>(primaryVarsMeaningGas));
623 priVars.setPrimaryVarsMeaningWater(static_cast<PVM_W>(primaryVarsMeaningWater));
624 priVars.setPrimaryVarsMeaningPressure(static_cast<PVM_P>(primaryVarsMeaningPressure));
625
626 priVars.setPvtRegionIndex(pvtRegionIdx);
627 }
628
636 template <class Restarter>
637 void deserialize(Restarter& res)
638 {
639 ParentType::deserialize(res);
640
641 // set the PVT indices of the primary variables. This is also done by writing
642 // them into the restart file and re-reading them, but it is better to calculate
643 // them from scratch because the input could have been changed in this regard...
644 ElementContext elemCtx(this->simulator_);
645 for (const auto& elem : elements(this->gridView())) {
646 elemCtx.updateStencil(elem);
647 for (unsigned dofIdx = 0; dofIdx < elemCtx.numPrimaryDof(/*timIdx=*/0); ++dofIdx) {
648 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, /*timIdx=*/0);
649 updatePvtRegionIndex_(this->solution(/*timeIdx=*/0)[globalDofIdx],
650 elemCtx,
651 dofIdx,
652 /*timeIdx=*/0);
653 }
654 }
655
656 this->solution(/*timeIdx=*/1) = this->solution(/*timeIdx=*/0);
657 }
658
659/*
660 // hack: this interferes with the static polymorphism trick
661protected:
662 friend ParentType;
663 friend Discretization;
664*/
665
666 template <class Context>
668 const Context& context,
669 unsigned dofIdx,
670 unsigned timeIdx)
671 { updatePvtRegionIndex_(priVars, context, dofIdx, timeIdx); }
672
674 {
676
677 // add the VTK output modules which make sense for the blackoil model
678 SolventModule::registerOutputModules(asImp_(), this->simulator_);
679 PolymerModule::registerOutputModules(asImp_(), this->simulator_);
680 EnergyModule::registerOutputModules(asImp_(), this->simulator_);
681 BioeffectsModule::registerOutputModules(asImp_(), this->simulator_);
682
683 this->addOutputModule(std::make_unique<VtkBlackOilModule<TypeTag>>(this->simulator_));
684 this->addOutputModule(std::make_unique<VtkCompositionModule<TypeTag>>(this->simulator_));
685
686 if constexpr (enableDiffusion) {
687 this->addOutputModule(std::make_unique<VtkDiffusionModule<TypeTag>>(this->simulator_));
688 }
689 }
690
691private:
692 std::vector<Scalar> eqWeights_;
693
694 Implementation& asImp_()
695 { return *static_cast<Implementation*>(this); }
696
697 const Implementation& asImp_() const
698 { return *static_cast<const Implementation*>(this); }
699
700 template <class Context>
701 void updatePvtRegionIndex_(PrimaryVariables& priVars,
702 const Context& context,
703 unsigned dofIdx,
704 unsigned timeIdx)
705 {
706 const unsigned regionIdx = context.problem().pvtRegionIndex(context, dofIdx, timeIdx);
707 priVars.setPvtRegionIndex(regionIdx);
708 }
709};
710
711} // namespace Opm
712
713#endif // OPM_BLACK_OIL_MODEL_HPP
Contains the classes required to extend the black-oil model by bioeffects.
Contains the classes required to extend the black-oil model by brine.
This file contains the default flux module of the blackoil model.
Classes required for molecular diffusion.
Classes required for mechanical dispersion.
Contains the classes required to extend the black-oil model by energy.
Contains the classes required to extend the black-oil model by solvent component. For details,...
Contains the classes required to extend the black-oil model to include the effects of foam.
Contains the classes required to extend the black-oil model by polymer.
Declares the properties required by the black oil model.
Contains the classes required to extend the black-oil model by solvents.
static void registerParameters()
Register all run-time parameters for the black-oil bioeffects module.
Definition: blackoilbioeffectsmodules.hh:141
static void registerOutputModules(Model &model, Simulator &simulator)
Register all bioeffects specific VTK and ECL output modules.
Definition: blackoilbioeffectsmodules.hh:150
Implements a boundary vector for the fully implicit black-oil model.
Definition: blackoilboundaryratevector.hh:48
static std::string eqName(unsigned eqIdx)
Definition: blackoilenergymodules.hh:142
static bool eqApplies(unsigned eqIdx)
Definition: blackoilenergymodules.hh:132
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilenergymodules.hh:117
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilenergymodules.hh:316
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilenergymodules.hh:124
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilenergymodules.hh:326
static void registerOutputModules(Model &model, Simulator &simulator)
Register all energy specific VTK and ECL output modules.
Definition: blackoilenergymodules.hh:99
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilenergymodules.hh:107
static void registerParameters()
Register all run-time parameters for the black-oil energy module.
Definition: blackoilenergymodules.hh:89
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilextbomodules.hh:107
static std::string eqName(unsigned eqIdx)
Definition: blackoilextbomodules.hh:142
static bool eqApplies(unsigned eqIdx)
Definition: blackoilextbomodules.hh:132
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilextbomodules.hh:124
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilextbomodules.hh:293
static void registerParameters()
Register all run-time parameters for the black-oil solvent module.
Definition: blackoilextbomodules.hh:97
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilextbomodules.hh:304
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilextbomodules.hh:117
This template class contains the data which is required to calculate the fluxes of the fluid phases o...
Definition: blackoilextensivequantities.hh:59
Contains the quantities which are are constant within a finite volume in the black-oil model.
Definition: blackoilintensivequantities.hh:87
Calculates the local residual of the black oil model.
Definition: blackoillocalresidual.hh:56
A fully-implicit black-oil flow model.
Definition: blackoilmodel.hh:340
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: blackoilmodel.hh:344
BlackOilModel(Simulator &simulator)
Definition: blackoilmodel.hh:374
std::string primaryVarName(int pvIdx) const
Given an primary variable index, return a human readable name.
Definition: blackoilmodel.hh:409
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: blackoilmodel.hh:467
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: blackoilmodel.hh:523
void supplementInitialSolution_(PrimaryVariables &priVars, const Context &context, unsigned dofIdx, unsigned timeIdx)
Definition: blackoilmodel.hh:667
void serializeEntity(std::ostream &outstream, const DofEntity &dof)
Write the current solution for a degree of freedom to a restart file.
Definition: blackoilmodel.hh:546
void registerOutputModules_()
Definition: blackoilmodel.hh:673
std::string eqName(int eqIdx) const
Given an equation index, return a human readable name.
Definition: blackoilmodel.hh:440
GetPropType< TypeTag, Properties::FluidSystem > FluidSystem
Definition: blackoilmodel.hh:343
static std::string name()
Definition: blackoilmodel.hh:403
GetPropType< TypeTag, Properties::Indices > Indices
Definition: blackoilmodel.hh:342
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:583
GetPropType< TypeTag, Properties::LocalResidual > LocalResidual
Definition: blackoilmodel.hh:372
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: blackoilmodel.hh:383
void setEqWeight(unsigned eqIdx, Scalar value)
Definition: blackoilmodel.hh:534
void deserialize(Restarter &res)
Deserializes the state of the model.
Definition: blackoilmodel.hh:637
A newton solver which is specific to the black oil model.
Definition: blackoilnewtonmethod.hpp:61
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:180
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilpolymermodules.hh:342
static void registerParameters()
Register all run-time parameters for the black-oil polymer module.
Definition: blackoilpolymermodules.hh:147
static bool eqApplies(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:200
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:192
static void registerOutputModules(Model &model, Simulator &simulator)
Register all polymer specific VTK and ECL output modules.
Definition: blackoilpolymermodules.hh:157
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilpolymermodules.hh:165
static std::string eqName(unsigned eqIdx)
Definition: blackoilpolymermodules.hh:215
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilpolymermodules.hh:353
Represents the primary variables used by the black-oil model.
Definition: blackoilprimaryvariables.hh:72
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:62
static Scalar primaryVarWeight(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:151
static std::string eqName(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:169
static void serializeEntity(const Model &model, std::ostream &outstream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:341
static void registerOutputModules(Model &model, Simulator &simulator)
Register all solvent specific VTK and ECL output modules.
Definition: blackoilsolventmodules.hh:126
static void registerParameters()
Register all run-time parameters for the black-oil solvent module.
Definition: blackoilsolventmodules.hh:116
static std::string primaryVarName(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:144
static bool eqApplies(unsigned eqIdx)
Definition: blackoilsolventmodules.hh:159
static void deserializeEntity(Model &model, std::istream &instream, const DofEntity &dof)
Definition: blackoilsolventmodules.hh:352
static bool primaryVarApplies(unsigned pvIdx)
Definition: blackoilsolventmodules.hh:134
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:80
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)
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:47
The type of the base class for all problems which use this model.
Definition: fvbaseproperties.hh:84
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:119
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:153
GetPropType< TypeTag, Properties::Evaluation > Evaluation
Definition: blackoilmodel.hh:160
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: blackoilmodel.hh:159
BlackOilFluidSystem< Scalar > type
Definition: blackoilmodel.hh:161
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:133
The type of the local residual function.
Definition: fvbaseproperties.hh:94
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:130
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:116
The discretization specific part of the intensive quantities.
Definition: fvbaseproperties.hh:137
The type tag for the black-oil problems.
Definition: blackoilmodel.hh:86
std::tuple< VtkBlackOilPolymer, MultiPhaseBaseModel > InheritsFrom
Definition: blackoilmodel.hh:86