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<VtkDiffusion,
76} // namespace TTag
77
79template<class TypeTag>
80struct LocalResidual<TypeTag, TTag::NcpModel> { using type = NcpLocalResidual<TypeTag>; };
81
83template<class TypeTag>
84struct NewtonMethod<TypeTag, TTag::NcpModel> { using type = NcpNewtonMethod<TypeTag>; };
85
87template<class TypeTag>
88struct Model<TypeTag, TTag::NcpModel> { using type = NcpModel<TypeTag>; };
89
91template<class TypeTag>
92struct BaseProblem<TypeTag, TTag::NcpModel> { using type = MultiPhaseBaseProblem<TypeTag>; };
93
95template<class TypeTag>
96struct EnableEnergy<TypeTag, TTag::NcpModel> { static constexpr bool value = false; };
97
99template<class TypeTag>
100struct EnableDiffusion<TypeTag, TTag::NcpModel> { static constexpr bool value = false; };
101
103template<class TypeTag>
104struct RateVector<TypeTag, TTag::NcpModel> { using type = NcpRateVector<TypeTag>; };
105
107template<class TypeTag>
109
111template<class TypeTag>
112struct PrimaryVariables<TypeTag, TTag::NcpModel> { using type = NcpPrimaryVariables<TypeTag>; };
113
115template<class TypeTag>
117
119template<class TypeTag>
121
123template<class TypeTag>
124struct Indices<TypeTag, TTag::NcpModel> { using type = NcpIndices<TypeTag, 0>; };
125
127template<class TypeTag>
128struct NcpPressureBaseWeight<TypeTag, TTag::NcpModel>
129{
131 static constexpr type value = 1.0;
132};
134template<class TypeTag>
135struct NcpSaturationsBaseWeight<TypeTag, TTag::NcpModel>
136{
138 static constexpr type value = 1.0;
139};
141template<class TypeTag>
142struct NcpFugacitiesBaseWeight<TypeTag, TTag::NcpModel>
143{
145 static constexpr type value = 1.0e-6;
146};
147
148} // namespace Opm::Properties
149
150namespace Opm {
151
226template <class TypeTag>
228 : public MultiPhaseBaseModel<TypeTag>
229{
230 using ParentType = MultiPhaseBaseModel<TypeTag>;
231
238
239 enum { numPhases = FluidSystem::numPhases };
240 enum { numComponents = FluidSystem::numComponents };
241 enum { fugacity0Idx = Indices::fugacity0Idx };
242 enum { pressure0Idx = Indices::pressure0Idx };
243 enum { saturation0Idx = Indices::saturation0Idx };
244 enum { conti0EqIdx = Indices::conti0EqIdx };
245 enum { ncp0EqIdx = Indices::ncp0EqIdx };
246 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
247 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
248
249 using ComponentVector = Dune::FieldVector<Scalar, numComponents>;
250
251 using Toolbox = MathToolbox<Evaluation>;
252
255
256public:
257 NcpModel(Simulator& simulator)
258 : ParentType(simulator)
259 {}
260
264 static void registerParameters()
265 {
267
268 DiffusionModule::registerParameters();
269 EnergyModule::registerParameters();
270
271 // register runtime parameters of the VTK output modules
273
274 if (enableDiffusion)
276
277 if (enableEnergy)
279 }
280
285 {
286 ParentType::finishInit();
287
288 minActivityCoeff_.resize(this->numGridDof());
289 std::fill(minActivityCoeff_.begin(), minActivityCoeff_.end(), 1.0);
290 }
291
293 {
294 ParentType::adaptGrid();
295 minActivityCoeff_.resize(this->numGridDof());
296 }
297
301 static std::string name()
302 { return "ncp"; }
303
307 std::string primaryVarName(unsigned pvIdx) const
308 {
309 std::string s;
310 if (!(s = EnergyModule::primaryVarName(pvIdx)).empty())
311 return s;
312
313 std::ostringstream oss;
314 if (pvIdx == pressure0Idx)
315 oss << "pressure_" << FluidSystem::phaseName(/*phaseIdx=*/0);
316 else if (saturation0Idx <= pvIdx && pvIdx < saturation0Idx + (numPhases - 1))
317 oss << "saturation_" << FluidSystem::phaseName(/*phaseIdx=*/pvIdx - saturation0Idx);
318 else if (fugacity0Idx <= pvIdx && pvIdx < fugacity0Idx + numComponents)
319 oss << "fugacity^" << FluidSystem::componentName(pvIdx - fugacity0Idx);
320 else
321 assert(false);
322
323 return oss.str();
324 }
325
329 std::string eqName(unsigned eqIdx) const
330 {
331 std::string s;
332 if (!(s = EnergyModule::eqName(eqIdx)).empty())
333 return s;
334
335 std::ostringstream oss;
336 if (conti0EqIdx <= eqIdx && eqIdx < conti0EqIdx + numComponents)
337 oss << "continuity^" << FluidSystem::componentName(eqIdx - conti0EqIdx);
338 else if (ncp0EqIdx <= eqIdx && eqIdx < ncp0EqIdx + numPhases)
339 oss << "ncp_" << FluidSystem::phaseName(/*phaseIdx=*/eqIdx - ncp0EqIdx);
340 else
341 assert(false);
342
343 return oss.str();
344 }
345
350 {
351 ParentType::updateBegin();
352
353 // find the a reference pressure. The first degree of freedom
354 // might correspond to non-interior entities which would lead
355 // to an undefined value, so we have to iterate...
356 for (unsigned dofIdx = 0; dofIdx < this->numGridDof(); ++ dofIdx) {
357 if (this->isLocalDof(dofIdx)) {
359 this->solution(/*timeIdx=*/0)[dofIdx][/*pvIdx=*/Indices::pressure0Idx];
360 break;
361 }
362 }
363 }
364
368 void updatePVWeights(const ElementContext& elemCtx) const
369 {
370 for (unsigned dofIdx = 0; dofIdx < elemCtx.numDof(/*timeIdx=*/0); ++dofIdx) {
371 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
372
373 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
374 minActivityCoeff_[globalIdx][compIdx] = 1e100;
375 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
376 const auto& fs = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0).fluidState();
377
378 minActivityCoeff_[globalIdx][compIdx] =
379 std::min(minActivityCoeff_[globalIdx][compIdx],
380 Toolbox::value(fs.fugacityCoefficient(phaseIdx, compIdx))
381 * Toolbox::value(fs.pressure(phaseIdx)));
382 Valgrind::CheckDefined(minActivityCoeff_[globalIdx][compIdx]);
383 }
384 if (minActivityCoeff_[globalIdx][compIdx] <= 0)
385 throw NumericalProblem("The minimum activity coefficient for component "+std::to_string(compIdx)
386 +" on DOF "+std::to_string(globalIdx)+" is negative or zero!");
387 }
388 }
389 }
390
394 Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
395 {
396 Scalar tmp = EnergyModule::primaryVarWeight(*this, globalDofIdx, pvIdx);
397 Scalar result;
398 if (tmp > 0)
399 // energy related quantity
400 result = tmp;
401 else if (fugacity0Idx <= pvIdx && pvIdx < fugacity0Idx + numComponents) {
402 // component fugacity
403 unsigned compIdx = pvIdx - fugacity0Idx;
404 assert(compIdx <= numComponents);
405
406 Valgrind::CheckDefined(minActivityCoeff_[globalDofIdx][compIdx]);
407 static const Scalar fugacityBaseWeight =
408 getPropValue<TypeTag, Properties::NcpFugacitiesBaseWeight>();
409 result = fugacityBaseWeight / minActivityCoeff_[globalDofIdx][compIdx];
410 }
411 else if (Indices::pressure0Idx == pvIdx) {
412 static const Scalar pressureBaseWeight = getPropValue<TypeTag, Properties::NcpPressureBaseWeight>();
413 result = pressureBaseWeight / referencePressure_;
414 }
415 else {
416#ifndef NDEBUG
417 unsigned phaseIdx = pvIdx - saturation0Idx;
418 assert(phaseIdx < numPhases - 1);
419#endif
420
421 // saturation
422 static const Scalar saturationsBaseWeight =
423 getPropValue<TypeTag, Properties::NcpSaturationsBaseWeight>();
424 result = saturationsBaseWeight;
425 }
426
427 assert(std::isfinite(result));
428 assert(result > 0);
429
430 return result;
431 }
432
436 Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
437 {
438 Scalar tmp = EnergyModule::eqWeight(*this, globalDofIdx, eqIdx);
439 if (tmp > 0)
440 // an energy related equation
441 return tmp;
442 // an NCP
443 else if (ncp0EqIdx <= eqIdx && eqIdx < Indices::ncp0EqIdx + numPhases)
444 return 1.0;
445
446 // a mass conservation equation
447 unsigned compIdx = eqIdx - Indices::conti0EqIdx;
448 assert(compIdx <= numComponents);
449
450 // make all kg equal
451 return FluidSystem::molarMass(compIdx);
452 }
453
461 Scalar minActivityCoeff(unsigned globalDofIdx, unsigned compIdx) const
462 { return minActivityCoeff_[globalDofIdx][compIdx]; }
463
468 {
470
471 this->addOutputModule(new VtkCompositionModule<TypeTag>(this->simulator_));
472 if (enableDiffusion)
473 this->addOutputModule(new VtkDiffusionModule<TypeTag>(this->simulator_));
474 if (enableEnergy)
475 this->addOutputModule(new VtkEnergyModule<TypeTag>(this->simulator_));
476 }
477
478 mutable Scalar referencePressure_;
479 mutable std::vector<ComponentVector> minActivityCoeff_;
480};
481
482} // namespace Opm
483
484#endif
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:47
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
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:57
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:229
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:368
std::string eqName(unsigned eqIdx) const
Given an equation index, return a human readable name.
Definition: ncpmodel.hh:329
std::string primaryVarName(unsigned pvIdx) const
Given an primary variable index, return a human readable name.
Definition: ncpmodel.hh:307
NcpModel(Simulator &simulator)
Definition: ncpmodel.hh:257
static std::string name()
Definition: ncpmodel.hh:301
Scalar referencePressure_
Definition: ncpmodel.hh:478
Scalar primaryVarWeight(unsigned globalDofIdx, unsigned pvIdx) const
Returns the relative weight of a primary variable for calculating relative errors.
Definition: ncpmodel.hh:394
void adaptGrid()
Definition: ncpmodel.hh:292
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:461
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: ncpmodel.hh:264
void finishInit()
Apply the initial conditions to the model.
Definition: ncpmodel.hh:284
void updateBegin()
Called by the update() method before it tries to apply the newton method. This is primary a hook whic...
Definition: ncpmodel.hh:349
Scalar eqWeight(unsigned globalDofIdx, unsigned eqIdx) const
Returns the relative weight of an equation.
Definition: ncpmodel.hh:436
void registerOutputModules_()
Definition: ncpmodel.hh:467
std::vector< ComponentVector > minActivityCoeff_
Definition: ncpmodel.hh:479
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: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
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 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:101
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
GetPropType< TypeTag, Scalar > type
Definition: ncpmodel.hh:144
The unmodified weight for the fugacity primary variables.
Definition: ncpproperties.hh:48
GetPropType< TypeTag, Scalar > type
Definition: ncpmodel.hh:130
The unmodified weight for the pressure primary variable.
Definition: ncpproperties.hh:42
GetPropType< TypeTag, Scalar > type
Definition: ncpmodel.hh:137
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:147
Vector containing volumetric or areal rates of quantities.
Definition: fvbaseproperties.hh:133
Definition: multiphasebasemodel.hh:57
Define the type tag for the compositional NCP model.
Definition: ncpmodel.hh:72
std::tuple< VtkDiffusion, VtkEnergy, VtkComposition, MultiPhaseBaseModel > InheritsFrom
Definition: ncpmodel.hh:75
Definition: vtkcompositionmodule.hh:43
Definition: vtkdiffusionmodule.hh:46
Definition: vtkenergymodule.hh:43