3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/compositional/volumevariables.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 *****************************************************************************/
25#ifndef DUMUX_FREEFLOW_NC_VOLUMEVARIABLES_HH
26#define DUMUX_FREEFLOW_NC_VOLUMEVARIABLES_HH
27
28#include <array>
29#include <dune/common/exceptions.hh>
31
32namespace Dumux {
33
38template <class Traits>
39class FreeflowNCVolumeVariables : public FreeFlowVolumeVariables< Traits, FreeflowNCVolumeVariables<Traits> >
40{
43
44 using Scalar = typename Traits::PrimaryVariables::value_type;
45
46 static constexpr bool useMoles = Traits::ModelTraits::useMoles();
47 static constexpr int numFluidComps = ParentType::numFluidComponents();
48 using DiffusionCoefficients = typename Traits::DiffusionType::DiffusionCoefficientsContainer;
49
50public:
52 using FluidSystem = typename Traits::FluidSystem;
54 using FluidState = typename Traits::FluidState;
56 using Indices = typename Traits::ModelTraits::Indices;
57
67 template<class ElementSolution, class Problem, class Element, class SubControlVolume>
68 void update(const ElementSolution &elemSol,
69 const Problem &problem,
70 const Element &element,
71 const SubControlVolume& scv)
72 {
73 ParentType::update(elemSol, problem, element, scv);
74
75 completeFluidState(elemSol, problem, element, scv, fluidState_);
76
77 typename FluidSystem::ParameterCache paramCache;
78 paramCache.updateAll(fluidState_);
79
80 auto getDiffusionCoefficient = [&](int phaseIdx, int compIIdx, int compJIdx)
81 {
82 return FluidSystem::binaryDiffusionCoefficient(this->fluidState_,
83 paramCache,
84 0,
85 compIIdx,
86 compJIdx);
87 };
88
89 diffCoefficient_.update(getDiffusionCoefficient);
90 }
91
95 template<class ElementSolution, class Problem, class Element, class SubControlVolume>
96 static void completeFluidState(const ElementSolution& elemSol,
97 const Problem& problem,
98 const Element& element,
99 const SubControlVolume& scv,
101 {
102 fluidState.setTemperature(ParentType::temperature(elemSol, problem, element, scv));
103 fluidState.setPressure(0, elemSol[0][Indices::pressureIdx]);
104
105 // saturation in a single phase is always 1 and thus redundant
106 // to set. But since we use the fluid state shared by the
107 // immiscible multi-phase models, so we have to set it here...
108 fluidState.setSaturation(0, 1.0);
109
110 Scalar sumFracMinorComp = 0.0;
111
112 for(int compIdx = 1; compIdx < ParentType::numFluidComponents(); ++compIdx)
113 {
114 // temporary add 1.0 to remove spurious differences in mole fractions
115 // which are below the numerical accuracy
116 Scalar moleOrMassFraction = elemSol[0][Indices::conti0EqIdx+compIdx] + 1.0;
117 moleOrMassFraction = moleOrMassFraction - 1.0;
118 if(useMoles)
119 fluidState.setMoleFraction(0, compIdx, moleOrMassFraction);
120 else
121 fluidState.setMassFraction(0, compIdx, moleOrMassFraction);
122 sumFracMinorComp += moleOrMassFraction;
123 }
124 if(useMoles)
125 fluidState.setMoleFraction(0, 0, 1.0 - sumFracMinorComp);
126 else
127 fluidState.setMassFraction(0, 0, 1.0 - sumFracMinorComp);
128
129 typename FluidSystem::ParameterCache paramCache;
130 paramCache.updateAll(fluidState);
131
132 Scalar value = FluidSystem::density(fluidState, paramCache, 0);
133 fluidState.setDensity(0, value);
134
135 value = FluidSystem::molarDensity(fluidState, paramCache, 0);
136 fluidState.setMolarDensity(0, value);
137
138 value = FluidSystem::viscosity(fluidState, paramCache, 0);
139 fluidState.setViscosity(0, value);
140
141 // compute and set the enthalpy
142 const Scalar h = ParentType::enthalpy(fluidState, paramCache);
143 fluidState.setEnthalpy(0, h);
144 }
145
150 Scalar pressure(int phaseIdx = 0) const
151 { return fluidState_.pressure(0); }
152
157 Scalar density(int phaseIdx = 0) const
158 { return fluidState_.density(0); }
159
167 Scalar temperature() const
168 { return fluidState_.temperature(); }
169
174 Scalar effectiveViscosity() const
175 { return viscosity(); }
176
181 Scalar viscosity(int phaseIdx = 0) const
182 { return fluidState_.viscosity(0); }
183
190 Scalar massFraction(int phaseIdx, int compIdx) const
191 {
192 return fluidState_.massFraction(0, compIdx);
193 }
194
200 Scalar massFraction(int compIdx) const
201 {
202 return fluidState_.massFraction(0, compIdx);
203 }
204
211 Scalar moleFraction(int phaseIdx, int compIdx) const
212 {
213 return fluidState_.moleFraction(0, compIdx);
214 }
215
221 Scalar moleFraction(int compIdx) const
222 {
223 return fluidState_.moleFraction(0, compIdx);
224 }
225
229 Scalar molarDensity(int phaseIdx = 0) const
230 {
231 return fluidState_.molarDensity(0);
232 }
233
239 Scalar molarMass(int compIdx) const
240 {
241 return FluidSystem::molarMass(compIdx);
242 }
243
249 Scalar averageMolarMass(const int phaseIdx = 0) const
250 { return fluidState_.averageMolarMass(phaseIdx); }
251
258 [[deprecated("Will be removed after release 3.2. Use diffusionCoefficient(phaseIdx, compIIdx, compJIdx)!")]]
259 Scalar diffusionCoefficient(int compIIdx, int compJIdx = 0) const
260 {
261 if (compIIdx == compJIdx)
262 DUNE_THROW(Dune::InvalidStateException, "Diffusion coefficient called for compIIdx = compJIdx");
263 return diffCoefficient_(0, compIIdx, compJIdx);
264 }
265
269 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
270 { return diffCoefficient_(0, compIIdx, compJIdx); }
271
278 [[deprecated("Signature deprecated. Use effectiveDiffusionCoefficient(phaseIdx, compIIdx, compJIdx)!")]]
279 Scalar effectiveDiffusivity(int compIIdx, int compJIdx = 0) const
280 { return diffusionCoefficient(compIIdx, compJIdx); }
281
285 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
286 { return diffusionCoefficient(0, compIIdx, compJIdx); }
287
291 const FluidState& fluidState() const
292 { return fluidState_; }
293
294protected:
296 // Binary diffusion coefficient
297 DiffusionCoefficients diffCoefficient_;
298};
299
300} // end namespace Dumux
301
302#endif
Definition: adapt.hh:29
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:83
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Volume variables for the single-phase, multi-component free-flow model.
Definition: freeflow/compositional/volumevariables.hh:40
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the mole fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:211
const FluidState & fluidState() const
Return the fluid state of the control volume.
Definition: freeflow/compositional/volumevariables.hh:291
Scalar diffusionCoefficient(int compIIdx, int compJIdx=0) const
Returns the diffusion coefficient .
Definition: freeflow/compositional/volumevariables.hh:259
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:190
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: freeflow/compositional/volumevariables.hh:269
Scalar effectiveViscosity() const
Return the effective dynamic viscosity of the fluid within the control volume.
Definition: freeflow/compositional/volumevariables.hh:174
FluidState fluidState_
Definition: freeflow/compositional/volumevariables.hh:295
typename Traits::ModelTraits::Indices Indices
export the indices type
Definition: freeflow/compositional/volumevariables.hh:56
Scalar temperature() const
Return temperature inside the sub-control volume.
Definition: freeflow/compositional/volumevariables.hh:167
static void completeFluidState(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv, FluidState &fluidState)
Update the fluid state.
Definition: freeflow/compositional/volumevariables.hh:96
typename Traits::FluidState FluidState
export the fluid state type
Definition: freeflow/compositional/volumevariables.hh:54
void update(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv)
Update all quantities for a given control volume.
Definition: freeflow/compositional/volumevariables.hh:68
Scalar viscosity(int phaseIdx=0) const
Return the dynamic viscosity of the fluid within the control volume.
Definition: freeflow/compositional/volumevariables.hh:181
Scalar pressure(int phaseIdx=0) const
Return the effective pressure of a given phase within the control volume.
Definition: freeflow/compositional/volumevariables.hh:150
Scalar molarMass(int compIdx) const
Returns the molar mass of a given component.
Definition: freeflow/compositional/volumevariables.hh:239
Scalar moleFraction(int compIdx) const
Returns the mole fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:221
Scalar density(int phaseIdx=0) const
Return the mass density of a given phase within the control volume.
Definition: freeflow/compositional/volumevariables.hh:157
Scalar massFraction(int compIdx) const
Returns the mass fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:200
DiffusionCoefficients diffCoefficient_
Definition: freeflow/compositional/volumevariables.hh:297
Scalar effectiveDiffusivity(int compIIdx, int compJIdx=0) const
Returns the effective diffusion coefficient .
Definition: freeflow/compositional/volumevariables.hh:279
Scalar molarDensity(int phaseIdx=0) const
Returns the mass density of a given phase .
Definition: freeflow/compositional/volumevariables.hh:229
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: freeflow/compositional/volumevariables.hh:285
typename Traits::FluidSystem FluidSystem
export the underlying fluid system
Definition: freeflow/compositional/volumevariables.hh:52
Scalar averageMolarMass(const int phaseIdx=0) const
Returns the average molar mass of the fluid phase.
Definition: freeflow/compositional/volumevariables.hh:249
Definition: freeflow/volumevariables.hh:34