problem.hh 13.8 KB
Newer Older
1
2
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
Bernd Flemisch's avatar
Bernd Flemisch committed
3
/*****************************************************************************
4
 *   See the file COPYING for full copying permissions.                      *
Bernd Flemisch's avatar
Bernd Flemisch committed
5
 *                                                                           *
6
 *   This program is free software: you can redistribute it and/or modify    *
Bernd Flemisch's avatar
Bernd Flemisch committed
7
 *   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
 *   (at your option) any later version.                                     *
Bernd Flemisch's avatar
Bernd Flemisch committed
10
 *                                                                           *
11
12
 *   This program is distributed in the hope that it will be useful,         *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of          *
13
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
14
15
16
17
 *   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/>.   *
Bernd Flemisch's avatar
Bernd Flemisch committed
18
 *****************************************************************************/
Melanie Darcis's avatar
doc    
Melanie Darcis committed
19
20
/*!
 * \file
Johannes Hommel's avatar
Johannes Hommel committed
21
 * \ingroup TwoPTwoCTests
22
 * \brief Problem where air is injected under a low permeable layer in a depth of 2700m.
Bernd Flemisch's avatar
Bernd Flemisch committed
23
 */
24

25
26
#ifndef DUMUX_INJECTION_PROBLEM_HH
#define DUMUX_INJECTION_PROBLEM_HH
Bernd Flemisch's avatar
Bernd Flemisch committed
27

28
29
#include <dune/grid/yaspgrid.hh>

30
#include <dumux/common/boundarytypes.hh>
31
#include <dumux/common/numeqvector.hh>
32

33
#include <dumux/discretization/cellcentered/mpfa/omethod/staticinteractionvolume.hh>
34
35
36
#include <dumux/discretization/ccmpfa.hh>
#include <dumux/discretization/cctpfa.hh>
#include <dumux/discretization/box.hh>
Timo Koch's avatar
Timo Koch committed
37

38
#include <dumux/porousmediumflow/problem.hh>
39
#include <dumux/porousmediumflow/2p2c/model.hh>
40
#include <dumux/material/fluidsystems/h2on2.hh>
Bernd Flemisch's avatar
Bernd Flemisch committed
41

42
#include "spatialparams.hh"
Bernd Flemisch's avatar
Bernd Flemisch committed
43

44
namespace Dumux {
Bernd Flemisch's avatar
Bernd Flemisch committed
45

46
47
48
49
#ifndef ENABLECACHING
#define ENABLECACHING 0
#endif

Bernd Flemisch's avatar
Bernd Flemisch committed
50
51
52
template <class TypeTag>
class InjectionProblem;

53
namespace Properties {
54
55
56
57
58
59
60
// Create new type tags
namespace TTag {
struct Injection { using InheritsFrom = std::tuple<TwoPTwoC>; };
struct InjectionBox { using InheritsFrom = std::tuple<Injection, BoxModel>; };
struct InjectionCCTpfa { using InheritsFrom = std::tuple<Injection, CCTpfaModel>; };
struct InjectionCCMpfa { using InheritsFrom = std::tuple<Injection, CCMpfaModel>; };
} // end namespace TTag
61

Bernd Flemisch's avatar
Bernd Flemisch committed
62
// Set the grid type
63
64
template<class TypeTag>
struct Grid<TypeTag, TTag::Injection> { using type = Dune::YaspGrid<2>; };
Bernd Flemisch's avatar
Bernd Flemisch committed
65
66

// Set the problem property
67
68
template<class TypeTag>
struct Problem<TypeTag, TTag::Injection> { using type = InjectionProblem<TypeTag>; };
Bernd Flemisch's avatar
Bernd Flemisch committed
69
70

// Set fluid configuration
71
72
73
74
75
76
template<class TypeTag>
struct FluidSystem<TypeTag, TTag::Injection>
{
    using type = FluidSystems::H2ON2<GetPropType<TypeTag, Properties::Scalar>,
                                     FluidSystems::H2ON2DefaultPolicy</*fastButSimplifiedRelations=*/true>>;
};
Bernd Flemisch's avatar
Bernd Flemisch committed
77

78
// Set the spatial parameters
79
80
template<class TypeTag>
struct SpatialParams<TypeTag, TTag::Injection>
81
{
82
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
83
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
84
    using type = InjectionSpatialParams<GridGeometry, Scalar>;
85
};
86

87
// Define whether mole(true) or mass (false) fractions are used
88
89
template<class TypeTag>
struct UseMoles<TypeTag, TTag::Injection> { static constexpr bool value = true; };
90
91

// Enable caching or not (reference solutions created without caching)
92
template<class TypeTag>
93
struct EnableGridGeometryCache<TypeTag, TTag::Injection> { static constexpr bool value = ENABLECACHING; };
94
95
96
97
template<class TypeTag>
struct EnableGridVolumeVariablesCache<TypeTag, TTag::Injection> { static constexpr bool value = ENABLECACHING; };
template<class TypeTag>
struct EnableGridFluxVariablesCache<TypeTag, TTag::Injection> { static constexpr bool value = ENABLECACHING; };
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115

// use the static interaction volume around interior vertices in the mpfa test
template<class TypeTag>
struct PrimaryInteractionVolume<TypeTag, TTag::InjectionCCMpfa>
{
private:
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using NodalIndexSet = GetPropType<TypeTag, Properties::DualGridNodalIndexSet>;

    // structured two-d grid
    static constexpr int numIvScvs = 4;
    static constexpr int numIvScvfs = 4;

    // use the default traits
    using Traits = CCMpfaODefaultStaticInteractionVolumeTraits< NodalIndexSet, Scalar, numIvScvs, numIvScvfs >;
public:
    using type = CCMpfaOStaticInteractionVolume< Traits >;
};
116
} // end namespace Properties
Bernd Flemisch's avatar
Bernd Flemisch committed
117
118

/*!
Johannes Hommel's avatar
Johannes Hommel committed
119
 * \ingroup TwoPTwoCTests
120
 * \brief Problem where air is injected under a low permeable layer in a depth of 2700m.
Bernd Flemisch's avatar
Bernd Flemisch committed
121
122
 *
 * The domain is sized 60m times 40m and consists of two layers, a moderately
123
124
 * permeable one (\f$ K=10e-12\f$) for \f$ y<22m\f$ and one with a lower
 * permeablility (\f$ K=10e-13\f$) in the rest of the domain.
Bernd Flemisch's avatar
Bernd Flemisch committed
125
 *
126
127
128
129
130
 * A mixture of Nitrogen and Water vapor, which is composed according to the
 * prevailing conditions (temperature, pressure) enters a water-filled aquifer.
 * This is realized with a solution-dependent Neumann boundary condition at the
 * right boundary (\f$ 5m<y<15m\f$). The aquifer is situated 2700m below sea level.
 * The injected fluid phase migrates upwards due to buoyancy.
131
 * It accumulates and partially enters the lower permeable aquitard.
132
 *
133
134
135
136
 * The model is able to use either mole or mass fractions. The property useMoles
 * can be set to either true or false in the problem file. Make sure that the
 * according units are used in the problem set-up.
 * The default setting for useMoles is true.
137
 *
138
 * This problem uses the \ref TwoPTwoCModel.
139
140
 *
 * To run the simulation execute the following line in shell:
141
142
 * <tt>./test_box2p2c</tt> or
 * <tt>./test_cc2p2c</tt>
Bernd Flemisch's avatar
Bernd Flemisch committed
143
 */
144
template <class TypeTag>
145
class InjectionProblem : public PorousMediumFlowProblem<TypeTag>
Bernd Flemisch's avatar
Bernd Flemisch committed
146
{
147
    using ParentType = PorousMediumFlowProblem<TypeTag>;
148
    using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
149
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
150
151
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
Bernd Flemisch's avatar
Bernd Flemisch committed
152

153
    using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>;
154
    using Indices = typename ModelTraits::Indices;
155

156
157
158
159
160
    // primary variable indices
    enum
    {
        pressureIdx = Indices::pressureIdx,
        switchIdx = Indices::switchIdx
Bernd Flemisch's avatar
Bernd Flemisch committed
161
162
    };

163
    // phase presence
164
    enum { wPhaseOnly = Indices::firstPhaseOnly };
Andreas Lauser's avatar
Andreas Lauser committed
165

166
167
168
    // equation indices
    enum
    {
169
170
        contiH2OEqIdx = Indices::conti0EqIdx + FluidSystem::H2OIdx,
        contiN2EqIdx = Indices::conti0EqIdx + FluidSystem::N2Idx,
Bernd Flemisch's avatar
Bernd Flemisch committed
171
172
    };

173
174
175
    // phase indices
    enum
    {
176
177
178
        gasPhaseIdx = FluidSystem::N2Idx,
        H2OIdx = FluidSystem::H2OIdx,
        N2Idx = FluidSystem::N2Idx
179
180
    };

181
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
182
    using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
183
184
185
    using GridVariables = GetPropType<TypeTag, Properties::GridVariables>;
    using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
    using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
186
    using BoundaryTypes = Dumux::BoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>;
187
    using Element = typename GridView::template Codim<0>::Entity;
188
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
189
    using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
190
    using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
Bernd Flemisch's avatar
Bernd Flemisch committed
191

192
    //! Property that defines whether mole or mass fractions are used
193
    static constexpr bool useMoles = ModelTraits::useMoles();
194

Bernd Flemisch's avatar
Bernd Flemisch committed
195
public:
196
197
    InjectionProblem(std::shared_ptr<const GridGeometry> gridGeometry)
    : ParentType(gridGeometry)
Bernd Flemisch's avatar
Bernd Flemisch committed
198
    {
199
200
201
202
203
204
205
206
207
        nTemperature_       = getParam<int>("Problem.NTemperature");
        nPressure_          = getParam<int>("Problem.NPressure");
        pressureLow_        = getParam<Scalar>("Problem.PressureLow");
        pressureHigh_       = getParam<Scalar>("Problem.PressureHigh");
        temperatureLow_     = getParam<Scalar>("Problem.TemperatureLow");
        temperatureHigh_    = getParam<Scalar>("Problem.TemperatureHigh");
        temperature_        = getParam<Scalar>("Problem.InitialTemperature");
        depthBOR_           = getParam<Scalar>("Problem.DepthBOR");
        name_               = getParam<std::string>("Problem.Name");
208

Bernd Flemisch's avatar
Bernd Flemisch committed
209
        // initialize the tables of the fluid system
210
211
212
213
214
215
        FluidSystem::init(/*Tmin=*/temperatureLow_,
                          /*Tmax=*/temperatureHigh_,
                          /*nT=*/nTemperature_,
                          /*pmin=*/pressureLow_,
                          /*pmax=*/pressureHigh_,
                          /*np=*/nPressure_);
216

217
        // stating in the console whether mole or mass fractions are used
218
219
        if(useMoles)
            std::cout<<"problem uses mole-fractions"<<std::endl;
220
        else
221
            std::cout<<"problem uses mass-fractions"<<std::endl;
Bernd Flemisch's avatar
Bernd Flemisch committed
222
223
224
225
226
227
228
229
    }

    /*!
     * \name Problem parameters
     */
    // \{

    /*!
230
     * \brief Returns the problem name
Bernd Flemisch's avatar
Bernd Flemisch committed
231
232
233
     *
     * This is used as a prefix for files generated by the simulation.
     */
234
    const std::string& name() const
235
    { return name_; }
Bernd Flemisch's avatar
Bernd Flemisch committed
236
237

    /*!
238
     * \brief Returns the temperature [K]
Bernd Flemisch's avatar
Bernd Flemisch committed
239
     */
240
    Scalar temperature() const
241
    { return temperature_; }
Bernd Flemisch's avatar
Bernd Flemisch committed
242
243
244
245
246
247
248
249
250
251

    // \}

    /*!
     * \name Boundary conditions
     */
    // \{

    /*!
     * \brief Specifies which kind of boundary condition should be
252
     *        used for which equation on a given boundary segment
Melanie Darcis's avatar
doc    
Melanie Darcis committed
253
     *
254
     * \param globalPos The global position
Bernd Flemisch's avatar
Bernd Flemisch committed
255
     */
256
    BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
257
    {
258
        BoundaryTypes bcTypes;
259
        if (globalPos[0] < eps_)
260
            bcTypes.setAllDirichlet();
Bernd Flemisch's avatar
Bernd Flemisch committed
261
        else
262
263
            bcTypes.setAllNeumann();
        return bcTypes;
Bernd Flemisch's avatar
Bernd Flemisch committed
264
265
266
    }

    /*!
267
     * \brief Evaluates the boundary conditions for a Dirichlet boundary segment
Melanie Darcis's avatar
doc    
Melanie Darcis committed
268
     *
269
     * \param globalPos The global position
Bernd Flemisch's avatar
Bernd Flemisch committed
270
     */
271
    PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
272
    {
273
        return initial_(globalPos);
Bernd Flemisch's avatar
Bernd Flemisch committed
274
275
276
    }

    /*!
277
     * \brief Evaluates the boundary conditions for a Neumann
278
279
     *        boundary segment in dependency on the current solution.
     *
Melanie Darcis's avatar
doc    
Melanie Darcis committed
280
     * \param element The finite element
281
     * \param fvGeometry The finite volume geometry of the element
282
     * \param elemVolVars All volume variables for the element
283
     * \param elemFluxVarsCache Flux variables caches for all faces in stencil
284
     * \param scvf The sub-control volume face
Melanie Darcis's avatar
doc    
Melanie Darcis committed
285
     *
286
287
288
289
     * This method is used for cases, when the Neumann condition depends on the
     * solution and requires some quantities that are specific to the fully-implicit method.
     * The \a values store the mass flux of each phase normal to the boundary.
     * Negative values indicate an inflow.
Bernd Flemisch's avatar
Bernd Flemisch committed
290
     */
291
    NumEqVector neumann(const Element& element,
292
293
                        const FVElementGeometry& fvGeometry,
                        const ElementVolumeVariables& elemVolVars,
294
                        const ElementFluxVariablesCache& elemFluxVarsCache,
295
                        const SubControlVolumeFace& scvf) const
Bernd Flemisch's avatar
Bernd Flemisch committed
296
    {
297
        NumEqVector values(0.0);
298
299
300
301

        const auto& globalPos = scvf.ipGlobal();

        Scalar injectedPhaseMass = 1e-3;
302
        Scalar moleFracW = elemVolVars[scvf.insideScvIdx()].moleFraction(gasPhaseIdx, H2OIdx);
Timo Koch's avatar
Timo Koch committed
303
        if (globalPos[1] < 14 - eps_ && globalPos[1] > 6.5 - eps_)
304
        {
305
306
            values[contiN2EqIdx] = -(1-moleFracW)*injectedPhaseMass/FluidSystem::molarMass(N2Idx); //mole/(m^2*s) -> kg/(s*m^2)
            values[contiH2OEqIdx] = -moleFracW*injectedPhaseMass/FluidSystem::molarMass(H2OIdx); //mole/(m^2*s) -> kg/(s*m^2)
307
308
        }
        return values;
Bernd Flemisch's avatar
Bernd Flemisch committed
309
310
311
312
313
314
315
316
317
318
    }

    // \}

    /*!
     * \name Volume terms
     */
    // \{

    /*!
319
     * \brief Evaluates the initial values for a control volume.
Melanie Darcis's avatar
doc    
Melanie Darcis committed
320
     *
321
     * \param globalPos The global position
Bernd Flemisch's avatar
Bernd Flemisch committed
322
     */
323
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
324
    { return initial_(globalPos); }
Bernd Flemisch's avatar
Bernd Flemisch committed
325
326
327
328

    // \}

private:
329
    /*!
330
     * \brief Evaluates the initial values for a control volume.
331
332
333
334
335
     *
     * The internal method for the initial condition
     *
     * \param globalPos The global position
     */
336
    PrimaryVariables initial_(const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
337
    {
338
        PrimaryVariables priVars(0.0);
339
        priVars.setState(wPhaseOnly);
340

341
342
        Scalar densityW = FluidSystem::H2O::liquidDensity(temperature_, 1e5);

Timo Koch's avatar
Timo Koch committed
343
        Scalar pl = 1e5 - densityW*this->spatialParams().gravity(globalPos)[1]*(depthBOR_ - globalPos[1]);
344
345
        Scalar moleFracLiquidN2 = pl*0.95/BinaryCoeff::H2O_N2::henry(temperature_);
        Scalar moleFracLiquidH2O = 1.0 - moleFracLiquidN2;
346

347
        Scalar meanM =
348
349
            FluidSystem::molarMass(H2OIdx)*moleFracLiquidH2O +
            FluidSystem::molarMass(N2Idx)*moleFracLiquidN2;
350
        if(useMoles)
351
        {
352
            //mole-fraction formulation
353
            priVars[switchIdx] = moleFracLiquidN2;
354
        }
355
        else
356
        {
357
            //mass fraction formulation
358
            Scalar massFracLiquidN2 = moleFracLiquidN2*FluidSystem::molarMass(N2Idx)/meanM;
359
            priVars[switchIdx] = massFracLiquidN2;
360
        }
361
        priVars[pressureIdx] = pl;
362
        return priVars;
Bernd Flemisch's avatar
Bernd Flemisch committed
363
364
    }

365
366
    Scalar temperature_;
    Scalar depthBOR_;
367
    static constexpr Scalar eps_ = 1e-6;
368
369
370
371
372

    int nTemperature_;
    int nPressure_;
    Scalar pressureLow_, pressureHigh_;
    Scalar temperatureLow_, temperatureHigh_;
Timo Koch's avatar
Timo Koch committed
373
374

    std::string name_;
Bernd Flemisch's avatar
Bernd Flemisch committed
375
};
376

Timo Koch's avatar
Timo Koch committed
377
} // end namespace Dumux
Bernd Flemisch's avatar
Bernd Flemisch committed
378
379

#endif