volumevariables.hh 7.34 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
// -*- 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    *
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (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

Timo Koch's avatar
Timo Koch committed
27
28
#include <dune/common/fvector.hh>

Sina Ackermann's avatar
Sina Ackermann committed
29
#include <dumux/common/properties.hh>
30
#include <dumux/porousmediumflow/volumevariables.hh>
31
32
33
34
35
36
37
38
39
40

namespace Dumux
{

/*!
 * \ingroup TracerModel
 * \brief Contains the quantities which are constant within a
 *        finite volume for the tracer model.
 */
template <class TypeTag>
41
class TracerVolumeVariables : public PorousMediumFlowVolumeVariables<TypeTag>
42
{
43
    using ParentType = PorousMediumFlowVolumeVariables<TypeTag>;
44
45
46
47
48

    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
    using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using Implementation = typename GET_PROP_TYPE(TypeTag, VolumeVariables);
49
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::LocalView;
50
    using SubControlVolume = typename FVElementGeometry::SubControlVolume;
51
52
53
54
    using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
    using ElementSolutionVector = typename GET_PROP_TYPE(TypeTag, ElementSolutionVector);
    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);

55
56
57
    static constexpr bool useMoles = GET_PROP_VALUE(TypeTag, UseMoles);
    static constexpr int dimWorld = GridView::dimensionworld;
    static constexpr int numComponents = GET_PROP_VALUE(TypeTag, NumComponents);
58
59
60
61
62
63
64

    using GlobalPosition = Dune::FieldVector<Scalar,dimWorld>;
    using Element = typename GridView::template Codim<0>::Entity;

public:

    /*!
65
66
67
68
69
70
71
     * \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
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
     */
    void update(const ElementSolutionVector &elemSol,
                const Problem &problem,
                const Element &element,
                const SubControlVolume &scv)
    {
        // update parent type sets primary variables
        ParentType::update(elemSol, problem, element, scv);

        porosity_ = problem.spatialParams().porosity(element, scv, elemSol);
        dispersivity_ = problem.spatialParams().dispersivity(element, scv, elemSol);

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

        for (int compIdx = 0; compIdx < numComponents; ++compIdx)
        {
90
            moleOrMassFraction_[compIdx] = this->priVars()[compIdx];
91
            diffCoeff_[compIdx] = FluidSystem::binaryDiffusionCoefficient(compIdx, problem, element, scv);
92
93
94
95
96
97
98
        }
    }

    /*!
     * \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).
99
100
     *
     * \param pIdx TODO docme!
101
102
103
104
105
106
107
108
109
     */
    Scalar density(int pIdx = 0) const
    { return fluidDensity_; }

    /*!
     * \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.
110
111
     *
     * \param pIdx TODO docme!
112
113
114
115
116
117
118
119
120
     */
    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
121
122
     *
     * \param pIdx TODO docme!
123
124
125
126
127
128
     */
    Scalar mobility(int pIdx = 0) const
    { return 1.0; }

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

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

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

    /*!
     * \brief Return concentration \f$\mathrm{[mol/m^3]}\f$  of a component in the phase.
155
156
     *
     * \param pIdx TODO docme!
157
158
159
160
161
162
163
     * \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.
164
165
166
     *
     * \param pIdx TODO docme!
     * \param compIdx The index of the component
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
     */
    Scalar diffusionCoefficient(int pIdx, int compIdx) const
    { return diffCoeff_[compIdx]; }

    /*!
     * \brief Returns the dispersivity of the fluid's streamlines.
     */
    const GlobalPosition &dispersivity() const
    { return dispersivity_; }

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

protected:
184
    Scalar porosity_;    // Effective porosity within the control volume
185
186
187
    Scalar fluidDensity_, fluidMolarMass_;
    GlobalPosition dispersivity_;
    std::array<Scalar, numComponents> diffCoeff_;
188
    std::array<Scalar, numComponents> moleOrMassFraction_;
189
190
191
192
193
};

} // end namespace Dumux

#endif