3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/tracer/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 *****************************************************************************/
24#ifndef DUMUX_TRACER_VOLUME_VARIABLES_HH
25#define DUMUX_TRACER_VOLUME_VARIABLES_HH
26
27#include <cassert>
28#include <array>
29#include <type_traits>
30
31#include <dune/common/std/type_traits.hh>
32
35
36namespace Dumux {
37
38namespace Detail {
39// helper structs and functions detecting if the user-defined spatial params class
40// has user-specified functions saturation() for multi-phase tracer.
41template <typename T, typename ...Ts>
42using SaturationDetector = decltype(std::declval<T>().spatialParams().saturation(std::declval<Ts>()...));
43
44template<class T, typename ...Args>
45static constexpr bool hasSaturation()
46{ return Dune::Std::is_detected<SaturationDetector, T, Args...>::value; }
47
48} // end namespace Detail
49
55template <class Traits>
58{
60 using Scalar = typename Traits::PrimaryVariables::value_type;
61 static constexpr bool useMoles = Traits::ModelTraits::useMoles();
62 using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
63 static constexpr int numFluidComps = ParentType::numFluidComponents();
64
65public:
67 using FluidSystem = typename Traits::FluidSystem;
69 using SolidState = typename Traits::SolidState;
71 using Indices = typename Traits::ModelTraits::Indices;
72
82 template<class ElemSol, class Problem, class Element, class Scv>
83 void update(const ElemSol &elemSol,
84 const Problem &problem,
85 const Element &element,
86 const Scv &scv)
87 {
88 // update parent type sets primary variables
89 ParentType::update(elemSol, problem, element, scv);
90
91 updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
92
93 // the spatial params special to the tracer model
94 fluidDensity_ = problem.spatialParams().fluidDensity(element, scv);
95 fluidMolarMass_ = problem.spatialParams().fluidMolarMass(element, scv);
96
97 if constexpr (Detail::hasSaturation<Problem, Element, Scv>())
98 fluidSaturation_ = problem.spatialParams().saturation(element, scv);
99
100 for (int compIdx = 0; compIdx < ParentType::numFluidComponents(); ++compIdx)
101 {
102 moleOrMassFraction_[compIdx] = this->priVars()[compIdx];
103 diffCoeff_[compIdx] = FluidSystem::binaryDiffusionCoefficient(compIdx, problem, element, scv);
104 effectiveDiffCoeff_[compIdx] = EffDiffModel::effectiveDiffusionCoefficient(*this, 0, 0, compIdx);
105 }
106 }
107
115 Scalar density(int phaseIdx = 0) const
116 { return fluidDensity_; }
117
123 Scalar averageMolarMass(int phaseIdx = 0) const
124 { return fluidMolarMass_; }
125
129 const SolidState &solidState() const
130 { return solidState_; }
131
142 Scalar saturation(int phaseIdx = 0) const
143 { return fluidSaturation_ ; }
144
153 Scalar mobility(int phaseIdx = 0) const
154 { return 1.0; }
155
161 Scalar molarDensity(int phaseIdx = 0) const
163
170 Scalar moleFraction(int phaseIdx, int compIdx) const
171 { return useMoles ? moleOrMassFraction_[compIdx] : moleOrMassFraction_[compIdx]/FluidSystem::molarMass(compIdx)*fluidMolarMass_; }
172
179 Scalar massFraction(int phaseIdx, int compIdx) const
180 { return useMoles ? moleOrMassFraction_[compIdx]*FluidSystem::molarMass(compIdx)/fluidMolarMass_ : moleOrMassFraction_[compIdx]; }
181
188 Scalar molarity(int phaseIdx, int compIdx) const
189 { return moleFraction(phaseIdx, compIdx)*molarDensity(); }
190
194 Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
195 {
196 if (phaseIdx != compIIdx) std::swap(compIIdx, compJIdx);
197 assert(phaseIdx == 0);
198 assert(phaseIdx == compIIdx);
199 return diffCoeff_[compJIdx]; }
200
204 Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
205 {
206 if (phaseIdx != compIIdx) std::swap(compIIdx, compJIdx);
207 assert(phaseIdx == 0);
208 assert(phaseIdx == compIIdx);
209 return effectiveDiffCoeff_[compJIdx];
210 }
211
215 Scalar porosity() const
216 { return solidState_.porosity(); }
217
218protected:
221 Scalar fluidSaturation_ = 1.0;
222
223 std::array<Scalar, numFluidComps> diffCoeff_;
224 std::array<Scalar, numFluidComps> effectiveDiffCoeff_;
225 std::array<Scalar, numFluidComps> moleOrMassFraction_;
226};
227
228} // end namespace Dumux
229
230#endif
Update the solid volume fractions (inert and reacitve) and set them in the solidstate.
void updateSolidVolumeFractions(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv, SolidState &solidState, const int solidVolFracOffset)
update the solid volume fractions (inert and reacitve) and set them in the solidstate
Definition: updatesolidvolumefractions.hh:36
Definition: adapt.hh:29
decltype(std::declval< T >().spatialParams().saturation(std::declval< Ts >()...)) SaturationDetector
Definition: porousmediumflow/tracer/volumevariables.hh:42
static constexpr bool hasSaturation()
Definition: porousmediumflow/tracer/volumevariables.hh:45
Contains the quantities which are constant within a finite volume for the tracer model.
Definition: porousmediumflow/tracer/volumevariables.hh:58
Scalar averageMolarMass(int phaseIdx=0) const
Returns the average molar mass of the fluid phase.
Definition: porousmediumflow/tracer/volumevariables.hh:123
typename Traits::FluidSystem FluidSystem
Export the fluid system type.
Definition: porousmediumflow/tracer/volumevariables.hh:67
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the mole fraction of a component in the phase.
Definition: porousmediumflow/tracer/volumevariables.hh:170
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of a component in the phase.
Definition: porousmediumflow/tracer/volumevariables.hh:179
std::array< Scalar, numFluidComps > effectiveDiffCoeff_
Definition: porousmediumflow/tracer/volumevariables.hh:224
Scalar molarDensity(int phaseIdx=0) const
Returns the molar density the of the fluid phase.
Definition: porousmediumflow/tracer/volumevariables.hh:161
std::array< Scalar, numFluidComps > moleOrMassFraction_
Definition: porousmediumflow/tracer/volumevariables.hh:225
typename Traits::SolidState SolidState
Export the solid state type.
Definition: porousmediumflow/tracer/volumevariables.hh:69
Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the effective diffusion coefficients for a phase in .
Definition: porousmediumflow/tracer/volumevariables.hh:204
Scalar fluidDensity_
Definition: porousmediumflow/tracer/volumevariables.hh:220
SolidState solidState_
Definition: porousmediumflow/tracer/volumevariables.hh:219
Scalar density(int phaseIdx=0) const
Returns the density the of the fluid phase.
Definition: porousmediumflow/tracer/volumevariables.hh:115
const SolidState & solidState() const
Returns the phase state for the control volume.
Definition: porousmediumflow/tracer/volumevariables.hh:129
Scalar mobility(int phaseIdx=0) const
Returns the mobility.
Definition: porousmediumflow/tracer/volumevariables.hh:153
Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
Returns the binary diffusion coefficients for a phase in .
Definition: porousmediumflow/tracer/volumevariables.hh:194
Scalar saturation(int phaseIdx=0) const
Returns the saturation.
Definition: porousmediumflow/tracer/volumevariables.hh:142
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/tracer/volumevariables.hh:83
Scalar porosity() const
Return the average porosity within the control volume.
Definition: porousmediumflow/tracer/volumevariables.hh:215
typename Traits::ModelTraits::Indices Indices
Export the indices.
Definition: porousmediumflow/tracer/volumevariables.hh:71
Scalar fluidSaturation_
Definition: porousmediumflow/tracer/volumevariables.hh:221
Scalar fluidMolarMass_
Definition: porousmediumflow/tracer/volumevariables.hh:220
Scalar molarity(int phaseIdx, int compIdx) const
Returns the concentration of a component in the phase.
Definition: porousmediumflow/tracer/volumevariables.hh:188
std::array< Scalar, numFluidComps > diffCoeff_
Definition: porousmediumflow/tracer/volumevariables.hh:223
The isothermal base class.
Definition: porousmediumflow/volumevariables.hh:40
static constexpr int numFluidComponents()
Return number of components considered by the model.
Definition: porousmediumflow/volumevariables.hh:52
const PrimaryVariables & priVars() const
Returns the vector of primary variables.
Definition: porousmediumflow/volumevariables.hh:76
void update(const ElemSol &elemSol, const Problem &problem, const Element &element, const Scv &scv)
Updates all quantities for a given control volume.
Definition: porousmediumflow/volumevariables.hh:64
Base class for the model specific class which provides access to all volume averaged quantities.