pvsmodel.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_PVS_MODEL_HH
29#define EWOMS_PVS_MODEL_HH
30
31#include <opm/common/Exceptions.hpp>
32
33#include <opm/material/densead/Math.hpp>
34#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
35#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
36
40
45
55
56#include <cassert>
57#include <cstddef>
58#include <iostream>
59#include <memory>
60#include <sstream>
61#include <stdexcept>
62#include <string>
63#include <tuple>
64#include <vector>
65
66namespace Opm {
67
68template <class TypeTag>
69class PvsModel;
70
71}
72
73namespace Opm::Properties {
74
75namespace TTag {
76
79{ using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
80
81} // namespace TTag
82
84template<class TypeTag>
85struct LocalResidual<TypeTag, TTag::PvsModel>
87
89template<class TypeTag>
90struct NewtonMethod<TypeTag, TTag::PvsModel>
92
94template<class TypeTag>
95struct Model<TypeTag, TTag::PvsModel>
97
99template<class TypeTag>
100struct PrimaryVariables<TypeTag, TTag::PvsModel>
102
104template<class TypeTag>
105struct RateVector<TypeTag, TTag::PvsModel>
107
109template<class TypeTag>
110struct BoundaryRateVector<TypeTag, TTag::PvsModel>
112
114template<class TypeTag>
115struct IntensiveQuantities<TypeTag, TTag::PvsModel>
117
119template<class TypeTag>
120struct ExtensiveQuantities<TypeTag, TTag::PvsModel>
122
124template<class TypeTag>
125struct Indices<TypeTag, TTag::PvsModel>
126{ using type = Opm::PvsIndices<TypeTag, /*PVIdx=*/0>; };
127
129template<class TypeTag>
130struct EnableEnergy<TypeTag, TTag::PvsModel>
131{ static constexpr bool value = false; };
132
133// disable molecular diffusion by default
134template<class TypeTag>
135struct EnableDiffusion<TypeTag, TTag::PvsModel>
136{ static constexpr bool value = false; };
137
139template<class TypeTag>
140struct PvsPressureBaseWeight<TypeTag, TTag::PvsModel>
141{
143 static constexpr type value = 1.0;
144};
145
147template<class TypeTag>
148struct PvsSaturationsBaseWeight<TypeTag, TTag::PvsModel>
149{
151 static constexpr type value = 1.0;
152};
153
155template<class TypeTag>
157{
159 static constexpr type value = 1.0;
160};
161
162} // namespace Opm::Properties
163
164namespace Opm::Parameters {
165
168{ static constexpr int value = 1; };
169
170} // namespace Opm::Parameters
171
172namespace Opm {
173
269template <class TypeTag>
271 : public MultiPhaseBaseModel<TypeTag>
272{
273 using ParentType = MultiPhaseBaseModel<TypeTag>;
274
279
284
285 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
286 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
287 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
288 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
289
290 using Element = typename GridView::template Codim<0>::Entity;
291 using ElementIterator = typename GridView::template Codim<0>::Iterator;
292
294
295public:
296 explicit PvsModel(Simulator& simulator)
297 : ParentType(simulator)
298 {
299 verbosity_ = Parameters::Get<Parameters::PvsVerbosity>();
300 numSwitched_ = 0;
301 }
302
306 static void registerParameters()
307 {
309
310 // register runtime parameters of the VTK output modules
313
314 if constexpr (enableDiffusion) {
316 }
317
318 if constexpr (enableEnergy) {
320 }
321
322 Parameters::Register<Parameters::PvsVerbosity>
323 ("The verbosity level of the primary variable "
324 "switching model");
325 }
326
330 static std::string name()
331 { return "pvs"; }
332
336 std::string primaryVarName(unsigned pvIdx) const
337 {
338 const std::string s = EnergyModule::primaryVarName(pvIdx);
339 if (!s.empty()) {
340 return s;
341 }
342
343 std::ostringstream oss;
344 if (pvIdx == Indices::pressure0Idx) {
345 oss << "pressure_" << FluidSystem::phaseName(/*phaseIdx=*/0);
346 }
347 else if (Indices::switch0Idx <= pvIdx &&
348 pvIdx < Indices::switch0Idx + numPhases - 1)
349 {
350 oss << "switch_" << pvIdx - Indices::switch0Idx;
351 }
352 else if (Indices::switch0Idx + numPhases - 1 <= pvIdx &&
353 pvIdx < Indices::switch0Idx + numComponents - 1)
354 {
355 oss << "auxMoleFrac^" << FluidSystem::componentName(pvIdx);
356 } else {
357 assert(false);
358 }
359
360 return oss.str();
361 }
362
366 std::string eqName(unsigned eqIdx) const
367 {
368 const std::string s = EnergyModule::eqName(eqIdx);
369 if (!s.empty()) {
370 return s;
371 }
372
373 std::ostringstream oss;
374 if (Indices::conti0EqIdx <= eqIdx &&
375 eqIdx < Indices::conti0EqIdx + numComponents)
376 {
377 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
378 oss << "continuity^" << FluidSystem::componentName(compIdx);
379 }
380 else {
381 assert(false);
382 }
383
384 return oss.str();
385 }
386
391 {
392 ParentType::updateFailed();
393 numSwitched_ = 0;
394 }
395
400 {
401 ParentType::updateBegin();
402
403 // find the a reference pressure. The first degree of freedom
404 // might correspond to non-interior entities which would lead
405 // to an undefined value, so we have to iterate...
406 const std::size_t nDof = this->numTotalDof();
407 for (unsigned dofIdx = 0; dofIdx < nDof; ++dofIdx) {
408 if (this->dofTotalVolume(dofIdx) > 0.0) {
410 this->solution(/*timeIdx=*/0)[dofIdx][/*pvIdx=*/Indices::pressure0Idx];
411 if (referencePressure_ > 0.0) {
412 break;
413 }
414 }
415 }
416 }
417
421 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
422 {
423 const Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
424 if (tmp > 0) {
425 // energy related quantity
426 return tmp;
427 }
428
429 if (Indices::pressure0Idx == pvIdx) {
430 return 10 / referencePressure_;
431 }
432
433 if (Indices::switch0Idx <= pvIdx &&
434 pvIdx < Indices::switch0Idx + numPhases - 1)
435 {
436 const unsigned phaseIdx = pvIdx - Indices::switch0Idx;
437
438 if (!this->solution(/*timeIdx=*/0)[globalDofIdx].phaseIsPresent(phaseIdx)) {
439 // for saturations, the weight is always 1
440 return 1;
441 }
442
443 // for saturations, the PvsMoleSaturationsBaseWeight
444 // property determines the weight
445 return getPropValue<TypeTag, Properties::PvsSaturationsBaseWeight>();
446 }
447
448 // for mole fractions, the PvsMoleFractionsBaseWeight
449 // property determines the weight
450 return getPropValue<TypeTag, Properties::PvsMoleFractionsBaseWeight>();
451 }
452
459 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
460 {
461 const Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
462 if (tmp > 0) {
463 // energy related equation
464 return tmp;
465 }
466
467 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
468 assert(compIdx <= numComponents);
469
470 // make all kg equal
471 return FluidSystem::molarMass(compIdx);
472 }
473
478 {
479 ParentType::advanceTimeLevel();
480 numSwitched_ = 0;
481 }
482
487 bool switched() const
488 { return numSwitched_ > 0; }
489
496 template <class DofEntity>
497 void serializeEntity(std::ostream& outstream, const DofEntity& dofEntity)
498 {
499 // write primary variables
500 ParentType::serializeEntity(outstream, dofEntity);
501
502 const unsigned dofIdx = static_cast<unsigned>(this->dofMapper().index(dofEntity));
503 if (!outstream.good()) {
504 throw std::runtime_error("Could not serialize DOF " + std::to_string(dofIdx));
505 }
506
507 outstream << this->solution(/*timeIdx=*/0)[dofIdx].phasePresence() << " ";
508 }
509
516 template <class DofEntity>
517 void deserializeEntity(std::istream& instream, const DofEntity& dofEntity)
518 {
519 // read primary variables
520 ParentType::deserializeEntity(instream, dofEntity);
521
522 // read phase presence
523 const unsigned dofIdx = static_cast<unsigned>(this->dofMapper().index(dofEntity));
524 if (!instream.good()) {
525 throw std::runtime_error("Could not deserialize DOF " + std::to_string(dofIdx));
526 }
527
528 short tmp;
529 instream >> tmp;
530 this->solution(/*timeIdx=*/0)[dofIdx].setPhasePresence(tmp);
531 this->solution(/*timeIdx=*/1)[dofIdx].setPhasePresence(tmp);
532 }
533
542 {
543 numSwitched_ = 0;
544
545 int succeeded;
546 try {
547 std::vector<bool> visited(this->numGridDof(), false);
548 ElementContext elemCtx(this->simulator_);
549
550 for (const auto& elem : elements(this->gridView_, Dune::Partitions::interior)) {
551 elemCtx.updateStencil(elem);
552
553 const std::size_t numLocalDof = elemCtx.stencil(/*timeIdx=*/0).numPrimaryDof();
554 for (unsigned dofIdx = 0; dofIdx < numLocalDof; ++dofIdx) {
555 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
556
557 if (visited[globalIdx]) {
558 continue;
559 }
560 visited[globalIdx] = true;
561
562 // compute the intensive quantities of the current degree of freedom
563 auto& priVars = this->solution(/*timeIdx=*/0)[globalIdx];
564 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
565 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
566
567 // evaluate primary variable switch
568 const short oldPhasePresence = priVars.phasePresence();
569
570 // set the primary variables and the new phase state
571 // from the current fluid state
572 priVars.assignNaive(intQuants.fluidState());
573
574 if (oldPhasePresence != priVars.phasePresence()) {
575 if (verbosity_ > 1) {
576 printSwitchedPhases_(elemCtx,
577 dofIdx,
578 intQuants.fluidState(),
579 oldPhasePresence,
580 priVars);
581 }
582 ++numSwitched_;
583 }
584 }
585 }
586
587 succeeded = 1;
588 }
589 catch (...)
590 {
591 std::cout << "rank " << this->simulator_.gridView().comm().rank()
592 << " caught an exception during primary variable switching"
593 << "\n" << std::flush;
594 succeeded = 0;
595 }
596 succeeded = this->simulator_.gridView().comm().min(succeeded);
597
598 if (!succeeded) {
599 throw NumericalProblem("A process did not succeed in adapting the primary variables");
600 }
601
602 // make sure that if there was a variable switch in an
603 // other partition we will also set the switch flag
604 // for our partition.
605 numSwitched_ = this->gridView_.comm().sum(numSwitched_);
606
607 if (verbosity_ > 0) {
608 this->simulator_.model().newtonMethod().endIterMsg()
609 << ", num switched=" << numSwitched_;
610 }
611 }
612
613 template <class FluidState>
614 void printSwitchedPhases_(const ElementContext& elemCtx,
615 unsigned dofIdx,
616 const FluidState& fs,
617 short oldPhasePresence,
618 const PrimaryVariables& newPv) const
619 {
620 using FsToolbox = MathToolbox<typename FluidState::ValueType>;
621
622 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
623 const bool oldPhasePresent = (oldPhasePresence & (1 << phaseIdx)) > 0;
624 const bool newPhasePresent = newPv.phaseIsPresent(phaseIdx);
625 if (oldPhasePresent == newPhasePresent) {
626 continue;
627 }
628
629 const auto& pos = elemCtx.pos(dofIdx, /*timeIdx=*/0);
630 if (oldPhasePresent) {
631 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
632 << "' phase disappears at position " << pos
633 << ". saturation=" << fs.saturation(phaseIdx)
634 << std::flush;
635 }
636 else {
637 Scalar sumx = 0;
638 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
639 sumx += FsToolbox::value(fs.moleFraction(phaseIdx, compIdx));
640 }
641
642 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
643 << "' phase appears at position " << pos
644 << " sum x = " << sumx << std::flush;
645 }
646 }
647
648 std::cout << ", new primary variables: ";
649 newPv.print(std::cout);
650 std::cout << "\n" << std::flush;
651 }
652
654 {
656
657 // add the VTK output modules which are meaningful for the model
658 this->addOutputModule(std::make_unique<VtkPhasePresenceModule<TypeTag>>(this->simulator_));
659 this->addOutputModule(std::make_unique<VtkCompositionModule<TypeTag>>(this->simulator_));
660 if constexpr (enableDiffusion) {
661 this->addOutputModule(std::make_unique<VtkDiffusionModule<TypeTag>>(this->simulator_));
662 }
663 if constexpr (enableEnergy) {
664 this->addOutputModule(std::make_unique<VtkEnergyModule<TypeTag>>(this->simulator_));
665 }
666 }
667
669
670 // number of switches of the phase state in the last Newton
671 // iteration
672 unsigned numSwitched_;
673
674 // verbosity of the model
676};
677
678} // namespace Opm
679
680#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
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
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:272
void advanceTimeLevel()
Called by the problem if a time integration was successful, post processing of the solution is done a...
Definition: pvsmodel.hh:477
Scalar referencePressure_
Definition: pvsmodel.hh:668
void serializeEntity(std::ostream &outstream, const DofEntity &dofEntity)
Write the current solution for a degree of freedom to a restart file.
Definition: pvsmodel.hh:497
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: pvsmodel.hh:366
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: pvsmodel.hh:336
void switchPrimaryVars_()
Definition: pvsmodel.hh:541
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: pvsmodel.hh:399
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: pvsmodel.hh:459
int verbosity_
Definition: pvsmodel.hh:675
static std::string name()
Definition: pvsmodel.hh:330
bool switched() const
Return true if the primary variables were switched for at least one vertex after the last timestep.
Definition: pvsmodel.hh:487
static void registerParameters()
Register all run-time parameters for the PVS compositional model.
Definition: pvsmodel.hh:306
void printSwitchedPhases_(const ElementContext &elemCtx, unsigned dofIdx, const FluidState &fs, short oldPhasePresence, const PrimaryVariables &newPv) const
Definition: pvsmodel.hh:614
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: pvsmodel.hh:421
void registerOutputModules_()
Definition: pvsmodel.hh:653
unsigned numSwitched_
Definition: pvsmodel.hh:672
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:517
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: pvsmodel.hh:390
PvsModel(Simulator &simulator)
Definition: pvsmodel.hh:296
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
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
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hpp:58
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hpp:87
VTK output module for the fluid composition.
Definition: vtkphasepresencemodule.hpp:48
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.
Definition: blackoilnewtonmethodparams.hpp:31
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)
Declares the properties required for the compositional multi-phase primary variable switching model.
The verbosity of the model (0 -> do not print anything, 2 -> spam stdout a lot)
Definition: pvsmodel.hh:168
static constexpr int value
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:153
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
GetPropType< TypeTag, Scalar > type
Definition: pvsmodel.hh:158
The basis value for the weight of the mole fraction primary variables.
Definition: pvsproperties.hh:45
GetPropType< TypeTag, Scalar > type
Definition: pvsmodel.hh:142
The basis value for the weight of the pressure primary variable.
Definition: pvsproperties.hh:39
GetPropType< TypeTag, Scalar > type
Definition: pvsmodel.hh:150
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
std::tuple< MultiPhaseBaseModel > InheritsFrom
Definition: pvsmodel.hh:79