opm-simulators
Loading...
Searching...
No Matches
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>
86{ using type = Opm::PvsLocalResidual<TypeTag>; };
87
89template<class TypeTag>
90struct NewtonMethod<TypeTag, TTag::PvsModel>
91{ using type = Opm::PvsNewtonMethod<TypeTag>; };
92
94template<class TypeTag>
95struct Model<TypeTag, TTag::PvsModel>
96{ using type = Opm::PvsModel<TypeTag>; };
97
99template<class TypeTag>
101{ using type = Opm::PvsPrimaryVariables<TypeTag>; };
102
104template<class TypeTag>
105struct RateVector<TypeTag, TTag::PvsModel>
106{ using type = Opm::PvsRateVector<TypeTag>; };
107
109template<class TypeTag>
112
114template<class TypeTag>
117
119template<class TypeTag>
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>
136{ static constexpr bool value = false; };
137
139template<class TypeTag>
141{
142 using type = GetPropType<TypeTag, Scalar>;
143 static constexpr type value = 1.0;
144};
145
147template<class TypeTag>
149{
150 using type = GetPropType<TypeTag, Scalar>;
151 static constexpr type value = 1.0;
152};
153
155template<class TypeTag>
157{
158 using type = GetPropType<TypeTag, Scalar>;
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>
271class PvsModel
272 : public MultiPhaseBaseModel<TypeTag>
273{
274 using ParentType = MultiPhaseBaseModel<TypeTag>;
275
280
285
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 {
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
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) {
410 referencePressure_ =
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
447 }
448
449 // for mole fractions, the PvsMoleFractionsBaseWeight
450 // property determines the weight
452 }
453
460 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
461 {
462 const Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
463 if (tmp > 0) {
464 // energy related equation
465 return tmp;
466 }
467
468 const unsigned compIdx = eqIdx - Indices::conti0EqIdx;
469 assert(compIdx <= numComponents);
470
471 // make all kg equal
472 return FluidSystem::molarMass(compIdx);
473 }
474
479 {
480 ParentType::advanceTimeLevel();
481 numSwitched_ = 0;
482 }
483
488 bool switched() const
489 { return numSwitched_ > 0; }
490
497 template <class DofEntity>
498 void serializeEntity(std::ostream& outstream, const DofEntity& dofEntity)
499 {
500 // write primary variables
501 ParentType::serializeEntity(outstream, dofEntity);
502
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));
506 }
507
508 outstream << this->solution(/*timeIdx=*/0)[dofIdx].phasePresence() << " ";
509 }
510
517 template <class DofEntity>
518 void deserializeEntity(std::istream& instream, const DofEntity& dofEntity)
519 {
520 // read primary variables
521 ParentType::deserializeEntity(instream, dofEntity);
522
523 // read phase presence
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));
527 }
528
529 short tmp;
530 instream >> tmp;
531 this->solution(/*timeIdx=*/0)[dofIdx].setPhasePresence(tmp);
532 this->solution(/*timeIdx=*/1)[dofIdx].setPhasePresence(tmp);
533 }
534
542 void switchPrimaryVars_()
543 {
544 numSwitched_ = 0;
545
546 int succeeded;
547 try {
548 std::vector<bool> visited(this->numGridDof(), false);
549 ElementContext elemCtx(this->simulator_);
550
551 for (const auto& elem : elements(this->gridView_, Dune::Partitions::interior)) {
552 elemCtx.updateStencil(elem);
553
554 const std::size_t numLocalDof = elemCtx.stencil(/*timeIdx=*/0).numPrimaryDof();
555 for (unsigned dofIdx = 0; dofIdx < numLocalDof; ++dofIdx) {
556 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
557
558 if (visited[globalIdx]) {
559 continue;
560 }
561 visited[globalIdx] = true;
562
563 // compute the intensive quantities of the current degree of freedom
564 auto& priVars = this->solution(/*timeIdx=*/0)[globalIdx];
565 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
566 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
567
568 // evaluate primary variable switch
569 const short oldPhasePresence = priVars.phasePresence();
570
571 // set the primary variables and the new phase state
572 // from the current fluid state
573 priVars.assignNaive(intQuants.fluidState());
574
575 if (oldPhasePresence != priVars.phasePresence()) {
576 if (verbosity_ > 1) {
577 printSwitchedPhases_(elemCtx,
578 dofIdx,
579 intQuants.fluidState(),
580 oldPhasePresence,
581 priVars);
582 }
583 ++numSwitched_;
584 }
585 }
586 }
587
588 succeeded = 1;
589 }
590 catch (...)
591 {
592 std::cout << "rank " << this->simulator_.gridView().comm().rank()
593 << " caught an exception during primary variable switching"
594 << "\n" << std::flush;
595 succeeded = 0;
596 }
597 succeeded = this->simulator_.gridView().comm().min(succeeded);
598
599 if (!succeeded) {
600 throw NumericalProblem("A process did not succeed in adapting the primary variables");
601 }
602
603 // make sure that if there was a variable switch in an
604 // other partition we will also set the switch flag
605 // for our partition.
606 numSwitched_ = this->gridView_.comm().sum(numSwitched_);
607
608 if (verbosity_ > 0) {
609 this->simulator_.model().newtonMethod().endIterMsg()
610 << ", num switched=" << numSwitched_;
611 }
612 }
613
614 template <class FluidState>
615 void printSwitchedPhases_(const ElementContext& elemCtx,
616 unsigned dofIdx,
617 const FluidState& fs,
618 short oldPhasePresence,
619 const PrimaryVariables& newPv) const
620 {
621 using FsToolbox = MathToolbox<typename FluidState::ValueType>;
622
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) {
627 continue;
628 }
629
630 const auto& pos = elemCtx.pos(dofIdx, /*timeIdx=*/0);
631 if (oldPhasePresent) {
632 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
633 << "' phase disappears at position " << pos
634 << ". saturation=" << fs.saturation(phaseIdx)
635 << std::flush;
636 }
637 else {
638 Scalar sumx = 0;
639 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
640 sumx += FsToolbox::value(fs.moleFraction(phaseIdx, compIdx));
641 }
642
643 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
644 << "' phase appears at position " << pos
645 << " sum x = " << sumx << std::flush;
646 }
647 }
648
649 std::cout << ", new primary variables: ";
650 newPv.print(std::cout);
651 std::cout << "\n" << std::flush;
652 }
653
654 void registerOutputModules_()
655 {
656 ParentType::registerOutputModules_();
657
658 // add the VTK output modules which are meaningful for the model
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_));
663 }
664 if constexpr (enableEnergy) {
665 this->addOutputModule(std::make_unique<VtkEnergyModule<TypeTag>>(this->simulator_));
666 }
667 }
668
669 Scalar referencePressure_{};
670
671 // number of switches of the phase state in the last Newton
672 // iteration
673 unsigned numSwitched_;
674
675 // verbosity of the model
676 int verbosity_;
677};
678
679} // namespace Opm
680
681#endif
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:197
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:82
This file contains a set of helper functions used by VFPProd / VFPInj.
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
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: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
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.