volumevariables.hh 10.8 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 TwoPModel
22
23
24
25
26
27
 * \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

28
#include <dumux/porousmediumflow/volumevariables.hh>
29
#include "indices.hh"
30

31
32
namespace Dumux {

33
34
35
36
37
/*!
 * \ingroup TwoPModel
 * \brief Contains the quantities which are are constant within a
 *        finite volume in the two-phase model.
 */
38
39
40
template <class Traits>
class TwoPVolumeVariables
: public PorousMediumFlowVolumeVariables<Traits, TwoPVolumeVariables<Traits> >
41
{
42
43
44
45
46
47
    using ParentType = PorousMediumFlowVolumeVariables<Traits, TwoPVolumeVariables<Traits> >;

    using PermeabilityType = typename Traits::PermeabilityType;
    using ModelTraits = typename Traits::ModelTraits;
    using Indices = typename ModelTraits::Indices;
    using Scalar = typename Traits::PrimaryVariables::value_type;
48
    using FS = typename Traits::FluidSystem;
49
50
51

    enum
    {
52
53
        pressureIdx = Indices::pressureIdx,
        saturationIdx = Indices::saturationIdx,
54
55
56

        phase0Idx = FS::phase0Idx,
        phase1Idx = FS::phase1Idx
57
58
    };

59
    static constexpr auto formulation = ModelTraits::priVarFormulation();
60

61
public:
62
63
64
65
    //! export type of fluid system
    using FluidSystem = typename Traits::FluidSystem;
    //! export type of fluid state
    using FluidState = typename Traits::FluidState;
66
67

    /*!
68
69
70
71
72
73
74
75
     * \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
    */
76
77
    template<class ElemSol, class Problem, class Element, class Scv>
    void update(const ElemSol &elemSol,
78
79
                const Problem &problem,
                const Element &element,
80
                const Scv& scv)
81
    {
82
        ParentType::update(elemSol, problem, element, scv);
83

84
        completeFluidState(elemSol, problem, element, scv, fluidState_);
85

86
        using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
87
        const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
88

89
90
91
        const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
        const int nPhaseIdx = 1 - wPhaseIdx;

92
93
94
95
96
97
98
99
100
        mobility_[wPhaseIdx] =
            MaterialLaw::krw(materialParams, fluidState_.saturation(wPhaseIdx))
            / fluidState_.viscosity(wPhaseIdx);

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

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

    }

    /*!
107
108
109
110
111
112
113
114
115
     * \brief Complete the fluid state
     *
     * \param elemSol A vector containing all primary variables connected to the element
     * \param problem The problem
     * \param element The element
     * \param scv The sub control volume
     * \param fluidState The fluid state
     *
     * Set temperature, saturations, capillary pressures, viscosities, densities and enthalpies.
116
     */
117
    template<class ElemSol, class Problem, class Element, class Scv>
118
119
120
121
122
    void completeFluidState(const ElemSol& elemSol,
                            const Problem& problem,
                            const Element& element,
                            const Scv& scv,
                            FluidState& fluidState)
123
    {
124
        Scalar t =  ParentType::temperature(elemSol, problem, element, scv);
125
126
        fluidState.setTemperature(t);

127
        using MaterialLaw = typename Problem::SpatialParams::MaterialLaw;
128
        const auto& materialParams = problem.spatialParams().materialLawParams(element, scv, elemSol);
129
        const auto& priVars = ParentType::extractDofPriVars(elemSol, scv);
130

131
132
133
134
        const int wPhaseIdx = problem.spatialParams().template wettingPhase<FluidSystem>(element, scv, elemSol);
        if (formulation == TwoPFormulation::p0s1)
        {
            fluidState.setPressure(phase0Idx, priVars[pressureIdx]);
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
            if (wPhaseIdx == phase1Idx)
            {
                fluidState.setSaturation(phase1Idx, priVars[saturationIdx]);
                fluidState.setSaturation(phase0Idx, 1 - priVars[saturationIdx]);
                pc_ = MaterialLaw::pc(materialParams, fluidState.saturation(wPhaseIdx));
                fluidState.setPressure(phase1Idx, priVars[pressureIdx] - pc_);
            }
            else
            {
                const auto Sn = Traits::SaturationReconstruction::reconstructSn(problem.spatialParams(), element,
                                                                                scv, elemSol, priVars[saturationIdx]);
                fluidState.setSaturation(phase1Idx, Sn);
                fluidState.setSaturation(phase0Idx, 1 - Sn);
                pc_ = MaterialLaw::pc(materialParams, fluidState.saturation(wPhaseIdx));
                fluidState.setPressure(phase1Idx, priVars[pressureIdx] + pc_);
            }
151
        }
152
153
154
        else if (formulation == TwoPFormulation::p1s0)
        {
            fluidState.setPressure(phase1Idx, priVars[pressureIdx]);
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
            if (wPhaseIdx == phase1Idx)
            {
                const auto Sn = Traits::SaturationReconstruction::reconstructSn(problem.spatialParams(), element,
                                                                                scv, elemSol, priVars[saturationIdx]);
                fluidState.setSaturation(phase0Idx, Sn);
                fluidState.setSaturation(phase1Idx, 1 - Sn);
                pc_ = MaterialLaw::pc(materialParams, fluidState.saturation(wPhaseIdx));
                fluidState.setPressure(phase0Idx, priVars[pressureIdx] + pc_);
            }
            else
            {
                fluidState.setSaturation(phase1Idx, priVars[saturationIdx]);
                fluidState.setSaturation(phase0Idx, 1 - priVars[saturationIdx]);
                pc_ = MaterialLaw::pc(materialParams, fluidState.saturation(wPhaseIdx));
                fluidState.setPressure(phase0Idx, priVars[pressureIdx] - pc_);
            }
171
172
173
174
175
        }

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

176
        for (int phaseIdx = 0; phaseIdx < ModelTraits::numPhases(); ++phaseIdx) {
177
178
179
180
181
182
183
184
185
            // 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
186
            Scalar h = ParentType::enthalpy(fluidState, paramCache, phaseIdx);
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
215
216
217
218
219
220
221
222
223
224
225
226
227
228
            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
229
    { return pc_; }
230
231
232
233
234
235
236
237
238
239
240
241

    /*!
     * \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); }

242
243
244
245
246
247
248
249
250
    /*!
     * \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); }

251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
    /*!
     * \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_; }

266
267
268
    /*!
     * \brief Returns the permeability within the control volume in \f$[m^2]\f$.
     */
269
    const PermeabilityType& permeability() const
270
271
    { return permeability_; }

272
273
protected:
    FluidState fluidState_;
274
275

private:
276
    Scalar pc_;
277
    Scalar porosity_;
278
    PermeabilityType permeability_;
279
    Scalar mobility_[ModelTraits::numPhases()];
280
281
};

282
} // end namespace Dumux
283
284

#endif