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
270template <class TypeTag>
272 : public MultiPhaseBaseModel<TypeTag>
273{
274 using ParentType = MultiPhaseBaseModel<TypeTag>;
275
280
285
286 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
287 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
288 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
289 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
290
291 using Element = typename GridView::template Codim<0>::Entity;
292 using ElementIterator = typename GridView::template Codim<0>::Iterator;
293
295
296public:
297 explicit PvsModel(Simulator& simulator)
298 : ParentType(simulator)
299 {
300 verbosity_ = Parameters::Get<Parameters::PvsVerbosity>();
301 numSwitched_ = 0;
302 }
303
307 static void registerParameters()
308 {
310
311 // register runtime parameters of the VTK output modules
314
315 if constexpr (enableDiffusion) {
317 }
318
319 if constexpr (enableEnergy) {
321 }
322
323 Parameters::Register<Parameters::PvsVerbosity>
324 ("The verbosity level of the primary variable "
325 "switching model");
326 }
327
331 static std::string name()
332 { return "pvs"; }
333
337 std::string primaryVarName(unsigned pvIdx) const
338 {
339 const std::string s = EnergyModule::primaryVarName(pvIdx);
340 if (!s.empty()) {
341 return s;
342 }
343
344 std::ostringstream oss;
345 if (pvIdx == Indices::pressure0Idx) {
346 oss << "pressure_" << FluidSystem::phaseName(/*phaseIdx=*/0);
347 }
348 else if (Indices::switch0Idx <= pvIdx &&
349 pvIdx < Indices::switch0Idx + numPhases - 1)
350 {
351 oss << "switch_" << pvIdx - Indices::switch0Idx;
352 }
353 else if (Indices::switch0Idx + numPhases - 1 <= pvIdx &&
354 pvIdx < Indices::switch0Idx + numComponents - 1)
355 {
356 oss << "auxMoleFrac^" << FluidSystem::componentName(pvIdx);
357 } else {
358 assert(false);
359 }
360
361 return oss.str();
362 }
363
367 std::string eqName(unsigned eqIdx) const
368 {
369 const std::string s = EnergyModule::eqName(eqIdx);
370 if (!s.empty()) {
371 return s;
372 }
373
374 std::ostringstream oss;
375 if (Indices::conti0EqIdx <= eqIdx &&
376 eqIdx < Indices::conti0EqIdx + numComponents)
377 {
378 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
379 oss << "continuity^" << FluidSystem::componentName(compIdx);
380 }
381 else {
382 assert(false);
383 }
384
385 return oss.str();
386 }
387
392 {
393 ParentType::updateFailed();
394 numSwitched_ = 0;
395 }
396
401 {
402 ParentType::updateBegin();
403
404 // find the a reference pressure. The first degree of freedom
405 // might correspond to non-interior entities which would lead
406 // to an undefined value, so we have to iterate...
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(/*timeIdx=*/0)[dofIdx][/*pvIdx=*/Indices::pressure0Idx];
412 if (referencePressure_ > 0.0) {
413 break;
414 }
415 }
416 }
417 }
418
422 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
423 {
424 const Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
425 if (tmp > 0) {
426 // energy related quantity
427 return tmp;
428 }
429
430 if (Indices::pressure0Idx == pvIdx) {
431 return 10 / referencePressure_;
432 }
433
434 if (Indices::switch0Idx <= pvIdx &&
435 pvIdx < Indices::switch0Idx + numPhases - 1)
436 {
437 const unsigned phaseIdx = pvIdx - Indices::switch0Idx;
438
439 if (!this->solution(/*timeIdx=*/0)[globalDofIdx].phaseIsPresent(phaseIdx)) {
440 // for saturations, the weight is always 1
441 return 1;
442 }
443
444 // for saturations, the PvsMoleSaturationsBaseWeight
445 // property determines the weight
446 return getPropValue<TypeTag, Properties::PvsSaturationsBaseWeight>();
447 }
448
449 // for mole fractions, the PvsMoleFractionsBaseWeight
450 // property determines the weight
451 return getPropValue<TypeTag, Properties::PvsMoleFractionsBaseWeight>();
452 }
453
457 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
458 {
459 const Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
460 if (tmp > 0) {
461 // energy related equation
462 return tmp;
463 }
464
465 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
466 assert(compIdx <= numComponents);
467
468 // make all kg equal
469 return FluidSystem::molarMass(compIdx);
470 }
471
476 {
477 ParentType::advanceTimeLevel();
478 numSwitched_ = 0;
479 }
480
485 bool switched() const
486 { return numSwitched_ > 0; }
487
491 template <class DofEntity>
492 void serializeEntity(std::ostream& outstream, const DofEntity& dofEntity)
493 {
494 // write primary variables
495 ParentType::serializeEntity(outstream, dofEntity);
496
497 const unsigned dofIdx = static_cast<unsigned>(this->dofMapper().index(dofEntity));
498 if (!outstream.good()) {
499 throw std::runtime_error("Could not serialize DOF " + std::to_string(dofIdx));
500 }
501
502 outstream << this->solution(/*timeIdx=*/0)[dofIdx].phasePresence() << " ";
503 }
504
508 template <class DofEntity>
509 void deserializeEntity(std::istream& instream, const DofEntity& dofEntity)
510 {
511 // read primary variables
512 ParentType::deserializeEntity(instream, dofEntity);
513
514 // read phase presence
515 const unsigned dofIdx = static_cast<unsigned>(this->dofMapper().index(dofEntity));
516 if (!instream.good()) {
517 throw std::runtime_error("Could not deserialize DOF " + std::to_string(dofIdx));
518 }
519
520 short tmp;
521 instream >> tmp;
522 this->solution(/*timeIdx=*/0)[dofIdx].setPhasePresence(tmp);
523 this->solution(/*timeIdx=*/1)[dofIdx].setPhasePresence(tmp);
524 }
525
534 {
535 numSwitched_ = 0;
536
537 int succeeded;
538 try {
539 std::vector<bool> visited(this->numGridDof(), false);
540 ElementContext elemCtx(this->simulator_);
541
542 for (const auto& elem : elements(this->gridView_, Dune::Partitions::interior)) {
543 elemCtx.updateStencil(elem);
544
545 const std::size_t numLocalDof = elemCtx.stencil(/*timeIdx=*/0).numPrimaryDof();
546 for (unsigned dofIdx = 0; dofIdx < numLocalDof; ++dofIdx) {
547 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
548
549 if (visited[globalIdx]) {
550 continue;
551 }
552 visited[globalIdx] = true;
553
554 // compute the intensive quantities of the current degree of freedom
555 auto& priVars = this->solution(/*timeIdx=*/0)[globalIdx];
556 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
557 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
558
559 // evaluate primary variable switch
560 const short oldPhasePresence = priVars.phasePresence();
561
562 // set the primary variables and the new phase state
563 // from the current fluid state
564 priVars.assignNaive(intQuants.fluidState());
565
566 if (oldPhasePresence != priVars.phasePresence()) {
567 if (verbosity_ > 1) {
568 printSwitchedPhases_(elemCtx,
569 dofIdx,
570 intQuants.fluidState(),
571 oldPhasePresence,
572 priVars);
573 }
574 ++numSwitched_;
575 }
576 }
577 }
578
579 succeeded = 1;
580 }
581 catch (...)
582 {
583 std::cout << "rank " << this->simulator_.gridView().comm().rank()
584 << " caught an exception during primary variable switching"
585 << "\n" << std::flush;
586 succeeded = 0;
587 }
588 succeeded = this->simulator_.gridView().comm().min(succeeded);
589
590 if (!succeeded) {
591 throw NumericalProblem("A process did not succeed in adapting the primary variables");
592 }
593
594 // make sure that if there was a variable switch in an
595 // other partition we will also set the switch flag
596 // for our partition.
597 numSwitched_ = this->gridView_.comm().sum(numSwitched_);
598
599 if (verbosity_ > 0) {
600 this->simulator_.model().newtonMethod().endIterMsg()
601 << ", num switched=" << numSwitched_;
602 }
603 }
604
605 template <class FluidState>
606 void printSwitchedPhases_(const ElementContext& elemCtx,
607 unsigned dofIdx,
608 const FluidState& fs,
609 short oldPhasePresence,
610 const PrimaryVariables& newPv) const
611 {
612 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
613
614 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
615 const bool oldPhasePresent = (oldPhasePresence & (1 << phaseIdx)) > 0;
616 const bool newPhasePresent = newPv.phaseIsPresent(phaseIdx);
617 if (oldPhasePresent == newPhasePresent) {
618 continue;
619 }
620
621 const auto& pos = elemCtx.pos(dofIdx, /*timeIdx=*/0);
622 if (oldPhasePresent) {
623 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
624 << "' phase disappears at position " << pos
625 << ". saturation=" << fs.saturation(phaseIdx)
626 << std::flush;
627 }
628 else {
629 Scalar sumx = 0;
630 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
631 sumx += FsToolbox::value(fs.moleFraction(phaseIdx, compIdx));
632 }
633
634 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
635 << "' phase appears at position " << pos
636 << " sum x = " << sumx << std::flush;
637 }
638 }
639
640 std::cout << ", new primary variables: ";
641 newPv.print(std::cout);
642 std::cout << "\n" << std::flush;
643 }
644
646 {
648
649 // add the VTK output modules which are meaningful for the model
650 this->addOutputModule(std::make_unique<VtkPhasePresenceModule<TypeTag>>(this->simulator_));
651 this->addOutputModule(std::make_unique<VtkCompositionModule<TypeTag>>(this->simulator_));
652 if constexpr (enableDiffusion) {
653 this->addOutputModule(std::make_unique<VtkDiffusionModule<TypeTag>>(this->simulator_));
654 }
655 if constexpr (enableEnergy) {
656 this->addOutputModule(std::make_unique<VtkEnergyModule<TypeTag>>(this->simulator_));
657 }
658 }
659
661
662 // number of switches of the phase state in the last Newton
663 // iteration
664 unsigned numSwitched_;
665
666 // verbosity of the model
668};
669
670} // namespace Opm
671
672#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:168
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:194
void registerOutputModules_()
Definition: multiphasebasemodel.hh:270
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:475
Scalar referencePressure_
Definition: pvsmodel.hh:660
void serializeEntity(std::ostream &outstream, const DofEntity &dofEntity)
Write the current solution for a degree of freedom to a restart file.
Definition: pvsmodel.hh:492
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 switchPrimaryVars_()
Definition: pvsmodel.hh:533
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: pvsmodel.hh:400
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: pvsmodel.hh:457
int verbosity_
Definition: pvsmodel.hh:667
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:485
static void registerParameters()
Register all run-time parameters for the PVS compositional model.
Definition: pvsmodel.hh:307
void printSwitchedPhases_(const ElementContext &elemCtx, unsigned dofIdx, const FluidState &fs, short oldPhasePresence, const PrimaryVariables &newPv) const
Definition: pvsmodel.hh:606
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: pvsmodel.hh:422
void registerOutputModules_()
Definition: pvsmodel.hh:645
unsigned numSwitched_
Definition: pvsmodel.hh:664
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:509
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: pvsmodel.hh:391
PvsModel(Simulator &simulator)
Definition: pvsmodel.hh:297
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:79
Definition: blackoilboundaryratevector.hh:39
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: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
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