3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
initialconditionhelper.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_MPNC_INITIALCONDITION_HELPER_HH
25#define DUMUX_MPNC_INITIALCONDITION_HELPER_HH
26
29
31
32namespace Dumux {
33
34namespace MPNCInitialConditions {
35
38
39} // namespace MPNCInitialConditions
40
41template<class PrimaryVariables, class FluidSystem, class ModelTraits>
43{
44 using Scalar = typename PrimaryVariables::value_type;
47
48public:
49 template<class FluidState>
50 void solve(FluidState& fs, const AllPhasesPresent& allPhases) const
51 {
52 typename FluidSystem::ParameterCache paramCache;
54 paramCache,
55 allPhases.refPhaseIdx);
56 }
57
58 template<class FluidState>
59 void solve(FluidState& fs, const NotAllPhasesPresent& notAllPhases) const
60 {
61 typename FluidSystem::ParameterCache paramCache;
63 paramCache,
64 notAllPhases.refPhaseIdx);
65 }
66
67 template<class FluidState>
68 PrimaryVariables getPrimaryVariables(const FluidState& fs) const
69 {
70 PrimaryVariables priVars(0.0);
71
72 static constexpr auto numComponents = FluidSystem::numComponents;
73 static constexpr auto numPhases = FluidSystem::numPhases;
74 static constexpr auto fug0Idx = ModelTraits::Indices::fug0Idx;
75 static constexpr auto s0Idx = ModelTraits::Indices::s0Idx;
76 static constexpr auto p0Idx = ModelTraits::Indices::p0Idx;
77
78 // all N component fugacities
79 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
80 priVars[fug0Idx + compIdx] = fs.fugacity(0, compIdx);
81
82 // first M - 1 saturations
83 for (int phaseIdx = 0; phaseIdx < numPhases - 1; ++phaseIdx)
84 priVars[s0Idx + phaseIdx] = fs.saturation(phaseIdx);
85
86 static constexpr auto pressureFormulation = ModelTraits::pressureFormulation();
87 if (pressureFormulation == MpNcPressureFormulation::mostWettingFirst)
88 priVars[p0Idx] = fs.pressure(/*phaseIdx=*/0);
89 else if (pressureFormulation == MpNcPressureFormulation::leastWettingFirst)
90 priVars[p0Idx] = fs.pressure(numPhases-1);
91 else
92 DUNE_THROW(Dune::InvalidStateException,"unknown pressure formulation");
93
94 return priVars;
95 }
96};
97
98} // end namespace Dumux
99
100#endif
Computes all quantities of a generic fluid state if a reference phase has been specified.
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Enumeration of the formulations accepted by the MpNc model.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
static void solve(FluidState &fluidState, ParameterCache &paramCache, int refPhaseIdx)
Computes all quantities of a generic fluid state if a reference phase has been specified.
Definition: computefromreferencephase.hh:111
static void solve(FluidState &fluidState, ParameterCache &paramCache, int knownPhaseIdx=0)
Computes the composition of all phases of a N-phase, N-component fluid system assuming that all N pha...
Definition: misciblemultiphasecomposition.hh:81
Definition: initialconditionhelper.hh:36
int refPhaseIdx
Definition: initialconditionhelper.hh:36
Definition: initialconditionhelper.hh:37
int refPhaseIdx
Definition: initialconditionhelper.hh:37
Definition: initialconditionhelper.hh:43
PrimaryVariables getPrimaryVariables(const FluidState &fs) const
Definition: initialconditionhelper.hh:68
void solve(FluidState &fs, const AllPhasesPresent &allPhases) const
Definition: initialconditionhelper.hh:50
void solve(FluidState &fs, const NotAllPhasesPresent &notAllPhases) const
Definition: initialconditionhelper.hh:59