volumevariables.hh 10.3 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
// -*- 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
 *
 * \brief Contains the quantities which are constant within a
 *        finite volume in the two-phase model.
 */
#ifndef DUMUX_2P_VOLUME_VARIABLES_HH
#define DUMUX_2P_VOLUME_VARIABLES_HH

#include "properties.hh"

30
#include <dumux/discretization/volumevariables.hh>
31
32
33
34
35
36
37
38
39
40
41
42
43
44

#include <dune/common/fvector.hh>

namespace Dumux
{
/*!
 * \ingroup TwoPModel
 * \ingroup ImplicitVolumeVariables
 * \brief Contains the quantities which are are constant within a
 *        finite volume in the two-phase model.
 */
template <class TypeTag>
class TwoPVolumeVariables : public ImplicitVolumeVariables<TypeTag>
{
45
46
47
48
49
    using ParentType = ImplicitVolumeVariables<TypeTag>;

    using Implementation = typename GET_PROP_TYPE(TypeTag, VolumeVariables);
    using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
    using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
50
51
    using SpatialParams = typename GET_PROP_TYPE(TypeTag, SpatialParams);
    using PermeabilityType = typename SpatialParams::PermeabilityType;
52
53
54
55
56
57
58
59
60
    using Indices = typename GET_PROP_TYPE(TypeTag, Indices);
    using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
    using MaterialLaw = typename GET_PROP_TYPE(TypeTag, MaterialLaw);
    using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
    using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
    using ElementSolutionVector = typename GET_PROP_TYPE(TypeTag, ElementSolutionVector);

    enum
    {
61
62
63
64
65
66
67
68
69
70
        pwsn = Indices::pwsn,
        pnsw = Indices::pnsw,
        pressureIdx = Indices::pressureIdx,
        saturationIdx = Indices::saturationIdx,
        wPhaseIdx = Indices::wPhaseIdx,
        nPhaseIdx = Indices::nPhaseIdx,
        numPhases = GET_PROP_VALUE(TypeTag, NumPhases),
        formulation = GET_PROP_VALUE(TypeTag, Formulation)
    };

71
72
73
74
    using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
    using Element = typename GridView::template Codim<0>::Entity;

    static constexpr bool isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox);
75
76
77

public:
    // export type of fluid state for non-isothermal models
78
    using FluidState = typename GET_PROP_TYPE(TypeTag, FluidState);
79
80
81
82

    /*!
     * \copydoc ImplicitVolumeVariables::update
     */
83
    void update(const ElementSolutionVector &elemSol,
84
85
                const Problem &problem,
                const Element &element,
86
                const SubControlVolume& scv)
87
    {
88
        ParentType::update(elemSol, problem, element, scv);
89

90
        completeFluidState(elemSol, problem, element, scv, fluidState_);
91

92
        const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
93
94
95
96
97
98
99
100
101
102

        mobility_[wPhaseIdx] =
            MaterialLaw::krw(materialParams, fluidState_.saturation(wPhaseIdx))
            / fluidState_.viscosity(wPhaseIdx);

        mobility_[nPhaseIdx] =
            MaterialLaw::krn(materialParams, fluidState_.saturation(wPhaseIdx))
            / fluidState_.viscosity(nPhaseIdx);

        // porosity
103
104
        porosity_ = problem.spatialParams().porosity(element, scv, elemSol);
        permeability_ = problem.spatialParams().permeability(element, scv, elemSol);
105
106

        // energy related quantities not belonging to the fluid state
107
        asImp_().updateEnergy_(elemSol, problem, element, scv);
108
109
110
111
112
    }

    /*!
     * \copydoc ImplicitModel::completeFluidState
     */
113
    static void completeFluidState(const ElementSolutionVector& elemSol,
114
115
                                   const Problem& problem,
                                   const Element& element,
116
                                   const SubControlVolume& scv,
117
118
                                   FluidState& fluidState)
    {
119
        Scalar t = Implementation::temperature_(elemSol, problem, element, scv);
120
121
        fluidState.setTemperature(t);

122
        const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
123
        const auto& priVars = ParentType::extractDofPriVars(elemSol, scv);
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214

        if (int(formulation) == pwsn) {
            Scalar sn = priVars[saturationIdx];
            fluidState.setSaturation(nPhaseIdx, sn);
            fluidState.setSaturation(wPhaseIdx, 1 - sn);

            Scalar pw = priVars[pressureIdx];
            fluidState.setPressure(wPhaseIdx, pw);
            fluidState.setPressure(nPhaseIdx,
                                   pw + MaterialLaw::pc(materialParams, 1 - sn));
        }
        else if (int(formulation) == pnsw) {
            Scalar sw = priVars[saturationIdx];
            fluidState.setSaturation(wPhaseIdx, sw);
            fluidState.setSaturation(nPhaseIdx, 1 - sw);

            Scalar pn = priVars[pressureIdx];
            fluidState.setPressure(nPhaseIdx, pn);
            fluidState.setPressure(wPhaseIdx,
                                   pn - MaterialLaw::pc(materialParams, sw));
        }

        typename FluidSystem::ParameterCache paramCache;
        paramCache.updateAll(fluidState);

        for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
            // compute and set the viscosity
            Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
            fluidState.setViscosity(phaseIdx, mu);

            // compute and set the density
            Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
            fluidState.setDensity(phaseIdx, rho);

            // compute and set the enthalpy
            Scalar h = Implementation::enthalpy_(fluidState, paramCache, phaseIdx);
            fluidState.setEnthalpy(phaseIdx, h);
        }
    }

    /*!
     * \brief Returns the phase state for the control volume.
     */
    const FluidState &fluidState() const
    { return fluidState_; }

    /*!
     * \brief Returns the saturation of a given phase within
     *        the control volume in \f$[-]\f$.
     *
     * \param phaseIdx The phase index
     */
    Scalar saturation(int phaseIdx) const
    { return fluidState_.saturation(phaseIdx); }

    /*!
     * \brief Returns the mass density of a given phase within the
     *        control volume in \f$[kg/m^3]\f$.
     *
     * \param phaseIdx The phase index
     */
    Scalar density(int phaseIdx) const
    { return fluidState_.density(phaseIdx); }

    /*!
     * \brief Returns the effective pressure of a given phase within
     *        the control volume in \f$[kg/(m*s^2)=N/m^2=Pa]\f$.
     *
     * \param phaseIdx The phase index
     */
    Scalar pressure(int phaseIdx) const
    { return fluidState_.pressure(phaseIdx); }

    /*!
     * \brief Returns the capillary pressure within the control volume
     * in \f$[kg/(m*s^2)=N/m^2=Pa]\f$.
     */
    Scalar capillaryPressure() const
    { return fluidState_.pressure(nPhaseIdx) - fluidState_.pressure(wPhaseIdx); }

    /*!
     * \brief Returns temperature inside the sub-control volume
     * in \f$[K]\f$.
     *
     * Note that we assume thermodynamic equilibrium, i.e. the
     * temperature of the rock matrix and of all fluid phases are
     * identical.
     */
    Scalar temperature() const
    { return fluidState_.temperature(/*phaseIdx=*/0); }

215
216
217
218
219
220
221
222
223
    /*!
     * \brief Returns the dynamic viscosity of the fluid within the
     *        control volume in \f$\mathrm{[Pa s]}\f$.
     *
     * \param phaseIdx The phase index
     */
    Scalar viscosity(int phaseIdx) const
    { return fluidState_.viscosity(phaseIdx); }

224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
    /*!
     * \brief Returns the effective mobility of a given phase within
     *        the control volume in \f$[s*m/kg]\f$.
     *
     * \param phaseIdx The phase index
     */
    Scalar mobility(int phaseIdx) const
    { return mobility_[phaseIdx]; }

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

239
240
241
242
243
244
    /*!
     * \brief Returns the permeability within the control volume in \f$[m^2]\f$.
     */
    PermeabilityType permeability() const
    { return permeability_; }

245
protected:
246
    static Scalar temperature_(const ElementSolutionVector &elemSol,
247
248
                               const Problem& problem,
                               const Element &element,
249
                               const SubControlVolume &scv)
250
    {
251
        return problem.temperatureAtPos(scv.dofPosition());
252
253
254
255
256
    }

    template<class ParameterCache>
    static Scalar enthalpy_(const FluidState& fluidState,
                            const ParameterCache& paramCache,
257
                            const int phaseIdx)
258
259
260
261
262
    {
        return 0;
    }

    /*!
263
     * \brief Called by update() to compute the energy related quantities.
264
     */
265
    void updateEnergy_(const ElementSolutionVector &elemSol,
266
267
                       const Problem &problem,
                       const Element &element,
268
                       const SubControlVolume& scv)
269
    {}
270
271
272

    FluidState fluidState_;
    Scalar porosity_;
273
    PermeabilityType permeability_;
274
275
276
277
278
279
280
281
282
283
284
285
286
    Scalar mobility_[numPhases];

private:
    Implementation &asImp_()
    { return *static_cast<Implementation*>(this); }

    const Implementation &asImp_() const
    { return *static_cast<const Implementation*>(this); }
};

}

#endif