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/material/densead/Math.hpp>
32
33#include "pvsproperties.hh"
34#include "pvslocalresidual.hh"
35#include "pvsnewtonmethod.hh"
37#include "pvsratevector.hh"
41#include "pvsindices.hh"
42
43#include <opm/common/Exceptions.hpp>
44
51
52#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
53#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
54
55#include <iostream>
56#include <sstream>
57#include <string>
58#include <vector>
59
60namespace Opm {
61template <class TypeTag>
62class PvsModel;
63}
64
65namespace Opm::Properties {
66
67namespace TTag {
69struct PvsModel { using InheritsFrom = std::tuple<VtkDiffusion,
74} // namespace TTag
75
77template<class TypeTag>
78struct LocalResidual<TypeTag, TTag::PvsModel> { using type = Opm::PvsLocalResidual<TypeTag>; };
79
81template<class TypeTag>
82struct NewtonMethod<TypeTag, TTag::PvsModel> { using type = Opm::PvsNewtonMethod<TypeTag>; };
83
85template<class TypeTag>
86struct Model<TypeTag, TTag::PvsModel> { using type = Opm::PvsModel<TypeTag>; };
87
89template<class TypeTag>
91
93template<class TypeTag>
94struct RateVector<TypeTag, TTag::PvsModel> { using type = Opm::PvsRateVector<TypeTag>; };
95
97template<class TypeTag>
99
101template<class TypeTag>
103
105template<class TypeTag>
107
109template<class TypeTag>
110struct Indices<TypeTag, TTag::PvsModel> { using type = Opm::PvsIndices<TypeTag, /*PVIdx=*/0>; };
111
112// set the model to a medium verbosity
113template<class TypeTag>
114struct PvsVerbosity<TypeTag, TTag::PvsModel> { static constexpr int value = 1; };
115
117template<class TypeTag>
118struct EnableEnergy<TypeTag, TTag::PvsModel> { static constexpr bool value = false; };
119
120// disable molecular diffusion by default
121template<class TypeTag>
122struct EnableDiffusion<TypeTag, TTag::PvsModel> { static constexpr bool value = false; };
123
125template<class TypeTag>
126struct PvsPressureBaseWeight<TypeTag, TTag::PvsModel>
127{
129 static constexpr type value = 1.0;
130};
131
133template<class TypeTag>
134struct PvsSaturationsBaseWeight<TypeTag, TTag::PvsModel>
135{
137 static constexpr type value = 1.0;
138};
139
141template<class TypeTag>
143{
145 static constexpr type value = 1.0;
146};
147
148} // namespace Opm::Properties
149
150namespace Opm {
151
248template <class TypeTag>
250 : public MultiPhaseBaseModel<TypeTag>
251{
252 using ParentType = MultiPhaseBaseModel<TypeTag>;
253
258
263
264 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
265 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
266 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
267 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
268
269 using Element = typename GridView::template Codim<0>::Entity;
270 using ElementIterator = typename GridView::template Codim<0>::Iterator;
271
273
274public:
275 PvsModel(Simulator& simulator)
276 : ParentType(simulator)
277 {
278 verbosity_ = Parameters::get<TypeTag, Properties::PvsVerbosity>();
279 numSwitched_ = 0;
280 }
281
285 static void registerParameters()
286 {
288
289 // register runtime parameters of the VTK output modules
292
293 if (enableDiffusion)
295
296 if (enableEnergy)
298
299 Parameters::registerParam<TypeTag, Properties::PvsVerbosity>
300 ("The verbosity level of the primary variable "
301 "switching model");
302 }
303
307 static std::string name()
308 { return "pvs"; }
309
313 std::string primaryVarName(unsigned pvIdx) const
314 {
315 std::string s;
316 if (!(s = EnergyModule::primaryVarName(pvIdx)).empty())
317 return s;
318
319 std::ostringstream oss;
320 if (pvIdx == Indices::pressure0Idx)
321 oss << "pressure_" << FluidSystem::phaseName(/*phaseIdx=*/0);
322 else if (Indices::switch0Idx <= pvIdx
323 && pvIdx < Indices::switch0Idx + numPhases - 1)
324 oss << "switch_" << pvIdx - Indices::switch0Idx;
325 else if (Indices::switch0Idx + numPhases - 1 <= pvIdx
326 && pvIdx < Indices::switch0Idx + numComponents - 1)
327 oss << "auxMoleFrac^" << FluidSystem::componentName(pvIdx);
328 else
329 assert(false);
330
331 return oss.str();
332 }
333
337 std::string eqName(unsigned eqIdx) const
338 {
339 std::string s;
340 if (!(s = EnergyModule::eqName(eqIdx)).empty())
341 return s;
342
343 std::ostringstream oss;
344 if (Indices::conti0EqIdx <= eqIdx && eqIdx < Indices::conti0EqIdx
345 + numComponents) {
346 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
347 oss << "continuity^" << FluidSystem::componentName(compIdx);
348 }
349 else
350 assert(false);
351
352 return oss.str();
353 }
354
359 {
360 ParentType::updateFailed();
361 numSwitched_ = 0;
362 }
363
368 {
369 ParentType::updateBegin();
370
371 // find the a reference pressure. The first degree of freedom
372 // might correspond to non-interior entities which would lead
373 // to an undefined value, so we have to iterate...
374 size_t nDof = this->numTotalDof();
375 for (unsigned dofIdx = 0; dofIdx < nDof; ++ dofIdx) {
376 if (this->dofTotalVolume(dofIdx) > 0.0) {
378 this->solution(/*timeIdx=*/0)[dofIdx][/*pvIdx=*/Indices::pressure0Idx];
379 if (referencePressure_ > 0.0)
380 break;
381 }
382 }
383 }
384
388 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
389 {
390 Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
391 if (tmp > 0)
392 // energy related quantity
393 return tmp;
394
395 if (Indices::pressure0Idx == pvIdx) {
396 return 10 / referencePressure_;
397 }
398
399 if (Indices::switch0Idx <= pvIdx && pvIdx < Indices::switch0Idx
400 + numPhases - 1) {
401 unsigned phaseIdx = pvIdx - Indices::switch0Idx;
402
403 if (!this->solution(/*timeIdx=*/0)[globalDofIdx].phaseIsPresent(phaseIdx))
404 // for saturations, the weight is always 1
405 return 1;
406
407 // for saturations, the PvsMoleSaturationsBaseWeight
408 // property determines the weight
409 return getPropValue<TypeTag, Properties::PvsSaturationsBaseWeight>();
410 }
411
412 // for mole fractions, the PvsMoleFractionsBaseWeight
413 // property determines the weight
414 return getPropValue<TypeTag, Properties::PvsMoleFractionsBaseWeight>();
415 }
416
420 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
421 {
422 Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
423 if (tmp > 0)
424 // energy related equation
425 return tmp;
426
427 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
428 assert(compIdx <= numComponents);
429
430 // make all kg equal
431 return FluidSystem::molarMass(compIdx);
432 }
433
438 {
439 ParentType::advanceTimeLevel();
440 numSwitched_ = 0;
441 }
442
447 bool switched() const
448 { return numSwitched_ > 0; }
449
453 template <class DofEntity>
454 void serializeEntity(std::ostream& outstream, const DofEntity& dofEntity)
455 {
456 // write primary variables
457 ParentType::serializeEntity(outstream, dofEntity);
458
459 unsigned dofIdx = static_cast<unsigned>(this->dofMapper().index(dofEntity));
460 if (!outstream.good())
461 throw std::runtime_error("Could not serialize DOF "+std::to_string(dofIdx));
462
463 outstream << this->solution(/*timeIdx=*/0)[dofIdx].phasePresence() << " ";
464 }
465
469 template <class DofEntity>
470 void deserializeEntity(std::istream& instream, const DofEntity& dofEntity)
471 {
472 // read primary variables
473 ParentType::deserializeEntity(instream, dofEntity);
474
475 // read phase presence
476 unsigned dofIdx = static_cast<unsigned>(this->dofMapper().index(dofEntity));
477 if (!instream.good())
478 throw std::runtime_error("Could not deserialize DOF "+std::to_string(dofIdx));
479
480 short tmp;
481 instream >> tmp;
482 this->solution(/*timeIdx=*/0)[dofIdx].setPhasePresence(tmp);
483 this->solution(/*timeIdx=*/1)[dofIdx].setPhasePresence(tmp);
484 }
485
494 {
495 numSwitched_ = 0;
496
497 int succeeded;
498 try {
499 std::vector<bool> visited(this->numGridDof(), false);
500 ElementContext elemCtx(this->simulator_);
501
502 for (const auto& elem : elements(this->gridView_, Dune::Partitions::interior)) {
503 elemCtx.updateStencil(elem);
504
505 size_t numLocalDof = elemCtx.stencil(/*timeIdx=*/0).numPrimaryDof();
506 for (unsigned dofIdx = 0; dofIdx < numLocalDof; ++dofIdx) {
507 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
508
509 if (visited[globalIdx])
510 continue;
511 visited[globalIdx] = true;
512
513 // compute the intensive quantities of the current degree of freedom
514 auto& priVars = this->solution(/*timeIdx=*/0)[globalIdx];
515 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
516 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
517
518 // evaluate primary variable switch
519 short oldPhasePresence = priVars.phasePresence();
520
521 // set the primary variables and the new phase state
522 // from the current fluid state
523 priVars.assignNaive(intQuants.fluidState());
524
525 if (oldPhasePresence != priVars.phasePresence()) {
526 if (verbosity_ > 1)
527 printSwitchedPhases_(elemCtx,
528 dofIdx,
529 intQuants.fluidState(),
530 oldPhasePresence,
531 priVars);
532 ++numSwitched_;
533 }
534 }
535 }
536
537 succeeded = 1;
538 }
539 catch (...)
540 {
541 std::cout << "rank " << this->simulator_.gridView().comm().rank()
542 << " caught an exception during primary variable switching"
543 << "\n" << std::flush;
544 succeeded = 0;
545 }
546 succeeded = this->simulator_.gridView().comm().min(succeeded);
547
548 if (!succeeded)
549 throw NumericalProblem("A process did not succeed in adapting the primary variables");
550
551 // make sure that if there was a variable switch in an
552 // other partition we will also set the switch flag
553 // for our partition.
554 numSwitched_ = this->gridView_.comm().sum(numSwitched_);
555
556 if (verbosity_ > 0)
557 this->simulator_.model().newtonMethod().endIterMsg()
558 << ", num switched=" << numSwitched_;
559 }
560
561 template <class FluidState>
562 void printSwitchedPhases_(const ElementContext& elemCtx,
563 unsigned dofIdx,
564 const FluidState& fs,
565 short oldPhasePresence,
566 const PrimaryVariables& newPv) const
567 {
568 using FsToolbox = Opm::MathToolbox<typename FluidState::Scalar>;
569
570 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
571 bool oldPhasePresent = (oldPhasePresence& (1 << phaseIdx)) > 0;
572 bool newPhasePresent = newPv.phaseIsPresent(phaseIdx);
573 if (oldPhasePresent == newPhasePresent)
574 continue;
575
576 const auto& pos = elemCtx.pos(dofIdx, /*timeIdx=*/0);
577 if (oldPhasePresent && !newPhasePresent) {
578 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
579 << "' phase disappears at position " << pos
580 << ". saturation=" << fs.saturation(phaseIdx)
581 << std::flush;
582 }
583 else {
584 Scalar sumx = 0;
585 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx)
586 sumx += FsToolbox::value(fs.moleFraction(phaseIdx, compIdx));
587
588 std::cout << "'" << FluidSystem::phaseName(phaseIdx)
589 << "' phase appears at position " << pos
590 << " sum x = " << sumx << std::flush;
591 }
592 }
593
594 std::cout << ", new primary variables: ";
595 newPv.print();
596 std::cout << "\n" << std::flush;
597 }
598
600 {
602
603 // add the VTK output modules which are meaningful for the model
604 this->addOutputModule(new Opm::VtkPhasePresenceModule<TypeTag>(this->simulator_));
605 this->addOutputModule(new Opm::VtkCompositionModule<TypeTag>(this->simulator_));
606 if (enableDiffusion)
607 this->addOutputModule(new Opm::VtkDiffusionModule<TypeTag>(this->simulator_));
608 if (enableEnergy)
609 this->addOutputModule(new Opm::VtkEnergyModule<TypeTag>(this->simulator_));
610 }
611
612 mutable Scalar referencePressure_;
613
614 // number of switches of the phase state in the last Newton
615 // iteration
616 unsigned numSwitched_;
617
618 // verbosity of the model
620};
621}
622
623#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:48
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:153
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: multiphasebasemodel.hh:179
void registerOutputModules_()
Definition: multiphasebasemodel.hh:254
Implements a rate vector on the boundary for the fully implicit compositional multi-phase primary var...
Definition: pvsboundaryratevector.hh:47
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:61
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Definition: pvslocalresidual.hh:48
A generic compositional multi-phase model using primary-variable switching.
Definition: pvsmodel.hh:251
void advanceTimeLevel()
Called by the problem if a time integration was successful, post processing of the solution is done a...
Definition: pvsmodel.hh:437
Scalar referencePressure_
Definition: pvsmodel.hh:612
void serializeEntity(std::ostream &outstream, const DofEntity &dofEntity)
Write the current solution for a degree of freedom to a restart file.
Definition: pvsmodel.hh:454
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: pvsmodel.hh:337
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: pvsmodel.hh:313
void switchPrimaryVars_()
Definition: pvsmodel.hh:493
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: pvsmodel.hh:367
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: pvsmodel.hh:420
int verbosity_
Definition: pvsmodel.hh:619
static std::string name()
Definition: pvsmodel.hh:307
bool switched() const
Return true if the primary variables were switched for at least one vertex after the last timestep.
Definition: pvsmodel.hh:447
static void registerParameters()
Register all run-time parameters for the PVS compositional model.
Definition: pvsmodel.hh:285
void printSwitchedPhases_(const ElementContext &elemCtx, unsigned dofIdx, const FluidState &fs, short oldPhasePresence, const PrimaryVariables &newPv) const
Definition: pvsmodel.hh:562
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: pvsmodel.hh:388
void registerOutputModules_()
Definition: pvsmodel.hh:599
unsigned numSwitched_
Definition: pvsmodel.hh:616
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:470
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: pvsmodel.hh:358
PvsModel(Simulator &simulator)
Definition: pvsmodel.hh:275
A newton solver which is specific to the compositional multi-phase PVS model.
Definition: pvsnewtonmethod.hh:52
Represents the primary variables used in the primary variable switching compositional model.
Definition: pvsprimaryvariables.hh:60
Implements a vector representing molar, mass or volumetric rates.
Definition: pvsratevector.hh:53
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hh:97
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:124
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Definition: vtkdiffusionmodule.hh:82
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:109
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hh:85
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:112
VTK output module for the fluid composition.
Definition: vtkphasepresencemodule.hh:62
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkphasepresencemodule.hh:84
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
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:242
Declares the properties required for the compositional multi-phase primary variable switching model.
Type of object for specifying boundary conditions.
Definition: fvbaseproperties.hh:136
Enable diffusive fluxes?
Definition: multiphasebaseproperties.hh:82
Specify whether energy should be considered as a conservation quantity or not.
Definition: multiphasebaseproperties.hh:76
Data required to calculate a flux over a face.
Definition: fvbaseproperties.hh:166
Enumerations used by the model.
Definition: multiphasebaseproperties.hh:48
The secondary variables within a sub-control volume.
Definition: fvbaseproperties.hh:150
The type of the local residual function.
Definition: fvbaseproperties.hh:111
The type of the model.
Definition: basicproperties.hh:81
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:147
GetPropType< TypeTag, Scalar > type
Definition: pvsmodel.hh:144
The basis value for the weight of the mole fraction primary variables.
Definition: pvsproperties.hh:54
GetPropType< TypeTag, Scalar > type
Definition: pvsmodel.hh:128
The basis value for the weight of the pressure primary variable.
Definition: pvsproperties.hh:48
GetPropType< TypeTag, Scalar > type
Definition: pvsmodel.hh:136
The basis value for the weight of the saturation primary variables.
Definition: pvsproperties.hh:51
The verbosity of the model (0 -> do not print anything, 2 -> spam stdout a lot)
Definition: pvsproperties.hh:45
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:133
Definition: multiphasebasemodel.hh:57
The type tag for the isothermal single phase problems.
Definition: pvsmodel.hh:69
std::tuple< VtkDiffusion, VtkEnergy, VtkComposition, VtkPhasePresence, MultiPhaseBaseModel > InheritsFrom
Definition: pvsmodel.hh:73
Definition: vtkcompositionmodule.hh:43
Definition: vtkdiffusionmodule.hh:46
Definition: vtkenergymodule.hh:43
Definition: vtkphasepresencemodule.hh:41