volumevariables.hh 7.06 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
27
#ifndef DUMUX_TRACER_VOLUME_VARIABLES_HH
#define DUMUX_TRACER_VOLUME_VARIABLES_HH

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

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

33
namespace Dumux {
34
35
36
37
38
39

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

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

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

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

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

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

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

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

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

    /*!
113
     * \brief Returns the mobility.
114
115
116
     *
     * 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
117
     *
118
     * \param phaseIdx The phase index
119
     */
120
    Scalar mobility(int phaseIdx = 0) const
121
122
123
    { return 1.0; }

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

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

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

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

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

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

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

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

} // end namespace Dumux

#endif