volumevariables.hh 8.71 KB
Newer Older
1
2
3
4
5
6
7
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
 *   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
10
11
12
13
14
15
16
17
18
19
20
 *   (at your option) any later version.                                     *
 *                                                                           *
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
 *   GNU General Public License for more details.                            *
 *                                                                           *
 *   You should have received a copy of the GNU General Public License       *
 *   along with this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
21
 * \ingroup TracerModel
22
 * \brief Quantities required by the tracer model in a control volume.
23
24
25
26
 */
#ifndef DUMUX_TRACER_VOLUME_VARIABLES_HH
#define DUMUX_TRACER_VOLUME_VARIABLES_HH

Dennis Gläser's avatar
Dennis Gläser committed
27
#include <cassert>
Bernd Flemisch's avatar
Bernd Flemisch committed
28
#include <array>
Dennis Gläser's avatar
Dennis Gläser committed
29
30
#include <type_traits>

31
#include <dune/common/std/type_traits.hh>
Bernd Flemisch's avatar
Bernd Flemisch committed
32

33
#include <dumux/porousmediumflow/volumevariables.hh>
Timo Koch's avatar
Timo Koch committed
34
#include <dumux/material/solidstates/updatesolidvolumefractions.hh>
35

36
namespace Dumux {
37

38
39
40
41
42
43
44
45
46
47
48
49
namespace Detail {
// helper structs and functions detecting if the user-defined spatial params class
// has user-specified functions saturation() for multi-phase tracer.
template <typename T, typename ...Ts>
using SaturationDetector = decltype(std::declval<T>().spatialParams().saturation(std::declval<Ts>()...));

template<class T, typename ...Args>
static constexpr bool hasSaturation()
{ return Dune::Std::is_detected<SaturationDetector, T, Args...>::value; }

} // end namespace Detail

50
51
52
53
54
/*!
 * \ingroup TracerModel
 * \brief Contains the quantities which are constant within a
 *        finite volume for the tracer model.
 */
55
56
template <class Traits>
class TracerVolumeVariables
57
: public PorousMediumFlowVolumeVariables<Traits>
58
{
59
    using ParentType = PorousMediumFlowVolumeVariables<Traits>;
60
61
    using Scalar = typename Traits::PrimaryVariables::value_type;
    static constexpr bool useMoles = Traits::ModelTraits::useMoles();
62
    using EffDiffModel = typename Traits::EffectiveDiffusivityModel;
63
    static constexpr int numFluidComps = ParentType::numFluidComponents();
64

65
public:
66
    //! Export the fluid system type
67
    using FluidSystem = typename Traits::FluidSystem;
68
    //! Export the solid state type
69
    using SolidState = typename Traits::SolidState;
70
71
    //! Export the indices
    using Indices = typename Traits::ModelTraits::Indices;
72
73

    /*!
74
     * \brief Updates all quantities for a given control volume.
75
76
77
78
79
80
     *
     * \param elemSol A vector containing all primary variables connected to the element
     * \param problem The object specifying the problem which ought to
     *                be simulated
     * \param element An element which contains part of the control volume
     * \param scv The sub-control volume
81
     */
82
83
    template<class ElemSol, class Problem, class Element, class Scv>
    void update(const ElemSol &elemSol,
84
85
                const Problem &problem,
                const Element &element,
86
                const Scv &scv)
87
88
89
90
    {
        // update parent type sets primary variables
        ParentType::update(elemSol, problem, element, scv);

91
        updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, numFluidComps);
92
93
94
95

        // the spatial params special to the tracer model
        fluidDensity_ = problem.spatialParams().fluidDensity(element, scv);
        fluidMolarMass_ = problem.spatialParams().fluidMolarMass(element, scv);
96
97
98

        if constexpr (Detail::hasSaturation<Problem, Element, Scv>())
            fluidSaturation_ = problem.spatialParams().saturation(element, scv);
99

100
        for (int compIdx = 0; compIdx < ParentType::numFluidComponents(); ++compIdx)
101
        {
102
            moleOrMassFraction_[compIdx] = this->priVars()[compIdx];
103
104
            diffCoeff_[compIdx] = FluidSystem::binaryDiffusionCoefficient(compIdx, problem, element, scv);
            effectiveDiffCoeff_[compIdx] = EffDiffModel::effectiveDiffusionCoefficient(*this, 0, 0, compIdx);
105
106
107
108
        }
    }

    /*!
109
     * \brief Returns the density \f$\mathrm{[kg/m^3]}\f$ the of the fluid phase.
110
111
     *
     * We always forward to the fluid state with the phaseIdx property (see class description).
112
     *
113
     * \param phaseIdx The phase index
114
     */
115
    Scalar density(int phaseIdx = 0) const
116
117
    { return fluidDensity_; }

118
    /*!
Kilian Weishaupt's avatar
Kilian Weishaupt committed
119
     * \brief Returns the average molar mass \f$\mathrm{[kg/mol]}\f$ of the fluid phase.
120
121
122
123
124
125
     *
     * \param phaseIdx The phase index
     */
    Scalar averageMolarMass(int phaseIdx = 0) const
    { return fluidMolarMass_; }

126
127
128
129
130
131
    /*!
     * \brief Returns the phase state for the control volume.
     */
    const SolidState &solidState() const
    { return solidState_; }

132
    /*!
133
     * \brief Returns the saturation.
134
135
     *
     * This method is here for compatibility reasons with other models. The saturation
136
137
138
139
     * is always 1.0 in a one-phasic context, if two-phases or richards are considered,
     * the spatialParams serve as way to pass the saturation from the main-file to the
     * volVars and then to the localresidual for the tracer model.

140
     * \param phaseIdx The phase index
141
     */
142
    Scalar saturation(int phaseIdx = 0) const
143
    { return fluidSaturation_ ; }
144
145

    /*!
146
     * \brief Returns the mobility.
147
148
149
     *
     * This method is here for compatibility reasons with other models. The mobility is always 1
     * for one-phasic models where the velocity field is given
150
     *
151
     * \param phaseIdx The phase index
152
     */
153
    Scalar mobility(int phaseIdx = 0) const
154
155
156
    { return 1.0; }

    /*!
157
     * \brief Returns the molar density \f$\mathrm{[mol/m^3]}\f$ the of the fluid phase.
158
     *
159
     * \param phaseIdx The phase index
160
     */
161
    Scalar molarDensity(int phaseIdx = 0) const
162
163
164
    { return fluidDensity_/fluidMolarMass_; }

    /*!
165
     * \brief Returns the mole fraction \f$\mathrm{[mol/mol]}\f$ of a component in the phase.
166
     *
167
     * \param phaseIdx The phase index
168
169
     * \param compIdx The index of the component
     */
170
    Scalar moleFraction(int phaseIdx, int compIdx) const
171
    { return useMoles ? moleOrMassFraction_[compIdx] : moleOrMassFraction_[compIdx]/FluidSystem::molarMass(compIdx)*fluidMolarMass_; }
172
173

    /*!
174
     * \brief Returns the mass fraction \f$\mathrm{[kg/kg]}\f$ of a component in the phase.
175
     *
176
     * \param phaseIdx The phase index
177
178
     * \param compIdx The index of the component
     */
179
    Scalar massFraction(int phaseIdx, int compIdx) const
180
    { return useMoles ? moleOrMassFraction_[compIdx]*FluidSystem::molarMass(compIdx)/fluidMolarMass_ : moleOrMassFraction_[compIdx]; }
181
182

    /*!
183
     * \brief Returns the concentration \f$\mathrm{[mol/m^3]}\f$  of a component in the phase.
184
     *
185
     * \param phaseIdx The phase index
186
187
     * \param compIdx The index of the component
     */
188
189
    Scalar molarity(int phaseIdx, int compIdx) const
    { return moleFraction(phaseIdx, compIdx)*molarDensity(); }
190

191
    /*!
192
     * \brief Returns the binary diffusion coefficients for a phase in \f$[m^2/s]\f$.
193
     */
194
    Scalar diffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
195
196
197
198
199
    {
        if (phaseIdx != compIIdx) std::swap(compIIdx, compJIdx);
        assert(phaseIdx == 0);
        assert(phaseIdx == compIIdx);
        return diffCoeff_[compJIdx]; }
200
201
202
203
204

    /*!
     * \brief Returns the effective diffusion coefficients for a phase in \f$[m^2/s]\f$.
     */
    Scalar effectiveDiffusionCoefficient(int phaseIdx, int compIIdx, int compJIdx) const
205
206
207
208
209
210
    {
        if (phaseIdx != compIIdx) std::swap(compIIdx, compJIdx);
        assert(phaseIdx == 0);
        assert(phaseIdx == compIIdx);
        return effectiveDiffCoeff_[compJIdx];
    }
211
212
213
214
215

    /*!
     * \brief Return the average porosity \f$\mathrm{[-]}\f$ within the control volume.
     */
    Scalar porosity() const
216
    { return solidState_.porosity(); }
217
218

protected:
219
    SolidState solidState_;
220
    Scalar fluidDensity_, fluidMolarMass_;
221
    Scalar fluidSaturation_ = 1.0;
222

223
224
225
    std::array<Scalar, numFluidComps> diffCoeff_;
    std::array<Scalar, numFluidComps> effectiveDiffCoeff_;
    std::array<Scalar, numFluidComps> moleOrMassFraction_;
226
227
228
229
230
};

} // end namespace Dumux

#endif