2pproblem.hh 12.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
// -*- 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 Tutorial problem for a fully coupled twophase box model.
 */
24
25
#ifndef DUMUX_EXERCISE_FLUIDSYSTEM_A_PROBLEM_HH
#define DUMUX_EXERCISE_FLUIDSYSTEM_A_PROBLEM_HH
26
27
28
29
30

// The numerical model
#include <dumux/porousmediumflow/2p/model.hh>

// The box discretization
31
#include <dumux/discretization/box.hh>
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61

// The grid managers
#if HAVE_DUNE_ALUGRID
#include <dune/alugrid/grid.hh>
#elif HAVE_UG
#include <dune/grid/uggrid.hh>
#else
#include <dune/grid/yaspgrid.hh>
#endif

// The porous media base problem
#include <dumux/porousmediumflow/problem.hh>

// Spatially dependent parameters
#include "spatialparams.hh"

// The water component
#include <dumux/material/components/tabulatedcomponent.hh>
#include <dumux/material/components/h2o.hh>

// The components that will be created in this exercise
#include "components/myincompressiblecomponent.hh"
#include "components/mycompressiblecomponent.hh"

// We will only have liquid phases here
#include <dumux/material/fluidsystems/1pliquid.hh>

// The two-phase immiscible fluid system
#include <dumux/material/fluidsystems/2pimmiscible.hh>

62
63
64
// The interface to create plots during simulation
#include <dumux/io/gnuplotinterface.hh>

65
66
namespace Dumux{
// Forward declaration of the problem class
67
template <class TypeTag> class ExerciseFluidsystemProblemTwoP;
68
69
70

namespace Properties {
// Create a new type tag for the problem
71
72
73
74
// Create new type tags
namespace TTag {
struct ExerciseFluidsystemTwoPTypeTag { using InheritsFrom = std::tuple<BoxModel, TwoP>; };
} // end namespace TTag
75
76

// Set the "Problem" property
77
78
template<class TypeTag>
struct Problem<TypeTag, TTag::ExerciseFluidsystemTwoPTypeTag> { using type = ExerciseFluidsystemProblemTwoP<TypeTag>; };
79
80

// Set the spatial parameters
81
SET_TYPE_PROP(ExerciseFluidsystemTwoPTypeTag, SpatialParams,
82
83
              ExerciseFluidsystemSpatialParams<GetPropType<TypeTag, Properties::FVGridGeometry>,
                                         GetPropType<TypeTag, Properties::Scalar>>);
84
85
86

// Set grid to be used
#if HAVE_DUNE_ALUGRID
87
88
template<class TypeTag>
struct Grid<TypeTag, TTag::ExerciseFluidsystemTwoPTypeTag> { using type = Dune::ALUGrid</*dim=*/2, 2, Dune::cube, Dune::nonconforming>; };
89
#elif HAVE_UG
90
91
template<class TypeTag>
struct Grid<TypeTag, TTag::ExerciseFluidsystemTwoPTypeTag> { using type = Dune::UGGrid<2>; };
92
#else
93
94
template<class TypeTag>
struct Grid<TypeTag, TTag::ExerciseFluidsystemTwoPTypeTag> { using type = Dune::YaspGrid<2>; };
95
96
97
#endif // HAVE_DUNE_ALUGRID

// we use the immiscible fluid system here
98
99
template<class TypeTag>
struct FluidSystem<TypeTag, TTag::ExerciseFluidsystemTwoPTypeTag>
100
101
{
private:
102
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
103
    using TabulatedH2O = Components::TabulatedComponent<Components::H2O<Scalar>>;
104
    using LiquidWater = typename FluidSystems::OnePLiquid<Scalar, TabulatedH2O>;
105
106
107
108
    /*!
     * Uncomment first line and comment second line for using the incompressible component
     * Uncomment second line and comment first line for using the compressible component
     */
Ned Coltman's avatar
Ned Coltman committed
109
110
    using LiquidMyComponentPhase = typename FluidSystems::OnePLiquid<Scalar, MyIncompressibleComponent<Scalar> >;
    // using LiquidMyComponentPhase = typename FluidSystems::OnePLiquid<Scalar, MyCompressibleComponent<Scalar> >;
111
112

public:
Ned Coltman's avatar
Ned Coltman committed
113
    using type = typename FluidSystems::TwoPImmiscible<Scalar, LiquidWater, LiquidMyComponentPhase>;
114
115
116
117
118
119
120
121
122
};

}

/*!
 * \ingroup TwoPBoxModel
 * \brief  Tutorial problem for a fully coupled twophase box model.
 */
template <class TypeTag>
123
class ExerciseFluidsystemProblemTwoP : public PorousMediumFlowProblem<TypeTag>
124
125
{
    using ParentType = PorousMediumFlowProblem<TypeTag>;
126
127
    using GridView = GetPropType<TypeTag, Properties::GridView>;
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
128
129
130
131
132
133
134
135
136

    // Grid dimension
    enum { dim = GridView::dimension,
           dimWorld = GridView::dimensionworld
    };
    using Element = typename GridView::template Codim<0>::Entity;
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;

    // Dumux specific types
137
138
139
140
141
142
143
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
    using FVGridGeometry = GetPropType<TypeTag, Properties::FVGridGeometry>;
    using FVElementGeometry = typename GetPropType<TypeTag, Properties::FVGridGeometry>::LocalView;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using FluidState = GetPropType<TypeTag, Properties::FluidState>;
144
145
146
147
148
149
150
151
152

    enum {
        waterPressureIdx = Indices::pressureIdx,
        naplSaturationIdx = Indices::saturationIdx,
        contiWEqIdx = Indices::conti0EqIdx + FluidSystem::comp0Idx, // water transport equation index
        contiNEqIdx = Indices::conti0EqIdx + FluidSystem::comp1Idx // napl transport equation index
    };

public:
153
    ExerciseFluidsystemProblemTwoP(std::shared_ptr<const FVGridGeometry> fvGridGeometry)
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
    : ParentType(fvGridGeometry)
    , eps_(3e-6)
    {
#if !(HAVE_DUNE_ALUGRID || HAVE_UG)
        std::cout << "If you want to use simplices instead of cubes, install and use dune-ALUGrid or UGGrid." << std::endl;
#endif // !(HAVE_DUNE_ALUGRID || HAVE_UG)

        // initialize the tables for the water properties
        std::cout << "Initializing the tables for the water properties" << std::endl;
        Components::TabulatedComponent<Components::H2O<Scalar>>::init(/*tempMin=*/273.15,
                                                                      /*tempMax=*/623.15,
                                                                      /*numTemp=*/100,
                                                                      /*pMin=*/0.0,
                                                                      /*pMax=*/20e6,
                                                                      /*numP=*/200);

        // set the depth of the bottom of the reservoir
        depthBOR_ = this->fvGridGeometry().bBoxMax()[dimWorld-1];
172
173
174
175

        // plot density over pressure of the phase consisting of your component
        if(getParam<bool>("Output.PlotDensity"))
            plotDensity_();
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
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
    }

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

    /*!
     * \brief Returns the temperature \f$ K \f$
     */
    Scalar temperature() const
    { return 283.15; }

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

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     *
     * \param bcTypes The boundary types for the conservation equations
     * \param globalPos The position for which the bc type should be evaluated
     */
    BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
    {
         BoundaryTypes bcTypes;

        if (globalPos[0] < eps_ || globalPos[0] > this->fvGridGeometry().bBoxMax()[0] - eps_)
           bcTypes.setAllDirichlet();
        else
            bcTypes.setAllNeumann();

        return bcTypes;
    }

    /*!
     * \brief Evaluates the boundary conditions for a Dirichlet
     *        boundary segment
     *
     * \param globalPos The global position
     */
    PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const
    {
        // use the initial values as Dirichlet values
        return initialAtPos(globalPos);
    }

    /*!
     * \brief Evaluate the boundary conditions for a neumann
     *        boundary segment.
     *
     * \param values Stores the Neumann values for the conservation equations in
     *               \f$ [ \textnormal{unit of conserved quantity} / (m^(dim-1) \cdot s )] \f$
     * \param globalPos The position of the integration point of the boundary segment.
     *
     * For this method, the \a values parameter stores the mass flux
     * in normal direction of each phase. Negative values mean influx.
     */
    PrimaryVariables neumannAtPos(const GlobalPosition& globalPos) const
    {
        // initialize values to zero, i.e. no-flow Neumann boundary conditions
        PrimaryVariables values(0.0);

        Scalar up = this->fvGridGeometry().bBoxMax()[dimWorld-1];

        // influx of oil (30 g/m/s) over a segment of the top boundary
        if (globalPos[dimWorld-1] > up - eps_ && globalPos[0] > 20 && globalPos[0] < 40)
        {
            values[contiWEqIdx] = 0;
            values[contiNEqIdx] = -3e-2;
        }
        else
        {
            // no-flow on the remaining Neumann-boundaries.
            values[contiWEqIdx] = 0;
            values[contiNEqIdx] = 0;
        }

        return values;
    }

    // \}

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

    /*!
     * \brief Evaluate the initial value for a control volume.
     *
     * \param globalPos The position for which the initial condition should be evaluated
     *
     * For this method, the \a values parameter stores primary
     * variables.
     */
    PrimaryVariables initialAtPos(const GlobalPosition& globalPos) const

    {
        PrimaryVariables values(0.0);

        // use hydrostatic pressure distribution with 2 bar at the top and zero saturation
        values[waterPressureIdx] = 200.0e3 + 9.81*1000*(depthBOR_ - globalPos[dimWorld-1]);
        values[naplSaturationIdx] = 0.0;

        return values;
    }
    // \}

    /*!
     * \brief Returns the source term
     *
     * \param values Stores the source values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} / (m^\textrm{dim} \cdot s )] \f$
     * \param globalPos The global position
     */
    PrimaryVariables sourceAtPos(const GlobalPosition& globalPos) const
    {
        // we define no source term here
        PrimaryVariables values(0.0);
        return values;
    }

private:
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
    void plotDensity_()
    {
        FluidState fluidState;
        fluidState.setTemperature(temperature());
        int numberOfPoints = 100;
        Scalar xMin = 1e4;
        Scalar xMax = 1e7;
        Scalar spacing = std::pow((xMax/xMin), 1.0/(numberOfPoints-1));
        std::vector<double> x(numberOfPoints);
        std::vector<double> y(numberOfPoints);
        for (int i=0; i<numberOfPoints; ++i)
            x[i] = xMin*std::pow(spacing,i);
        for (int i=0; i<numberOfPoints; ++i)
        {
            fluidState.setPressure(FluidSystem::phase1Idx, x[i]);
            y[i] = FluidSystem::density(fluidState, FluidSystem::phase1Idx);
        }
        gnuplot_.resetPlot();
        gnuplot_.setXRange(xMin, xMax);
        gnuplot_.setOption("set logscale x 10");
        gnuplot_.setOption("set key left top");
        gnuplot_.setYRange(1440, 1480);
        gnuplot_.setXlabel("pressure [Pa]");
        gnuplot_.setYlabel("density [kg/m^3]");
        gnuplot_.addDataSetToPlot(x, y, "YourComponentPhase_density.dat", "w l t 'Phase consisting of your component'");
        gnuplot_.plot("YourComponentPhase_density");
    }

330
331
    Scalar eps_; //! small epsilon value
    Scalar depthBOR_; //! depth at the bottom of the reservoir
332
    Dumux::GnuplotInterface<double> gnuplot_; //! collects data for plotting
333
334
335
336
};
}

#endif