28#ifndef EWOMS_PVS_MODEL_HH
29#define EWOMS_PVS_MODEL_HH
31#include <opm/common/Exceptions.hpp>
33#include <opm/material/densead/Math.hpp>
34#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
35#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
68template <
class TypeTag>
73namespace Opm::Properties {
79{
using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
84template<
class TypeTag>
89template<
class TypeTag>
94template<
class TypeTag>
99template<
class TypeTag>
104template<
class TypeTag>
109template<
class TypeTag>
114template<
class TypeTag>
119template<
class TypeTag>
124template<
class TypeTag>
129template<
class TypeTag>
131{
static constexpr bool value =
false; };
134template<
class TypeTag>
136{
static constexpr bool value =
false; };
139template<
class TypeTag>
143 static constexpr type value = 1.0;
147template<
class TypeTag>
151 static constexpr type value = 1.0;
155template<
class TypeTag>
159 static constexpr type value = 1.0;
164namespace Opm::Parameters {
168{
static constexpr int value = 1; };
270template <
class TypeTag>
272 :
public MultiPhaseBaseModel<TypeTag>
274 using ParentType = MultiPhaseBaseModel<TypeTag>;
291 using Element =
typename GridView::template Codim<0>::Entity;
292 using ElementIterator =
typename GridView::template Codim<0>::Iterator;
297 explicit PvsModel(Simulator& simulator)
298 : ParentType(simulator)
315 if constexpr (enableDiffusion) {
319 if constexpr (enableEnergy) {
324 (
"The verbosity level of the primary variable "
339 const std::string s = EnergyModule::primaryVarName(pvIdx);
344 std::ostringstream oss;
345 if (pvIdx == Indices::pressure0Idx) {
346 oss <<
"pressure_" << FluidSystem::phaseName(0);
348 else if (Indices::switch0Idx <= pvIdx &&
349 pvIdx < Indices::switch0Idx + numPhases - 1)
351 oss <<
"switch_" << pvIdx - Indices::switch0Idx;
353 else if (Indices::switch0Idx + numPhases - 1 <= pvIdx &&
354 pvIdx < Indices::switch0Idx + numComponents - 1)
356 oss <<
"auxMoleFrac^" << FluidSystem::componentName(pvIdx);
369 const std::string s = EnergyModule::eqName(eqIdx);
374 std::ostringstream oss;
375 if (Indices::conti0EqIdx <= eqIdx &&
376 eqIdx < Indices::conti0EqIdx + numComponents)
378 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
379 oss <<
"continuity^" << FluidSystem::componentName(compIdx);
393 ParentType::updateFailed();
402 ParentType::updateBegin();
407 const std::size_t nDof = this->numTotalDof();
408 for (
unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx) {
409 if (this->dofTotalVolume(dofIdx) > 0.0) {
411 this->solution(0)[dofIdx][Indices::pressure0Idx];
412 if (referencePressure_ > 0.0) {
424 const Scalar tmp = EnergyModule::primaryVarWeight(*
this, globalDofIdx, pvIdx);
430 if (Indices::pressure0Idx == pvIdx) {
431 return 10 / referencePressure_;
434 if (Indices::switch0Idx <= pvIdx &&
435 pvIdx < Indices::switch0Idx + numPhases - 1)
437 const unsigned phaseIdx = pvIdx - Indices::switch0Idx;
439 if (!this->solution(0)[globalDofIdx].phaseIsPresent(phaseIdx)) {
460 Scalar
eqWeight(
unsigned globalDofIdx,
unsigned eqIdx)
const
462 const Scalar tmp = EnergyModule::eqWeight(*
this, globalDofIdx, eqIdx);
468 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
469 assert(compIdx <= numComponents);
472 return FluidSystem::molarMass(compIdx);
480 ParentType::advanceTimeLevel();
489 {
return numSwitched_ > 0; }
497 template <
class DofEntity>
501 ParentType::serializeEntity(outstream, dofEntity);
503 const unsigned dofIdx =
static_cast<unsigned>(this->dofMapper().index(dofEntity));
504 if (!outstream.good()) {
505 throw std::runtime_error(
"Could not serialize DOF " + std::to_string(dofIdx));
508 outstream << this->solution(0)[dofIdx].phasePresence() <<
" ";
517 template <
class DofEntity>
521 ParentType::deserializeEntity(instream, dofEntity);
524 const unsigned dofIdx =
static_cast<unsigned>(this->dofMapper().index(dofEntity));
525 if (!instream.good()) {
526 throw std::runtime_error(
"Could not deserialize DOF " + std::to_string(dofIdx));
531 this->solution(0)[dofIdx].setPhasePresence(tmp);
532 this->solution(1)[dofIdx].setPhasePresence(tmp);
542 void switchPrimaryVars_()
548 std::vector<bool> visited(this->numGridDof(),
false);
549 ElementContext elemCtx(this->simulator_);
551 for (
const auto& elem : elements(this->gridView_, Dune::Partitions::interior)) {
552 elemCtx.updateStencil(elem);
554 const std::size_t numLocalDof = elemCtx.stencil(0).numPrimaryDof();
555 for (
unsigned dofIdx = 0; dofIdx < numLocalDof; ++dofIdx) {
556 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
558 if (visited[globalIdx]) {
561 visited[globalIdx] =
true;
564 auto& priVars = this->solution(0)[globalIdx];
565 elemCtx.updateIntensiveQuantities(priVars, dofIdx, 0);
566 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, 0);
569 const short oldPhasePresence = priVars.phasePresence();
573 priVars.assignNaive(intQuants.fluidState());
575 if (oldPhasePresence != priVars.phasePresence()) {
576 if (verbosity_ > 1) {
577 printSwitchedPhases_(elemCtx,
579 intQuants.fluidState(),
592 std::cout <<
"rank " << this->simulator_.gridView().comm().rank()
593 <<
" caught an exception during primary variable switching"
594 <<
"\n" << std::flush;
597 succeeded = this->simulator_.gridView().comm().min(succeeded);
600 throw NumericalProblem(
"A process did not succeed in adapting the primary variables");
606 numSwitched_ = this->gridView_.comm().sum(numSwitched_);
608 if (verbosity_ > 0) {
609 this->simulator_.model().newtonMethod().endIterMsg()
610 <<
", num switched=" << numSwitched_;
614 template <
class Flu
idState>
615 void printSwitchedPhases_(
const ElementContext& elemCtx,
617 const FluidState& fs,
618 short oldPhasePresence,
619 const PrimaryVariables& newPv)
const
621 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
623 for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
624 const bool oldPhasePresent = (oldPhasePresence & (1 << phaseIdx)) > 0;
625 const bool newPhasePresent = newPv.phaseIsPresent(phaseIdx);
626 if (oldPhasePresent == newPhasePresent) {
630 const auto& pos = elemCtx.pos(dofIdx, 0);
631 if (oldPhasePresent) {
632 std::cout <<
"'" << FluidSystem::phaseName(phaseIdx)
633 <<
"' phase disappears at position " << pos
634 <<
". saturation=" << fs.saturation(phaseIdx)
639 for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
640 sumx += FsToolbox::value(fs.moleFraction(phaseIdx, compIdx));
643 std::cout <<
"'" << FluidSystem::phaseName(phaseIdx)
644 <<
"' phase appears at position " << pos
645 <<
" sum x = " << sumx << std::flush;
649 std::cout <<
", new primary variables: ";
650 newPv.print(std::cout);
651 std::cout <<
"\n" << std::flush;
654 void registerOutputModules_()
656 ParentType::registerOutputModules_();
659 this->addOutputModule(std::make_unique<VtkPhasePresenceModule<TypeTag>>(this->simulator_));
660 this->addOutputModule(std::make_unique<VtkCompositionModule<TypeTag>>(this->simulator_));
661 if constexpr (enableDiffusion) {
662 this->addOutputModule(std::make_unique<VtkDiffusionModule<TypeTag>>(this->simulator_));
664 if constexpr (enableEnergy) {
665 this->addOutputModule(std::make_unique<VtkEnergyModule<TypeTag>>(this->simulator_));
669 Scalar referencePressure_{};
673 unsigned numSwitched_;
Provides the auxiliary methods required for consideration of the energy equation.
Definition energymodule.hh:54
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition multiphasebasemodel.hh:190
Implements a rate vector on the boundary for the fully implicit compositional multi-phase primary var...
Definition pvsboundaryratevector.hh:51
Contains all data which is required to calculate all fluxes at a flux integration point for the prima...
Definition pvsextensivequantities.hh:54
The indices for the compositional multi-phase primary variable switching model.
Definition pvsindices.hh:48
Contains the quantities which are are constant within a finite volume in the compositional multi-phas...
Definition pvsintensivequantities.hh:62
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Definition pvslocalresidual.hh:49
A generic compositional multi-phase model using primary-variable switching.
Definition pvsmodel.hh:273
void advanceTimeLevel()
Called by the problem if a time integration was successful, post processing of the solution is done a...
Definition pvsmodel.hh:478
void serializeEntity(std::ostream &outstream, const DofEntity &dofEntity)
Write the current solution for a degree of freedom to a restart file.
Definition pvsmodel.hh:498
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition pvsmodel.hh:367
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition pvsmodel.hh:337
void updateBegin()
Called by the update() method before it tries to apply the newton method.
Definition pvsmodel.hh:400
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition pvsmodel.hh:460
static std::string name()
Definition pvsmodel.hh:331
bool switched() const
Return true if the primary variables were switched for at least one vertex after the last timestep.
Definition pvsmodel.hh:488
static void registerParameters()
Register all run-time parameters for the PVS compositional model.
Definition pvsmodel.hh:307
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition pvsmodel.hh:422
void deserializeEntity(std::istream &instream, const DofEntity &dofEntity)
Reads the current solution variables for a degree of freedom from a restart file.
Definition pvsmodel.hh:518
void updateFailed()
Called by the update() method if it was unsuccessful.
Definition pvsmodel.hh:391
A newton solver which is specific to the compositional multi-phase PVS model.
Definition pvsnewtonmethod.hh:54
Represents the primary variables used in the primary variable switching compositional model.
Definition pvsprimaryvariables.hh:62
Implements a vector representing molar, mass or volumetric rates.
Definition pvsratevector.hh:53
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkcompositionmodule.hpp:87
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkdiffusionmodule.hpp:88
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkenergymodule.hpp:87
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtkphasepresencemodule.hpp:71
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:81
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
void Register(const char *usageString)
Register a run-time parameter.
Definition parametersystem.hpp:292
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187
Implements a rate vector on the boundary for the fully implicit compositional multi-phase primary var...
Contains all data which is required to calculate all fluxes at a flux integration point for the prima...
The indices for the compositional multi-phase primary variable switching model.
Contains the quantities which are are constant within a finite volume in the compositional multi-phas...
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
A newton solver which is specific to the compositional multi-phase PVS model.
Represents the primary variables used in the primary variable switching compositional model.
Declares the properties required for the compositional multi-phase primary variable switching model.
Implements a vector representing molar, mass or volumetric rates.
The verbosity of the model (0 -> do not print anything, 2 -> spam stdout a lot).
Definition pvsmodel.hh:168
Type of object for specifying boundary conditions.
Definition fvbaseproperties.hh:119
Enable diffusive fluxes?
Definition multiphasebaseproperties.hh:91
Specify whether energy should be considered as a conservation quantity or not.
Definition multiphasebaseproperties.hh:87
Data required to calculate a flux over a face.
Definition fvbaseproperties.hh:149
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:88
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
The basis value for the weight of the mole fraction primary variables.
Definition pvsproperties.hh:45
The basis value for the weight of the pressure primary variable.
Definition pvsproperties.hh:39
The basis value for the weight of the saturation primary variables.
Definition pvsproperties.hh:42
Vector containing volumetric or areal rates of quantities.
Definition fvbaseproperties.hh:116
The type tag for the isothermal single phase problems.
Definition pvsmodel.hh:79
VTK output module for the fluid composition.
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
VTK output module for quantities which make sense for models which assume thermal equilibrium.
VTK output module for the fluid composition.