volumevariables.hh 6.92 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
23
24
25
26
 * \brief Quantities required by the tracer model in a control volume
 */
#ifndef DUMUX_TRACER_VOLUME_VARIABLES_HH
#define DUMUX_TRACER_VOLUME_VARIABLES_HH

Bernd Flemisch's avatar
Bernd Flemisch committed
27
28
#include <array>

29
#include <dumux/porousmediumflow/volumevariables.hh>
Timo Koch's avatar
Timo Koch committed
30
#include <dumux/material/solidstates/updatesolidvolumefractions.hh>
31

32
namespace Dumux {
33
34
35
36
37
38

/*!
 * \ingroup TracerModel
 * \brief Contains the quantities which are constant within a
 *        finite volume for the tracer model.
 */
39
40
template <class Traits>
class TracerVolumeVariables
41
: public PorousMediumFlowVolumeVariables<Traits>
42
{
43
    using ParentType = PorousMediumFlowVolumeVariables<Traits>;
44
45
    using Scalar = typename Traits::PrimaryVariables::value_type;
    static constexpr bool useMoles = Traits::ModelTraits::useMoles();
46
47

public:
48
49
    //! export fluid system type
    using FluidSystem = typename Traits::FluidSystem;
50
    using SolidState = typename Traits::SolidState;
51
52

    /*!
53
54
55
56
57
58
59
     * \brief Update all quantities for a given control volume
     *
     * \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
60
     */
61
62
    template<class ElemSol, class Problem, class Element, class Scv>
    void update(const ElemSol &elemSol,
63
64
                const Problem &problem,
                const Element &element,
65
                const Scv &scv)
66
67
68
69
    {
        // update parent type sets primary variables
        ParentType::update(elemSol, problem, element, scv);

70
        updateSolidVolumeFractions(elemSol, problem, element, scv, solidState_, ParentType::numFluidComponents());
71
        // dispersivity_ = problem.spatialParams().dispersivity(element, scv, elemSol);
72
73
74
75
76

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

77
        for (int compIdx = 0; compIdx < ParentType::numFluidComponents(); ++compIdx)
78
        {
79
            moleOrMassFraction_[compIdx] = this->priVars()[compIdx];
80
            diffCoeff_[compIdx] = FluidSystem::binaryDiffusionCoefficient(compIdx, problem, element, scv);
81
82
83
84
85
86
87
        }
    }

    /*!
     * \brief Return density \f$\mathrm{[kg/m^3]}\f$ the of the fluid phase.
     *
     * We always forward to the fluid state with the phaseIdx property (see class description).
88
89
     *
     * \param pIdx TODO docme!
90
91
92
93
     */
    Scalar density(int pIdx = 0) const
    { return fluidDensity_; }

94
95
96
97
98
99
    /*!
     * \brief Returns the phase state for the control volume.
     */
    const SolidState &solidState() const
    { return solidState_; }

100
101
102
103
104
    /*!
     * \brief Return the saturation
     *
     * This method is here for compatibility reasons with other models. The saturation
     * is always 1.0 in a one-phasic context.
105
106
     *
     * \param pIdx TODO docme!
107
108
109
110
111
112
113
114
115
     */
    Scalar saturation(int pIdx = 0) const
    { return 1.0; }

    /*!
     * \brief Return the mobility
     *
     * 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
116
117
     *
     * \param pIdx TODO docme!
118
119
120
121
122
123
     */
    Scalar mobility(int pIdx = 0) const
    { return 1.0; }

    /*!
     * \brief Return molar density \f$\mathrm{[mol/m^3]}\f$ the of the fluid phase.
124
125
     *
     * \param pIdx TODO docme!
126
127
128
129
130
131
     */
    Scalar molarDensity(int pIdx = 0) const
    { return fluidDensity_/fluidMolarMass_; }

    /*!
     * \brief Return mole fraction \f$\mathrm{[mol/mol]}\f$ of a component in the phase.
132
133
     *
     * \param pIdx TODO docme!
134
135
136
     * \param compIdx The index of the component
     */
    Scalar moleFraction(int pIdx, int compIdx) const
137
    { return useMoles ? moleOrMassFraction_[compIdx] : moleOrMassFraction_[compIdx]/FluidSystem::molarMass(compIdx)*fluidMolarMass_; }
138
139
140

    /*!
     * \brief Return mass fraction \f$\mathrm{[kg/kg]}\f$ of a component in the phase.
141
142
     *
     * \param pIdx TODO docme!
143
144
145
     * \param compIdx The index of the component
     */
    Scalar massFraction(int pIdx, int compIdx) const
146
    { return useMoles ? moleOrMassFraction_[compIdx]*FluidSystem::molarMass(compIdx)/fluidMolarMass_ : moleOrMassFraction_[compIdx]; }
147
148
149

    /*!
     * \brief Return concentration \f$\mathrm{[mol/m^3]}\f$  of a component in the phase.
150
151
     *
     * \param pIdx TODO docme!
152
153
154
155
156
157
158
     * \param compIdx The index of the component
     */
    Scalar molarity(int pIdx, int compIdx) const
    { return moleFraction(pIdx, compIdx)*molarDensity(); }

    /*!
     * \brief Return the binary diffusion coefficient \f$\mathrm{[m^2/s]}\f$ in the fluid.
159
160
161
     *
     * \param pIdx TODO docme!
     * \param compIdx The index of the component
162
163
164
165
     */
    Scalar diffusionCoefficient(int pIdx, int compIdx) const
    { return diffCoeff_[compIdx]; }

166
167
168
169
170
171
    // /*!
    //  * \brief Returns the dispersivity of the fluid's streamlines.
    //  * \todo implement me
    //  */
    // const DispersivityType &dispersivity() const
    // { return dispersivity_; }
172
173
174
175
176

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

protected:
180
    SolidState solidState_;
181
    Scalar fluidDensity_, fluidMolarMass_;
182
    // DispersivityType dispersivity_;
183
184
    std::array<Scalar, ParentType::numFluidComponents()> diffCoeff_;
    std::array<Scalar, ParentType::numFluidComponents()> moleOrMassFraction_;
185
186
187
188
189
};

} // end namespace Dumux

#endif