injectionproblem2pni.hh 13 KB
Newer Older
1
// $Id$
Bernd Flemisch's avatar
Bernd Flemisch committed
2
3
4
5
6
7
8
/*****************************************************************************
 *   Copyright (C) 2009 by Melanie Darcis                                    *
 *   Copyright (C) 2009 by Andreas Lauser                                    *
 *   Institute of Hydraulic Engineering                                      *
 *   University of Stuttgart, Germany                                        *
 *   email: <givenname>.<name>@iws.uni-stuttgart.de                          *
 *                                                                           *
9
 *   This program is free software: you can redistribute it and/or modify    *
Bernd Flemisch's avatar
Bernd Flemisch committed
10
 *   it under the terms of the GNU General Public License as published by    *
11
12
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (at your option) any later version.                                     *
Bernd Flemisch's avatar
Bernd Flemisch committed
13
 *                                                                           *
14
15
16
17
18
19
20
 *   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/>.   *
Bernd Flemisch's avatar
Bernd Flemisch committed
21
 *****************************************************************************/
Melanie Darcis's avatar
doc    
Melanie Darcis committed
22
23
24
25
26
27
28
29
/*!
 * \file
 *
 * \brief Nonisothermal gas injection problem where a gas (e.g. air) is injected into a fully
 *        water saturated medium. During buoyancy driven upward migration the gas
 *        passes a high temperature area.
 */

Bernd Flemisch's avatar
Bernd Flemisch committed
30
31
32
33
34
35
36
#ifndef DUMUX_INJECTIONPROBLEM2PNI_HH
#define DUMUX_INJECTIONPROBLEM2PNI_HH

#include <dune/grid/io/file/dgfparser/dgfug.hh>
#include <dune/grid/io/file/dgfparser/dgfs.hh>
#include <dune/grid/io/file/dgfparser/dgfyasp.hh>

37
#include <dumux/boxmodels/2pni/2pnimodel.hh>
Bernd Flemisch's avatar
Bernd Flemisch committed
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70

#include <dumux/material/fluidsystems/h2o_n2_system.hh>

// use the same spatial parameters as the injection problem of the
// 2p2c test program
#include "../2p2c/injectionspatialparameters.hh"

#define ISOTHERMAL 0

namespace Dumux {

template <class TypeTag>
class InjectionProblem2PNI;

namespace Properties
{
#if !ISOTHERMAL
NEW_TYPE_TAG(InjectionProblem2PNI, INHERITS_FROM(BoxTwoPNI));
#else
NEW_TYPE_TAG(InjectionProblem2PNI, INHERITS_FROM(BoxTwoP));
#endif

// Set the grid type
SET_PROP(InjectionProblem2PNI, Grid)
{
#if HAVE_UG
    typedef Dune::UGGrid<2> type;
#else
    typedef Dune::SGrid<2, 2> type;
    //typedef Dune::YaspGrid<2> type;
#endif
};

71
#if HAVE_DUNE_PDELAB
Bernd Flemisch's avatar
Bernd Flemisch committed
72
73
74
75
76
77
78
SET_PROP(InjectionProblem2PNI, LocalFEMSpace)
{
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
    enum{dim = GridView::dimension};

public:
79
80
    typedef Dune::PDELab::Q1LocalFiniteElementMap<Scalar,Scalar,dim> type; // for cubes
//    typedef Dune::PDELab::P1LocalFiniteElementMap<Scalar,Scalar,dim> type; // for simplices
Bernd Flemisch's avatar
Bernd Flemisch committed
81
};
82
#endif // HAVE_DUNE_PDELAB
Bernd Flemisch's avatar
Bernd Flemisch committed
83
84
85
86
87
88
89

// Set the problem property
SET_PROP(InjectionProblem2PNI, Problem)
{
    typedef Dumux::InjectionProblem2PNI<TypeTag> type;
};

90
// Set the spatial parameters. we use the same spatial parameters as the
Bernd Flemisch's avatar
Bernd Flemisch committed
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
// 2p2c injection problem
SET_PROP(InjectionProblem2PNI, SpatialParameters)
{ typedef InjectionSpatialParameters<TypeTag> type; };

#if 1
// Use the same fluid system as the 2p2c injection problem
SET_PROP(InjectionProblem2PNI, FluidSystem)
{
    typedef H2O_N2_System<TypeTag> type;
};
#else
// Set the wetting phase
SET_PROP(InjectionProblem2PNI, WettingPhase)
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
public:
    typedef Dumux::LiquidPhase<Scalar, Dumux::SimpleH2O<Scalar> > type;
};

// Set the non-wetting phase
SET_PROP(InjectionProblem2PNI, NonwettingPhase)
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
public:
    typedef Dumux::GasPhase<Scalar, Dumux::N2<Scalar> > type;
};
#endif

// Enable gravity
SET_BOOL_PROP(InjectionProblem2PNI, EnableGravity, true);

// write convergence behaviour to disk?
125
SET_BOOL_PROP(InjectionProblem2PNI, NewtonWriteConvergence, false);
Bernd Flemisch's avatar
Bernd Flemisch committed
126
127
128
}

/*!
129
 * \ingroup TwoPNIBoxModel
Bernd Flemisch's avatar
Bernd Flemisch committed
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
 * \brief Nonisothermal gas injection problem where a gas (e.g. air) is injected into a fully
 *        water saturated medium. During buoyancy driven upward migration the gas
 *        passes a high temperature area.
 *
 * The domain is sized 40 m times 40 m. The rectangular area with the increased temperature (380 K)
 * starts at (20 m, 5 m) and ends at (30 m, 35 m)
 *
 * For the mass conservation equation neumann boundary conditions are used on
 * the top and on the bottom of the domain, while dirichlet conditions
 * apply on the left and the right boundary.
 * For the energy conservation equation dirichlet boundary conditions are applied
 * on all boundaries.
 *
 * Gas is injected at the bottom boundary from 15 m to 25 m at a rate of
 * 0.001 kg/(s m), the remaining neumann boundaries are no-flow
 * boundaries.
 *
 * At the dirichlet boundaries a hydrostatic pressure, a gas saturation of zero and
 * a geothermal temperature gradient of 0.03 K/m are applied.
 *
 * This problem uses the \ref TwoPNIBoxModel.
 *
 * This problem should typically be simulated for 300000 seconds.
 * A good choice for the initial time step size is 1000 seconds.
 *
 * To run the simulation execute the following line in shell:
 * <tt>./test_2pni ./grids/test_2pni.dgf 300000 1000</tt>
 */
template<class TypeTag>
class InjectionProblem2PNI
#if !ISOTHERMAL
161
  : public TwoPNIProblem<TypeTag>
Bernd Flemisch's avatar
Bernd Flemisch committed
162
#else
163
  : public TwoPNIProblem<TypeTag>
Bernd Flemisch's avatar
Bernd Flemisch committed
164
165
#endif
{
166
167
168
169
    typedef InjectionProblem2PNI<TypeTag> ThisType;
    typedef TwoPNIProblem<TypeTag> ParentType;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
Bernd Flemisch's avatar
Bernd Flemisch committed
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193

#if ISOTHERMAL
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(TwoPIndices)) Indices;
#else
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(TwoPNIIndices)) Indices;
#endif
    enum {
        numEq = GET_PROP_VALUE(TypeTag, PTAG(NumEq)),
        pressureIdx = Indices::pressureIdx,
        saturationIdx = Indices::saturationIdx,

        contiWEqIdx = Indices::contiWEqIdx,
        contiNEqIdx = Indices::contiNEqIdx,

#if !ISOTHERMAL
        temperatureIdx = Indices::temperatureIdx,
        energyEqIdx = Indices::energyEqIdx,
#endif

        // Grid and world dimension
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld,
    };

194

195
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(PrimaryVariables)) PrimaryVariables;
196
197
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(BoundaryTypes)) BoundaryTypes;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(TimeManager)) TimeManager;
Bernd Flemisch's avatar
Bernd Flemisch committed
198
199
200
201
202
203
204
205
206
207
208
209

    typedef typename GridView::template Codim<0>::Entity Element;
    typedef typename GridView::template Codim<dim>::Entity Vertex;
    typedef typename GridView::Intersection Intersection;

    typedef typename GET_PROP_TYPE(TypeTag, PTAG(FVElementGeometry)) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(FluidSystem)) FluidSystem;

    typedef Dune::FieldVector<Scalar, dim> LocalPosition;
    typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition;

public:
Melanie Darcis's avatar
doc    
Melanie Darcis committed
210
211
212
213
214
215
    /*!
     * \brief The constructor
     *
     * \param timeManager The time manager
     * \param gridView The grid view
     */
216
217
    InjectionProblem2PNI(TimeManager &timeManager, const GridView &gridView)
        : ParentType(timeManager, gridView)
Bernd Flemisch's avatar
Bernd Flemisch committed
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
    {
        // initialize the tables of the fluid system
        FluidSystem::init();
    }

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

    /*!
     * \brief The problem name.
     *
     * This is used as a prefix for files generated by the simulation.
     */
    const char *name() const
    { return "injection2pni"; }

Melanie Darcis's avatar
doc    
Melanie Darcis committed
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
#if ISOTHERMAL
    /*!
     * \brief Returns the temperature within the domain.
     *
     * \param element The element
     * \param fvElemGeom The finite-volume geometry in the box scheme
     * \param scvIdx The local vertex index (SCV index)
     *
     * This problem assumes a temperature of 10 degrees Celsius.
     */
    Scalar temperature(const Element &element,
                       const FVElementGeometry &fvElemGeom,
                       int scvIdx) const
    {
        return 273.15 + 30; // [K]
    };
#endif

Bernd Flemisch's avatar
Bernd Flemisch committed
254
255
256
257
258
259
260
261
262
263
    // \}

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

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
Melanie Darcis's avatar
doc    
Melanie Darcis committed
264
265
     *
     * \param values The boundary types for the conservation equations
266
     * \param vertex The vertex for which the boundary type is set
Bernd Flemisch's avatar
Bernd Flemisch committed
267
     */
268
    void boundaryTypes(BoundaryTypes &values, const Vertex &vertex) const
Bernd Flemisch's avatar
Bernd Flemisch committed
269
    {
270
        const GlobalPosition globalPos = vertex.geometry().center();
Bernd Flemisch's avatar
Bernd Flemisch committed
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287

        if (globalPos[0] < eps_)
            values.setAllDirichlet();
        else
            values.setAllNeumann();

#if !ISOTHERMAL
        // set a dirichlet value for the temperature, use the energy
        // equation to set the value
        values.setDirichlet(temperatureIdx, energyEqIdx);
#endif
    }

    /*!
     * \brief Evaluate the boundary conditions for a dirichlet
     *        boundary segment.
     *
Melanie Darcis's avatar
doc    
Melanie Darcis committed
288
     * \param values The dirichlet values for the primary variables
289
     * \param vertex The vertex for which the boundary type is set
Melanie Darcis's avatar
doc    
Melanie Darcis committed
290
     *
Bernd Flemisch's avatar
Bernd Flemisch committed
291
292
     * For this method, the \a values parameter stores primary variables.
     */
293
    void dirichlet(PrimaryVariables &values, const Vertex &vertex) const
Bernd Flemisch's avatar
Bernd Flemisch committed
294
    {
295
        const GlobalPosition globalPos = vertex.geometry().center();
Bernd Flemisch's avatar
Bernd Flemisch committed
296
297
298
299
300
301
302
303
304
305
306
307
308

        Scalar densityW = 1000.0;
        values[pressureIdx] = 1e5 + (depthBOR_ - globalPos[1])*densityW*9.81;
        values[saturationIdx] = 0.0;
#if !ISOTHERMAL
        values[temperatureIdx] = 283.0 + (depthBOR_ - globalPos[1])*0.03;
#endif
    }

    /*!
     * \brief Evaluate the boundary conditions for a neumann
     *        boundary segment.
     *
Melanie Darcis's avatar
doc    
Melanie Darcis committed
309
310
311
312
313
314
315
     * \param values The neumann values for the conservation equations
     * \param element The finite element
     * \param fvElemGeom The finite-volume geometry in the box scheme
     * \param is The intersection between element and boundary
     * \param scvIdx The local vertex index
     * \param boundaryFaceIdx The index of the boundary face
     *
Bernd Flemisch's avatar
Bernd Flemisch committed
316
317
318
     * For this method, the \a values parameter stores the mass flux
     * in normal direction of each phase. Negative values mean influx.
     */
319
    void neumann(PrimaryVariables &values,
Bernd Flemisch's avatar
Bernd Flemisch committed
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
                 const Element &element,
                 const FVElementGeometry &fvElemGeom,
                 const Intersection &is,
                 int scvIdx,
                 int boundaryFaceIdx) const
    {
        const GlobalPosition &globalPos = element.geometry().corner(scvIdx);

        values = 0;
        if (globalPos[1] < 15 && globalPos[1] > 5) {
            // inject air. negative values mean injection
            values[contiNEqIdx] = -1e-3; // kg/(s*m^2)
        }
    }

    // \}

Melanie Darcis's avatar
doc    
Melanie Darcis committed
337

Bernd Flemisch's avatar
Bernd Flemisch committed
338
339
340
341
342
343
344
345
346
    /*!
     * \name Volume terms
     */
    // \{

    /*!
     * \brief Evaluate the source term for all phases within a given
     *        sub-control-volume.
     *
Melanie Darcis's avatar
doc    
Melanie Darcis committed
347
348
349
     * \param values The source and sink values for the conservation equations
     * \param element The finite element
     * \param fvElemGeom The finite-volume geometry in the box scheme
Bernd Flemisch's avatar
Bernd Flemisch committed
350
     * \param scvIdx The local vertex index
Melanie Darcis's avatar
doc    
Melanie Darcis committed
351
     *
Bernd Flemisch's avatar
Bernd Flemisch committed
352
353
354
355
     * For this method, the \a values parameter stores the rate mass
     * generated or annihilate per volume unit. Positive values mean
     * that mass is created, negative ones mean that it vanishes.
     */
356
    void source(PrimaryVariables &values,
Bernd Flemisch's avatar
Bernd Flemisch committed
357
                const Element &element,
Bernd Flemisch's avatar
Bernd Flemisch committed
358
359
                const FVElementGeometry &fvElemGeom,
                int scvIdx) const
Bernd Flemisch's avatar
Bernd Flemisch committed
360
361
362
363
364
365
366
    {
           values = Scalar(0.0);
    }

    /*!
     * \brief Evaluate the initial value for a control volume.
     *
Melanie Darcis's avatar
doc    
Melanie Darcis committed
367
368
369
370
371
     * \param values The initial values for the primary variables
     * \param element The finite element
     * \param fvElemGeom The finite-volume geometry in the box scheme
     * \param scvIdx The local vertex index
     *
Bernd Flemisch's avatar
Bernd Flemisch committed
372
373
374
     * For this method, the \a values parameter stores primary
     * variables.
     */
375
    void initial(PrimaryVariables &values,
Bernd Flemisch's avatar
Bernd Flemisch committed
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
                 const Element &element,
                 const FVElementGeometry &fvElemGeom,
                 int scvIdx) const
    {
        const GlobalPosition &globalPos = element.geometry().corner(scvIdx);

        Scalar densityW = 1000.0;
        values[pressureIdx] = 1e5 + (depthBOR_ - globalPos[1])*densityW*9.81;
        values[saturationIdx] = 0.0;

#if !ISOTHERMAL
        values[temperatureIdx] = 283.0 + (depthBOR_ - globalPos[1])*0.03;
        if (globalPos[0] > 20 && globalPos[0] < 30 && globalPos[1] > 5 && globalPos[1] < 35)
            values[temperatureIdx] = 380;
#endif // !ISOTHERMAL
    }
    // \}

private:
    static const Scalar depthBOR_ = 2700.0; // [m]
    static const Scalar eps_ = 1e-6;
};
} //end namespace

#endif