kuevetteproblem.hh 12.9 KB
Newer Older
Bernd Flemisch's avatar
Bernd Flemisch committed
1
2
3
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
4
 *   See the file COPYING for full copying permissions.                      *
Bernd Flemisch's avatar
Bernd Flemisch committed
5
6
7
8
9
10
11
12
 *                                                                           *
 *   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          *
13
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
Bernd Flemisch's avatar
Bernd Flemisch committed
14
15
16
17
18
19
20
21
 *   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
 *
22
23
 * \brief Non-isothermal gas injection problem where a gas (e.g. steam/air)
 *        is injected into a unsaturated porous medium with a residually
Bernd Flemisch's avatar
Bernd Flemisch committed
24
25
26
27
28
 *        trapped NAPL contamination.
 */
#ifndef DUMUX_KUEVETTE3P3CNIPROBLEM_HH
#define DUMUX_KUEVETTE3P3CNIPROBLEM_HH

29
#include <dune/common/version.hh>
Bernd Flemisch's avatar
Bernd Flemisch committed
30
31
32
33
#include <dune/grid/io/file/dgfparser/dgfyasp.hh>

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

34
#include <dumux/implicit/3p3c/3p3cmodel.hh>
35
#include <dumux/implicit/common/implicitporousmediaproblem.hh>
Bernd Flemisch's avatar
Bernd Flemisch committed
36

37
#include "kuevettespatialparams.hh"
Bernd Flemisch's avatar
Bernd Flemisch committed
38
39
40
41
42
43
44
45
46
47

#define ISOTHERMAL 0

namespace Dumux
{
template <class TypeTag>
class KuevetteProblem;

namespace Properties
{
48
49
50
51
NEW_TYPE_TAG(KuevetteProblem, INHERITS_FROM(ThreePThreeCNI, KuevetteSpatialParams));
NEW_TYPE_TAG(KuevetteBoxProblem, INHERITS_FROM(BoxModel, KuevetteProblem));
NEW_TYPE_TAG(KuevetteCCProblem, INHERITS_FROM(CCModel, KuevetteProblem));
    
Bernd Flemisch's avatar
Bernd Flemisch committed
52
// Set the grid type
Thomas Fetzer's avatar
Thomas Fetzer committed
53
SET_TYPE_PROP(KuevetteProblem, Grid, Dune::YaspGrid<2>);
Bernd Flemisch's avatar
Bernd Flemisch committed
54
55

// Set the problem property
Thomas Fetzer's avatar
Thomas Fetzer committed
56
SET_TYPE_PROP(KuevetteProblem, Problem, Dumux::KuevetteProblem<TypeTag>);
Bernd Flemisch's avatar
Bernd Flemisch committed
57
58

// Set the fluid system
59
SET_TYPE_PROP(KuevetteProblem, 
Bernd Flemisch's avatar
Bernd Flemisch committed
60
61
62
63
              FluidSystem,
              Dumux::FluidSystems::H2OAirMesitylene<typename GET_PROP_TYPE(TypeTag, Scalar)>);

// Enable gravity
64
SET_BOOL_PROP(KuevetteProblem, ProblemEnableGravity, true);
Bernd Flemisch's avatar
Bernd Flemisch committed
65

Holger Class's avatar
Holger Class committed
66
// Use central differences (backward -1, forward +1)
67
SET_INT_PROP(KuevetteProblem, ImplicitNumericDifferenceMethod, 0);
Bernd Flemisch's avatar
Bernd Flemisch committed
68

Holger Class's avatar
Holger Class committed
69
// Set the maximum time step
70
SET_SCALAR_PROP(KuevetteProblem, TimeManagerMaxTimeStepSize, 60.);
Holger Class's avatar
Holger Class committed
71
72
73

// set newton relative tolerance
SET_SCALAR_PROP(KuevetteProblem, NewtonRelTolerance, 1e-6);
Bernd Flemisch's avatar
Bernd Flemisch committed
74
75
76
77
}


/*!
Thomas Fetzer's avatar
Thomas Fetzer committed
78
 * \ingroup ThreePThreeCModel
79
 * \ingroup ImplicitTestProblems
80
81
 * \brief Non-isothermal gas injection problem where a gas (e.g. steam/air)
 *        is injected into a unsaturated porous medium with a residually
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
 *        trapped NAPL contamination.
 *
 * The domain is a quasi-two-dimensional container (kuevette). Its dimensions
 * are 1.5 m x 0.74 m. The top and bottom boundaries are closed, the right
 * boundary is a Dirichlet boundary allowing fluids to escape. From the left,
 * an injection of a hot water-air mixture is applied (Neumann boundary condition
 * for the mass components and the enthalpy), aimed at remediating an initial 
 * NAPL (Non-Aquoeus Phase Liquid) contamination in the heterogeneous domain.
 * The contamination is initially placed partly into the coarse sand
 * and partly into a fine sand lense.
 *
 * This simulation can be varied through assigning different boundary conditions
 * at the left boundary as described in Class (2001):
 * Theorie und numerische Modellierung nichtisothermer Mehrphasenprozesse in
 * NAPL-kontaminierten por"osen Medien, Dissertation, Eigenverlag des Instituts
 * f"ur Wasserbau
 *
 * This problem uses the \ref ThreePThreeCNIModel
 *
 * To see the basic effect and the differences to scenarios with pure steam or
 * pure air injection, it is sufficient to simulated for about 2-3 hours (10000 s).
 * Complete remediation of the domain requires much longer (about 10 days simulated time).
 * To adjust the simulation time it is necessary to edit the file test_3p3cni.input
 *
 * To run the simulation execute:
Thomas Fetzer's avatar
Thomas Fetzer committed
107
108
 * <tt>./test_box3p3cni -parameterFile test_box3p3cni.input</tt> or
 * <tt>./test_cc3p3cni -parameterFile test_cc3p3cni.input</tt>
Bernd Flemisch's avatar
Bernd Flemisch committed
109
110
 *  */
template <class TypeTag >
111
class KuevetteProblem : public ImplicitPorousMediaProblem<TypeTag>
Bernd Flemisch's avatar
Bernd Flemisch committed
112
113
114
115
116
{
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GridView::Grid Grid;

117
    typedef ImplicitPorousMediaProblem<TypeTag> ParentType;
Bernd Flemisch's avatar
Bernd Flemisch committed
118
119

    // copy some indices for convenience
120
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
Bernd Flemisch's avatar
Bernd Flemisch committed
121
122
123
124
125
126
127
128
129
130
131
132
133
134
    enum {

        pressureIdx = Indices::pressureIdx,
        switch1Idx = Indices::switch1Idx,
        switch2Idx = Indices::switch2Idx,
        temperatureIdx = Indices::temperatureIdx,
        energyEqIdx = Indices::energyEqIdx,

        // Phase State
        threePhases = Indices::threePhases,
        wgPhaseOnly = Indices::wgPhaseOnly,

        // Grid and world dimension
        dim = GridView::dimension,
Andreas Lauser's avatar
Andreas Lauser committed
135
        dimWorld = GridView::dimensionworld
Bernd Flemisch's avatar
Bernd Flemisch committed
136
137
138
139
140
141
142
143
    };


    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
    typedef typename GET_PROP_TYPE(TypeTag, BoundaryTypes) BoundaryTypes;
    typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager;

    typedef typename GridView::template Codim<0>::Entity Element;
144
    typedef typename GridView::template Codim<0>::Iterator ElementIterator;
Bernd Flemisch's avatar
Bernd Flemisch committed
145
146
147
148
149
150
151
152
    typedef typename GridView::template Codim<dim>::Entity Vertex;
    typedef typename GridView::Intersection Intersection;

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

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

153
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };
154
    enum { dofCodim = isBox ? dim : 0 };
155

Bernd Flemisch's avatar
Bernd Flemisch committed
156
157
public:
    /*!
158
     * \brief The constructor.
Bernd Flemisch's avatar
Bernd Flemisch committed
159
160
161
162
     *
     * \param timeManager The time manager
     * \param gridView The grid view
     */
163
    KuevetteProblem(TimeManager &timeManager, const GridView &gridView)
164
        : ParentType(timeManager, gridView), eps_(1e-6)
Bernd Flemisch's avatar
Bernd Flemisch committed
165
166
    {
        FluidSystem::init();
167

168
        name_ = GET_RUNTIME_PARAM(TypeTag, std::string, Problem.Name);
Bernd Flemisch's avatar
Bernd Flemisch committed
169
170
171
172
173
174
175
176
177
178
179
180
    }

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

    /*!
     * \brief The problem name.
     *
     * This is used as a prefix for files generated by the simulation.
     */
181
182
    const std::string name() const
    { return name_; }
183
184
185
186
187
188
189

    /*!
     * \brief Returns the source term at specific position in the domain.
     *
     * \param values The source values for the primary variables
     * \param globalPos The position of the center of the finite volume
     */
Bernd Flemisch's avatar
Bernd Flemisch committed
190
191
192
193
194
195
    void sourceAtPos(PrimaryVariables &values,
                     const GlobalPosition &globalPos) const
    {
        values = 0;
    }

196

Bernd Flemisch's avatar
Bernd Flemisch committed
197
198
199
200
201
202
203
204
205
206
207
208
    // \}

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

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     *
     * \param values The boundary types for the conservation equations
209
     * \param globalPos The position for which the bc type should be evaluated
Bernd Flemisch's avatar
Bernd Flemisch committed
210
     */
211
212
    void boundaryTypesAtPos(BoundaryTypes &values, 
                            const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
213
214
215
216
217
218
219
220
221
222
223
224
225
    {
        if(globalPos[0] > 1.5 - eps_)
            values.setAllDirichlet();
        else
            values.setAllNeumann();

    }

    /*!
     * \brief Evaluate the boundary conditions for a dirichlet
     *        boundary segment.
     *
     * \param values The dirichlet values for the primary variables
226
     * \param globalPos The position for which the bc type should be evaluated
Bernd Flemisch's avatar
Bernd Flemisch committed
227
228
229
     *
     * For this method, the \a values parameter stores primary variables.
     */
230
    void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
231
232
233
234
235
236
237
238
239
240
    {
        initial_(values, globalPos);
    }

    /*!
     * \brief Evaluate the boundary conditions for a neumann
     *        boundary segment.
     *
     * \param values The neumann values for the conservation equations
     * \param element The finite element
241
     * \param fvGeomtry The finite-volume geometry in the box scheme
242
     * \param intersection The intersection between element and boundary
Bernd Flemisch's avatar
Bernd Flemisch committed
243
244
245
246
247
248
249
250
     * \param scvIdx The local vertex index
     * \param boundaryFaceIdx The index of the boundary face
     *
     * For this method, the \a values parameter stores the mass flux
     * in normal direction of each phase. Negative values mean influx.
     */
    void neumann(PrimaryVariables &values,
                 const Element &element,
251
                 const FVElementGeometry &fvGeomtry,
252
                 const Intersection &intersection,
253
                 const int scvIdx,
Bernd Flemisch's avatar
Bernd Flemisch committed
254
255
256
                 int boundaryFaceIdx) const
    {
        values = 0;
257

258
259
260
261
        GlobalPosition globalPos;
        if (isBox)
            globalPos = element.geometry().corner(scvIdx);
        else 
262
            globalPos = intersection.geometry().center();
263

Bernd Flemisch's avatar
Bernd Flemisch committed
264
265
266
        // negative values for injection
        if (globalPos[0] < eps_)
        {
Holger Class's avatar
Holger Class committed
267
            values[Indices::contiWEqIdx] = -0.1435; // 0.3435 [mol/(s m)] in total
268
269
            values[Indices::contiGEqIdx] = -0.2;
            values[Indices::contiNEqIdx] =  0.0;
Holger Class's avatar
Holger Class committed
270
            values[Indices::energyEqIdx] = -6929.;
Bernd Flemisch's avatar
Bernd Flemisch committed
271
272
273
274
275
276
277
278
279
280
281
282
283
284
        }
    }

    // \}

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

    /*!
     * \brief Evaluate the initial value for a control volume.
     *
     * \param values The initial values for the primary variables
285
     * \param globalPos The position for which the initial condition should be evaluated
Bernd Flemisch's avatar
Bernd Flemisch committed
286
287
288
289
     *
     * For this method, the \a values parameter stores primary
     * variables.
     */
290
    void initialAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
291
292
293
294
295
296
297
    {
        initial_(values, globalPos);
    }

    /*!
     * \brief Return the initial phase state inside a control volume.
     *
298
     * \param vertex The vertex
299
     * \param vIdxGlobal The global index of the vertex
Bernd Flemisch's avatar
Bernd Flemisch committed
300
301
     * \param globalPos The global position
     */
302
    int initialPhasePresence(const Vertex &vertex,
303
                             const int &vIdxGlobal,
Bernd Flemisch's avatar
Bernd Flemisch committed
304
305
306
307
308
309
310
311
312
                             const GlobalPosition &globalPos) const
    {
        if((globalPos[0] >= 0.20) && (globalPos[0] <= 0.80) && (globalPos[1] >= 0.4) && (globalPos[1] <= 0.65))
        {
            return threePhases;
        }
        else return wgPhaseOnly;
    }

313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
       /*!
     * \brief Append all quantities of interest which can be derived
     *        from the solution of the current time step to the VTK
     *        writer. Adjust this in case of anisotropic permeabilities.
     */
    void addOutputVtkFields()
    {
        // get the number of degrees of freedom
        unsigned numDofs = this->model().numDofs();

        // create the scalar field required for the permeabilities
        typedef Dune::BlockVector<Dune::FieldVector<double, 1> > ScalarField;
        ScalarField *Kxx = this->resultWriter().allocateManagedBuffer(numDofs);

        FVElementGeometry fvGeometry;

        ElementIterator eIt = this->gridView().template begin<0>();
        ElementIterator eEndIt = this->gridView().template end<0>();
        for (; eIt != eEndIt; ++eIt)
        {
            fvGeometry.update(this->gridView(), *eIt);

            for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx)
            {
337
338
339
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
                int dofIdxGlobal = this->model().dofMapper().subIndex(*eIt, scvIdx, dofCodim);
#else
340
                int dofIdxGlobal = this->model().dofMapper().map(*eIt, scvIdx, dofCodim);
341
#endif
342
                (*Kxx)[dofIdxGlobal] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, scvIdx);
343
344
345
346
347
348
            }
        }

        this->resultWriter().attachDofData(*Kxx, "permeability", isBox);
    }

Bernd Flemisch's avatar
Bernd Flemisch committed
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366

private:
    // internal method for the initial condition (reused for the
    // dirichlet conditions!)
    void initial_(PrimaryVariables &values,
                  const GlobalPosition &globalPos) const
    {
        values[pressureIdx] = 1e5 ;
        values[switch1Idx] = 0.12;
        values[switch2Idx] = 1.e-6;
        values[temperatureIdx] = 293.0;

        if((globalPos[0] >= 0.20) && (globalPos[0] <= 0.80) && (globalPos[1] >= 0.4) && (globalPos[1] <= 0.65))
        {
            values[switch2Idx] = 0.07;
        }
    }

367
    const Scalar eps_;
368
    std::string name_;
Bernd Flemisch's avatar
Bernd Flemisch committed
369
370
371
372
};
} //end namespace

#endif