ncpmodel.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_NCP_MODEL_HH
29#define EWOMS_NCP_MODEL_HH
30
31#include <opm/material/densead/Math.hpp>
32
33#include "ncpproperties.hh"
34#include "ncplocalresidual.hh"
38#include "ncpratevector.hh"
40#include "ncpnewtonmethod.hh"
41#include "ncpindices.hh"
42
43#include <opm/common/Exceptions.hpp>
44
51
52#include <opm/material/common/Valgrind.hpp>
53
54#include <dune/common/fvector.hh>
55
56#include <sstream>
57#include <string>
58#include <vector>
59#include <array>
60
61namespace Opm {
62template <class TypeTag>
63class NcpModel;
64}
65
66namespace Opm::Properties {
67
68namespace TTag {
72struct NcpModel { using InheritsFrom = std::tuple<MultiPhaseBaseModel>; };
73} // namespace TTag
74
76template<class TypeTag>
77struct LocalResidual<TypeTag, TTag::NcpModel> { using type = NcpLocalResidual<TypeTag>; };
78
80template<class TypeTag>
81struct NewtonMethod<TypeTag, TTag::NcpModel> { using type = NcpNewtonMethod<TypeTag>; };
82
84template<class TypeTag>
85struct Model<TypeTag, TTag::NcpModel> { using type = NcpModel<TypeTag>; };
86
88template<class TypeTag>
89struct BaseProblem<TypeTag, TTag::NcpModel> { using type = MultiPhaseBaseProblem<TypeTag>; };
90
92template<class TypeTag>
93struct EnableEnergy<TypeTag, TTag::NcpModel> { static constexpr bool value = false; };
94
96template<class TypeTag>
97struct EnableDiffusion<TypeTag, TTag::NcpModel> { static constexpr bool value = false; };
98
100template<class TypeTag>
101struct RateVector<TypeTag, TTag::NcpModel> { using type = NcpRateVector<TypeTag>; };
102
104template<class TypeTag>
106
108template<class TypeTag>
109struct PrimaryVariables<TypeTag, TTag::NcpModel> { using type = NcpPrimaryVariables<TypeTag>; };
110
112template<class TypeTag>
114
116template<class TypeTag>
118
120template<class TypeTag>
121struct Indices<TypeTag, TTag::NcpModel> { using type = NcpIndices<TypeTag, 0>; };
122
124template<class TypeTag>
125struct NcpPressureBaseWeight<TypeTag, TTag::NcpModel>
126{
128 static constexpr type value = 1.0;
129};
131template<class TypeTag>
132struct NcpSaturationsBaseWeight<TypeTag, TTag::NcpModel>
133{
135 static constexpr type value = 1.0;
136};
138template<class TypeTag>
139struct NcpFugacitiesBaseWeight<TypeTag, TTag::NcpModel>
140{
142 static constexpr type value = 1.0e-6;
143};
144
145} // namespace Opm::Properties
146
147namespace Opm {
148
223template <class TypeTag>
225 : public MultiPhaseBaseModel<TypeTag>
226{
227 using ParentType = MultiPhaseBaseModel<TypeTag>;
228
235
236 enum { numPhases = FluidSystem::numPhases };
237 enum { numComponents = FluidSystem::numComponents };
238 enum { fugacity0Idx = Indices::fugacity0Idx };
239 enum { pressure0Idx = Indices::pressure0Idx };
240 enum { saturation0Idx = Indices::saturation0Idx };
241 enum { conti0EqIdx = Indices::conti0EqIdx };
242 enum { ncp0EqIdx = Indices::ncp0EqIdx };
243 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
244 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
245
246 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
247
248 using Toolbox = MathToolbox<Evaluation>;
249
252
253public:
254 NcpModel(Simulator& simulator)
255 : ParentType(simulator)
256 {}
257
261 static void registerParameters()
262 {
264
265 DiffusionModule::registerParameters();
266 EnergyModule::registerParameters();
267
268 // register runtime parameters of the VTK output modules
270
271 if (enableDiffusion)
273
274 if (enableEnergy)
276 }
277
282 {
283 ParentType::finishInit();
284
285 minActivityCoeff_.resize(this->numGridDof());
286 std::fill(minActivityCoeff_.begin(), minActivityCoeff_.end(), 1.0);
287 }
288
290 {
291 ParentType::adaptGrid();
292 minActivityCoeff_.resize(this->numGridDof());
293 }
294
298 static std::string name()
299 { return "ncp"; }
300
304 std::string primaryVarName(unsigned pvIdx) const
305 {
306 std::string s;
307 if (!(s = EnergyModule::primaryVarName(pvIdx)).empty())
308 return s;
309
310 std::ostringstream oss;
311 if (pvIdx == pressure0Idx)
312 oss << "pressure_" << FluidSystem::phaseName(/*phaseIdx=*/0);
313 else if (saturation0Idx <= pvIdx && pvIdx < saturation0Idx + (numPhases - 1))
314 oss << "saturation_" << FluidSystem::phaseName(/*phaseIdx=*/pvIdx - saturation0Idx);
315 else if (fugacity0Idx <= pvIdx && pvIdx < fugacity0Idx + numComponents)
316 oss << "fugacity^" << FluidSystem::componentName(pvIdx - fugacity0Idx);
317 else
318 assert(false);
319
320 return oss.str();
321 }
322
326 std::string eqName(unsigned eqIdx) const
327 {
328 std::string s;
329 if (!(s = EnergyModule::eqName(eqIdx)).empty())
330 return s;
331
332 std::ostringstream oss;
333 if (conti0EqIdx <= eqIdx && eqIdx < conti0EqIdx + numComponents)
334 oss << "continuity^" << FluidSystem::componentName(eqIdx - conti0EqIdx);
335 else if (ncp0EqIdx <= eqIdx && eqIdx < ncp0EqIdx + numPhases)
336 oss << "ncp_" << FluidSystem::phaseName(/*phaseIdx=*/eqIdx - ncp0EqIdx);
337 else
338 assert(false);
339
340 return oss.str();
341 }
342
347 {
348 ParentType::updateBegin();
349
350 // find the a reference pressure. The first degree of freedom
351 // might correspond to non-interior entities which would lead
352 // to an undefined value, so we have to iterate...
353 for (unsigned dofIdx = 0; dofIdx < this->numGridDof(); ++ dofIdx) {
354 if (this->isLocalDof(dofIdx)) {
356 this->solution(/*timeIdx=*/0)[dofIdx][/*pvIdx=*/Indices::pressure0Idx];
357 break;
358 }
359 }
360 }
361
365 void updatePVWeights(const ElementContext& elemCtx) const
366 {
367 for (unsigned dofIdx = 0; dofIdx < elemCtx.numDof(/*timeIdx=*/0); ++dofIdx) {
368 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
369
370 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
371 minActivityCoeff_[globalIdx][compIdx] = 1e100;
372 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
373 const auto& fs = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0).fluidState();
374
375 minActivityCoeff_[globalIdx][compIdx] =
376 std::min(minActivityCoeff_[globalIdx][compIdx],
377 Toolbox::value(fs.fugacityCoefficient(phaseIdx, compIdx))
378 * Toolbox::value(fs.pressure(phaseIdx)));
379 Valgrind::CheckDefined(minActivityCoeff_[globalIdx][compIdx]);
380 }
381 if (minActivityCoeff_[globalIdx][compIdx] <= 0)
382 throw NumericalProblem("The minimum activity coefficient for component "+std::to_string(compIdx)
383 +" on DOF "+std::to_string(globalIdx)+" is negative or zero!");
384 }
385 }
386 }
387
391 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
392 {
393 Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
394 Scalar result;
395 if (tmp > 0)
396 // energy related quantity
397 result = tmp;
398 else if (fugacity0Idx <= pvIdx && pvIdx < fugacity0Idx + numComponents) {
399 // component fugacity
400 unsigned compIdx = pvIdx - fugacity0Idx;
401 assert(compIdx <= numComponents);
402
403 Valgrind::CheckDefined(minActivityCoeff_[globalDofIdx][compIdx]);
404 static const Scalar fugacityBaseWeight =
405 getPropValue<TypeTag, Properties::NcpFugacitiesBaseWeight>();
406 result = fugacityBaseWeight / minActivityCoeff_[globalDofIdx][compIdx];
407 }
408 else if (Indices::pressure0Idx == pvIdx) {
409 static const Scalar pressureBaseWeight = getPropValue<TypeTag, Properties::NcpPressureBaseWeight>();
410 result = pressureBaseWeight / referencePressure_;
411 }
412 else {
413#ifndef NDEBUG
414 unsigned phaseIdx = pvIdx - saturation0Idx;
415 assert(phaseIdx < numPhases - 1);
416#endif
417
418 // saturation
419 static const Scalar saturationsBaseWeight =
420 getPropValue<TypeTag, Properties::NcpSaturationsBaseWeight>();
421 result = saturationsBaseWeight;
422 }
423
424 assert(std::isfinite(result));
425 assert(result > 0);
426
427 return result;
428 }
429
433 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
434 {
435 Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
436 if (tmp > 0)
437 // an energy related equation
438 return tmp;
439 // an NCP
440 else if (ncp0EqIdx <= eqIdx && eqIdx < Indices::ncp0EqIdx + numPhases)
441 return 1.0;
442
443 // a mass conservation equation
444 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
445 assert(compIdx <= numComponents);
446
447 // make all kg equal
448 return FluidSystem::molarMass(compIdx);
449 }
450
458 Scalar minActivityCoeff(unsigned globalDofIdx, unsigned compIdx) const
459 { return minActivityCoeff_[globalDofIdx][compIdx]; }
460
465 {
467
468 this->addOutputModule(new VtkCompositionModule<TypeTag>(this->simulator_));
469 if (enableDiffusion)
470 this->addOutputModule(new VtkDiffusionModule<TypeTag>(this->simulator_));
471 if (enableEnergy)
472 this->addOutputModule(new VtkEnergyModule<TypeTag>(this->simulator_));
473 }
474
475 mutable Scalar referencePressure_;
476 mutable std::vector<ComponentVector> minActivityCoeff_;
477};
478
479} // namespace Opm
480
481#endif
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:48
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
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
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:60
Implements a boundary vector for the fully implicit compositional multi-phase NCP model.
Definition: ncpboundaryratevector.hh:46
This template class represents the extensive quantities of the compositional NCP model.
Definition: ncpextensivequantities.hh:51
Contains the quantities which are are constant within a finite volume in the compositional multi-phas...
Definition: ncpintensivequantities.hh:58
Details needed to calculate the local residual in the compositional multi-phase NCP-model .
Definition: ncplocalresidual.hh:47
A compositional multi-phase model based on non-linear complementarity functions.
Definition: ncpmodel.hh:226
void updatePVWeights(const ElementContext &elemCtx) const
Update the weights of all primary variables within an element given the complete set of intensive qua...
Definition: ncpmodel.hh:365
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: ncpmodel.hh:326
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: ncpmodel.hh:304
NcpModel(Simulator &simulator)
Definition: ncpmodel.hh:254
static std::string name()
Definition: ncpmodel.hh:298
Scalar referencePressure_
Definition: ncpmodel.hh:475
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: ncpmodel.hh:391
void adaptGrid()
Definition: ncpmodel.hh:289
Scalar minActivityCoeff(unsigned globalDofIdx, unsigned compIdx) const
Returns the smallest activity coefficient of a component for the most current solution at a vertex.
Definition: ncpmodel.hh:458
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: ncpmodel.hh:261
void finishInit()
Apply the initial conditions to the model.
Definition: ncpmodel.hh:281
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: ncpmodel.hh:346
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: ncpmodel.hh:433
void registerOutputModules_()
Definition: ncpmodel.hh:464
std::vector< ComponentVector > minActivityCoeff_
Definition: ncpmodel.hh:476
A Newton solver specific to the NCP model.
Definition: ncpnewtonmethod.hh:55
Represents the primary variables used by the compositional multi-phase NCP model.
Definition: ncpprimaryvariables.hh:55
Implements a vector representing mass, molar or volumetric rates.
Definition: ncpratevector.hh:51
VTK output module for the fluid composition.
Definition: vtkcompositionmodule.hh:69
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkcompositionmodule.hh:96
VTK output module for quantities which make sense for models which incorperate molecular diffusion.
Definition: vtkdiffusionmodule.hh:65
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkdiffusionmodule.hh:92
VTK output module for quantities which make sense for models which assume thermal equilibrium.
Definition: vtkenergymodule.hh:66
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition: vtkenergymodule.hh:93
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:235
Declares the properties required for the NCP compositional multi-phase model.
The primary variable and equation indices for the compositional multi-phase NCP model.
Definition: ncpindices.hh:48
The type of the base class for all problems which use this model.
Definition: fvbaseproperties.hh:84
Type of object for specifying boundary conditions.
Definition: fvbaseproperties.hh:119
Enable diffusive fluxes?
Definition: multiphasebaseproperties.hh:79
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:149
Enumerations used by the model.
Definition: multiphasebaseproperties.hh:48
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
GetPropType< TypeTag, Scalar > type
Definition: ncpmodel.hh:141
The unmodified weight for the fugacity primary variables.
Definition: ncpproperties.hh:48
GetPropType< TypeTag, Scalar > type
Definition: ncpmodel.hh:127
The unmodified weight for the pressure primary variable.
Definition: ncpproperties.hh:42
GetPropType< TypeTag, Scalar > type
Definition: ncpmodel.hh:134
The weight for the saturation primary variables.
Definition: ncpproperties.hh:45
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
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:116
Define the type tag for the compositional NCP model.
Definition: ncpmodel.hh:72
std::tuple< MultiPhaseBaseModel > InheritsFrom
Definition: ncpmodel.hh:72