3.1-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
48public:
50 using FluidSystem = typename Traits::FluidSystem;
52 using FluidState = typename Traits::FluidState;
55
65 template<class ElementSolution, class Problem, class Element, class SubControlVolume>
66 void update(const ElementSolution &elemSol,
67 const Problem &problem,
68 const Element &element,
69 const SubControlVolume& scv)
70 {
71 ParentType::update(elemSol, problem, element, scv);
72
73 completeFluidState(elemSol, problem, element, scv, fluidState_);
74
75 typename FluidSystem::ParameterCache paramCache;
76 paramCache.updateAll(fluidState_);
77
78 for (unsigned int compIIdx = 0; compIIdx < ParentType::numFluidComponents(); ++compIIdx)
79 {
80 for (unsigned int compJIdx = 0; compJIdx < ParentType::numFluidComponents(); ++compJIdx)
81 {
82 // binary diffusion coefficients
83 if(compIIdx != compJIdx)
84 {
85 diffCoefficient_[compIIdx][compJIdx]
86 = FluidSystem::binaryDiffusionCoefficient(fluidState_,
87 paramCache,
88 0,
89 compIIdx,
90 compJIdx);
91 }
92 }
93 }
94 }
95
99 template<class ElementSolution, class Problem, class Element, class SubControlVolume>
100 static void completeFluidState(const ElementSolution& elemSol,
101 const Problem& problem,
102 const Element& element,
103 const SubControlVolume& scv,
105 {
106 fluidState.setTemperature(ParentType::temperature(elemSol, problem, element, scv));
107 fluidState.setPressure(0, elemSol[0][Indices::pressureIdx]);
108
109 // saturation in a single phase is always 1 and thus redundant
110 // to set. But since we use the fluid state shared by the
111 // immiscible multi-phase models, so we have to set it here...
112 fluidState.setSaturation(0, 1.0);
113
114 Scalar sumFracMinorComp = 0.0;
115
116 for(int compIdx = 1; compIdx < ParentType::numFluidComponents(); ++compIdx)
117 {
118 // temporary add 1.0 to remove spurious differences in mole fractions
119 // which are below the numerical accuracy
120 Scalar moleOrMassFraction = elemSol[0][Indices::conti0EqIdx+compIdx] + 1.0;
121 moleOrMassFraction = moleOrMassFraction - 1.0;
122 if(useMoles)
123 fluidState.setMoleFraction(0, compIdx, moleOrMassFraction);
124 else
125 fluidState.setMassFraction(0, compIdx, moleOrMassFraction);
126 sumFracMinorComp += moleOrMassFraction;
127 }
128 if(useMoles)
129 fluidState.setMoleFraction(0, 0, 1.0 - sumFracMinorComp);
130 else
131 fluidState.setMassFraction(0, 0, 1.0 - sumFracMinorComp);
132
133 typename FluidSystem::ParameterCache paramCache;
134 paramCache.updateAll(fluidState);
135
136 Scalar value = FluidSystem::density(fluidState, paramCache, 0);
137 fluidState.setDensity(0, value);
138
139 value = FluidSystem::molarDensity(fluidState, paramCache, 0);
140 fluidState.setMolarDensity(0, value);
141
142 value = FluidSystem::viscosity(fluidState, paramCache, 0);
143 fluidState.setViscosity(0, value);
144
145 // compute and set the enthalpy
146 const Scalar h = ParentType::enthalpy(fluidState, paramCache);
147 fluidState.setEnthalpy(0, h);
148 }
149
154 Scalar pressure(int phaseIdx = 0) const
155 { return fluidState_.pressure(0); }
156
161 Scalar density(int phaseIdx = 0) const
162 { return fluidState_.density(0); }
163
171 Scalar temperature() const
172 { return fluidState_.temperature(); }
173
178 Scalar effectiveViscosity() const
179 { return viscosity(); }
180
185 Scalar viscosity(int phaseIdx = 0) const
186 { return fluidState_.viscosity(0); }
187
194 Scalar massFraction(int phaseIdx, int compIdx) const
195 {
196 return fluidState_.massFraction(0, compIdx);
197 }
198
204 Scalar massFraction(int compIdx) const
205 {
206 return fluidState_.massFraction(0, compIdx);
207 }
208
215 Scalar moleFraction(int phaseIdx, int compIdx) const
216 {
217 return fluidState_.moleFraction(0, compIdx);
218 }
219
225 Scalar moleFraction(int compIdx) const
226 {
227 return fluidState_.moleFraction(0, compIdx);
228 }
229
233 Scalar molarDensity(int phaseIdx = 0) const
234 {
235 return fluidState_.molarDensity(0);
236 }
237
243 Scalar molarMass(int compIdx) const
244 {
245 return FluidSystem::molarMass(compIdx);
246 }
247
253 Scalar averageMolarMass(const int phaseIdx = 0) const
254 { return fluidState_.averageMolarMass(phaseIdx); }
255
262 Scalar diffusionCoefficient(int compIIdx, int compJIdx = 0) const
263 {
264 if (compIIdx == compJIdx)
265 DUNE_THROW(Dune::InvalidStateException, "Diffusion coefficient called for compIIdx = compJIdx");
266 return diffCoefficient_[compIIdx][compJIdx];
267 }
268
275 Scalar effectiveDiffusivity(int compIIdx, int compJIdx = 0) const
276 {
277 return diffusionCoefficient(compIIdx, compJIdx);
278 }
279
280
284 const FluidState& fluidState() const
285 { return fluidState_; }
286
287protected:
289 std::array<std::array<Scalar, ParentType::numFluidComponents()>, ParentType::numFluidComponents()> diffCoefficient_ = {};
290};
291
292} // end namespace Dumux
293
294#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
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:215
const FluidState & fluidState() const
Return the fluid state of the control volume.
Definition: freeflow/compositional/volumevariables.hh:284
Scalar diffusionCoefficient(int compIIdx, int compJIdx=0) const
Returns the diffusion coefficient .
Definition: freeflow/compositional/volumevariables.hh:262
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:194
std::array< std::array< Scalar, ParentType::numFluidComponents()>, ParentType::numFluidComponents()> diffCoefficient_
Definition: freeflow/compositional/volumevariables.hh:289
Scalar effectiveViscosity() const
Return the effective dynamic viscosity of the fluid within the control volume.
Definition: freeflow/compositional/volumevariables.hh:178
FluidState fluidState_
Definition: freeflow/compositional/volumevariables.hh:288
typename Traits::ModelTraits::Indices Indices
export the indices type
Definition: freeflow/compositional/volumevariables.hh:54
Scalar temperature() const
Return temperature inside the sub-control volume.
Definition: freeflow/compositional/volumevariables.hh:171
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:100
typename Traits::FluidState FluidState
export the fluid state type
Definition: freeflow/compositional/volumevariables.hh:52
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:66
Scalar viscosity(int phaseIdx=0) const
Return the dynamic viscosity of the fluid within the control volume.
Definition: freeflow/compositional/volumevariables.hh:185
Scalar pressure(int phaseIdx=0) const
Return the effective pressure of a given phase within the control volume.
Definition: freeflow/compositional/volumevariables.hh:154
Scalar molarMass(int compIdx) const
Returns the molar mass of a given component.
Definition: freeflow/compositional/volumevariables.hh:243
Scalar moleFraction(int compIdx) const
Returns the mole fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:225
Scalar density(int phaseIdx=0) const
Return the mass density of a given phase within the control volume.
Definition: freeflow/compositional/volumevariables.hh:161
Scalar massFraction(int compIdx) const
Returns the mass fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:204
Scalar effectiveDiffusivity(int compIIdx, int compJIdx=0) const
Returns the effective diffusion coefficient .
Definition: freeflow/compositional/volumevariables.hh:275
Scalar molarDensity(int phaseIdx=0) const
Returns the mass density of a given phase .
Definition: freeflow/compositional/volumevariables.hh:233
typename Traits::FluidSystem FluidSystem
export the underlying fluid system
Definition: freeflow/compositional/volumevariables.hh:50
Scalar averageMolarMass(const int phaseIdx=0) const
Returns the average molar mass of the fluid phase.
Definition: freeflow/compositional/volumevariables.hh:253
Definition: freeflow/volumevariables.hh:34