heterogeneousproblem.hh 19.4 KB
Newer Older
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.                      *
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            *
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
 * \brief Definition of a problem, where CO2 is injected in a reservoir.
23
24
25
26
27
28
 */
#ifndef DUMUX_HETEROGENEOUS_PROBLEM_HH
#define DUMUX_HETEROGENEOUS_PROBLEM_HH

#if HAVE_ALUGRID
#include <dune/grid/alugrid/2d/alugrid.hh>
29
30
#elif HAVE_DUNE_ALUGRID
#include <dune/alugrid/grid.hh>
31
#else
Christoph Grueninger's avatar
[co2ni]    
Christoph Grueninger committed
32
#warning ALUGrid is necessary for this test.
33
#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
34
35
#endif

36
37
#include <dumux/implicit/co2/co2model.hh>
#include <dumux/implicit/co2/co2volumevariables.hh>
38
#include <dumux/material/fluidsystems/brineco2fluidsystem.hh>
39
#include <dumux/implicit/common/implicitporousmediaproblem.hh>
40
#include <dumux/implicit/box/intersectiontovertexbc.hh>
41

42

43
44
45
46
47
48
49
50
51
52
53
#include "heterogeneousspatialparameters.hh"
#include "heterogeneousco2tables.hh"

namespace Dumux
{

template <class TypeTag>
class HeterogeneousProblem;

namespace Properties
{
54
55
56
NEW_TYPE_TAG(HeterogeneousProblem, INHERITS_FROM(TwoPTwoC, HeterogeneousSpatialParams));
NEW_TYPE_TAG(HeterogeneousBoxProblem, INHERITS_FROM(BoxModel, HeterogeneousProblem));
NEW_TYPE_TAG(HeterogeneousCCProblem, INHERITS_FROM(CCModel, HeterogeneousProblem));
57
58

// Set the grid type
59
#if HAVE_ALUGRID || HAVE_DUNE_ALUGRID
60
61
62
63
64
65
SET_TYPE_PROP(HeterogeneousProblem, Grid, Dune::ALUGrid<2, 2, Dune::cube, Dune::nonconforming>);
#else
SET_TYPE_PROP(HeterogeneousProblem, Grid, Dune::YaspGrid<2>);
#endif

// Set the problem property
66
SET_TYPE_PROP(HeterogeneousProblem, Problem, Dumux::HeterogeneousProblem<TypeTag>);
67
68

// Set fluid configuration
69
SET_TYPE_PROP(HeterogeneousProblem, FluidSystem, Dumux::BrineCO2FluidSystem<TypeTag>);
70
71
72

// Set the CO2 table to be used; in this case not the the default table
SET_TYPE_PROP(HeterogeneousProblem, CO2Table, Dumux::HeterogeneousCO2Tables::CO2Tables);
73

74
// Set the salinity mass fraction of the brine in the reservoir
75
SET_SCALAR_PROP(HeterogeneousProblem, ProblemSalinity, 1e-1);
76
77
78
79
80

//! the CO2 Model and VolumeVariables properties
SET_TYPE_PROP(HeterogeneousProblem, Model, CO2Model<TypeTag>);
SET_TYPE_PROP(HeterogeneousProblem, VolumeVariables, CO2VolumeVariables<TypeTag>);

81
// Use Moles
82
SET_BOOL_PROP(HeterogeneousProblem, UseMoles, false);
83
84
85
86
87
}


/*!
 * \ingroup CO2Model
88
 * \ingroup ImplicitTestProblems
89
 * \brief Definition of a problem, where CO2 is injected in a reservoir.
90
91
92
93
94
95
96
97
98
99
100
101
 *
 * The domain is sized 200m times 100m and consists of four layers, a
 * permeable reservoir layer at the bottom, a barrier rock layer with reduced permeability, another reservoir layer
 * and at the top a barrier rock layer with a very low permeablility.
 *
 * CO2 is injected at the permeable bottom layer
 * from the left side. The domain is initially filled with brine.
 *
 * The grid is unstructered and permeability and porosity for the elements are read in from the grid file. The grid file
 * also contains so-called boundary ids which can be used assigned during the grid creation in order to differentiate
 * between different parts of the boundary.
 * These boundary ids can be imported into the problem where the boundary conditions can then be assigned accordingly.
Thomas Fetzer's avatar
Thomas Fetzer committed
102
 * 
103
104
105
 * 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 setup. The default setting for useMoles is false.
 *
106
 * To run the simulation execute the following line in shell (works with the box and cell centered spatial discretization method):
107
 * <tt>./test_ccco2 </tt> or <tt>./test_boxco2 </tt>
108
 */
109
template <class TypeTag >
110
class HeterogeneousProblem : public ImplicitPorousMediaProblem<TypeTag>
111
112
{

113
    typedef ImplicitPorousMediaProblem<TypeTag> ParentType;
114
115
116
117
118
119
120
121
122
123
124
125
126
127
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, Grid) Grid;
    typedef Dune::GridPtr<Grid> GridPointer;
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) VolumeVariables;

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

    // copy some indices for convenience
128
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
129
130
131
132
    enum {
        lPhaseIdx = Indices::wPhaseIdx,
        gPhaseIdx = Indices::nPhaseIdx,

133
134
        wCompIdx = FluidSystem::wCompIdx,
        nCompIdx = FluidSystem::nCompIdx,
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155

        BrineIdx = FluidSystem::BrineIdx,
        CO2Idx = FluidSystem::CO2Idx,

        conti0EqIdx = Indices::conti0EqIdx,
        contiCO2EqIdx = conti0EqIdx + CO2Idx
    };

    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;
    typedef typename GridView::template Codim<0>::Iterator ElementIterator;
    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, GridCreator) GridCreator;

    typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition;
156
    typedef typename GET_PROP_TYPE(TypeTag, CO2Table) CO2Table;
157
    typedef Dumux::CO2<Scalar, CO2Table> CO2;
158
159
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };
    enum { dofCodim = isBox ? dim : 0 };
160
    //! property that defines whether mole or mass fractions are used
161
    static const bool useMoles = GET_PROP_VALUE(TypeTag, UseMoles);
162
163
164
165
166
167
168
169
170
171

public:
    /*!
     * \brief The constructor
     *
     * \param timeManager The time manager
     * \param gridView The grid view
     */
    HeterogeneousProblem(TimeManager &timeManager,
                     const GridView &gridView)
172
        : ParentType(timeManager, GridCreator::grid().leafGridView()),
173
174
175
176
177
178
          //Boundary Id Setup:
          injectionTop_(1),
          injectionBottom_(2),
          dirichletBoundary_(3),
          noFlowBoundary_(4),
          intersectionToVertexBC_(*this)
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
    {
            nTemperature_       = GET_RUNTIME_PARAM(TypeTag, int, FluidSystem.NTemperature);
            nPressure_          = GET_RUNTIME_PARAM(TypeTag, int, FluidSystem.NPressure);
            pressureLow_        = GET_RUNTIME_PARAM(TypeTag, Scalar, FluidSystem.PressureLow);
            pressureHigh_       = GET_RUNTIME_PARAM(TypeTag, Scalar, FluidSystem.PressureHigh);
            temperatureLow_     = GET_RUNTIME_PARAM(TypeTag, Scalar, FluidSystem.TemperatureLow);
            temperatureHigh_    = GET_RUNTIME_PARAM(TypeTag, Scalar, FluidSystem.TemperatureHigh);
            depthBOR_           = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.DepthBOR);
            name_               = GET_RUNTIME_PARAM(TypeTag, std::string, Problem.Name);
            injectionRate_      = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.InjectionRate);

        /* Alternative syntax:
         * typedef typename GET_PROP(TypeTag, ParameterTree) ParameterTree;
         * const Dune::ParameterTree &tree = ParameterTree::tree();
         * nTemperature_       = tree.template get<int>("FluidSystem.nTemperature");
         *
         * + We see what we do
         * - Reporting whether it was used does not work
         * - Overwriting on command line not possible
        */
        GridPointer *gridPtr = &GridCreator::gridPtr();
        this->spatialParams().setParams(gridPtr);

        eps_ = 1e-6;

        // initialize the tables of the fluid system
        FluidSystem::init(/*Tmin=*/temperatureLow_,
                          /*Tmax=*/temperatureHigh_,
                          /*nT=*/nTemperature_,
                          /*pmin=*/pressureLow_,
                          /*pmax=*/pressureHigh_,
                          /*np=*/nPressure_);
211
212
213
214
215
216
217
218
219
220

        //stateing in the console whether mole or mass fractions are used
        if(!useMoles)
        {
        	std::cout<<"problem uses mass-fractions"<<std::endl;
        }
        else
        {
        	std::cout<<"problem uses mole-fractions"<<std::endl;
        }
221
222
223
    }

    /*!
224
225
226
     * \brief User defined output after the time integration
     *
     * Will be called diretly after the time integration.
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
     */
    void postTimeStep()
    {
        // Calculate storage terms
        PrimaryVariables storageL, storageG;
        this->model().globalPhaseStorage(storageL, lPhaseIdx);
        this->model().globalPhaseStorage(storageG, gPhaseIdx);

        // Write mass balance information for rank 0
        if (this->gridView().comm().rank() == 0) {
            std::cout<<"Storage: liquid=[" << storageL << "]"
                     << " gas=[" << storageG << "]\n";
        }
    }

242
243
244
245
246
    /*!
     * \brief Append all quantities of interest which can be derived
     *        from the solution of the current time step to the VTK
     *        writer.
     */
247
248
    void addOutputVtkFields()
        {
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
        typedef Dune::BlockVector<Dune::FieldVector<double, 1> > ScalarField;

         // get the number of degrees of freedom
         unsigned numDofs = this->model().numDofs();
         unsigned numElements = this->gridView().size(0);

         //create required scalar fields
         ScalarField *Kxx = this->resultWriter().allocateManagedBuffer(numElements);
         ScalarField *cellPorosity = this->resultWriter().allocateManagedBuffer(numElements);
         ScalarField *boxVolume = this->resultWriter().allocateManagedBuffer(numDofs);
         (*boxVolume) = 0;

         //Fill the scalar fields with values
         ScalarField *rank = this->resultWriter().allocateManagedBuffer(numElements);

         FVElementGeometry fvGeometry;
         VolumeVariables volVars;

267
268
269
         ElementIterator eIt = this->gridView().template begin<0>();
         ElementIterator eEndIt = this->gridView().template end<0>();
         for (; eIt != eEndIt; ++eIt)
270
         {
271
272
             int eIdx = this->elementMapper().map(*eIt);
             (*rank)[eIdx] = this->gridView().comm().rank();
273
             fvGeometry.update(this->gridView(), *eIt);
274
275
276
277


             for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx)
             {
278
279
                 int dofIdxGlobal = this->model().dofMapper().map(*eIt, scvIdx, dofCodim);
                 volVars.update(this->model().curSol()[dofIdxGlobal],
280
                                *this,
281
                                *eIt,
282
283
284
                                fvGeometry,
                                scvIdx,
                                false);
285
                 (*boxVolume)[dofIdxGlobal] += fvGeometry.subContVol[scvIdx].volume;
286
             }
287
288
             (*Kxx)[eIdx] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, /*element data*/ 0);
             (*cellPorosity)[eIdx] = this->spatialParams().porosity(*eIt, fvGeometry, /*element data*/ 0);
289
290
291
292
293
294
         }

         //pass the scalar fields to the vtkwriter
         this->resultWriter().attachDofData(*Kxx, "Kxx", false); //element data
         this->resultWriter().attachDofData(*cellPorosity, "cellwisePorosity", false); //element data
         this->resultWriter().attachDofData(*boxVolume, "boxVolume", isBox);
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
        }
    /*!
     * \name Problem parameters
     */
    // \{

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

    /*!
     * \brief Returns the temperature within the domain.
     *
312
313
314
315
     * \param globalPos The position
     *
     * This problem assumes a geothermal gradient with 
     * a surface temperature of 10 degrees Celsius.
316
     */
317
    Scalar temperatureAtPos(const GlobalPosition &globalPos) const
318
319
    {
        return temperature_(globalPos);
320
    }
321

322
    /*!
323
     * \brief Returns the source term
324
     *
325
326
     * \param values Stores the source values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} / (m^\textrm{dim} \cdot s )] \f$
327
328
     * \param globalPos The global position
     */
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
    void sourceAtPos(PrimaryVariables &values,
                const GlobalPosition &globalPos) const
    {
        values = 0;
    }

    // \}

    /*!
     * \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
     * \param vertex The vertex for which the boundary type is set
     */
    void boundaryTypes(BoundaryTypes &values, const Vertex &vertex) const
    {
        intersectionToVertexBC_.boundaryTypes(values, vertex);
    }

    /*!
     * \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
359
     * \param intersection specifies the intersection at which boundary
360
     *           condition is to set
361
     */
362
    void boundaryTypes(BoundaryTypes &values, const Intersection &intersection) const
363
    {
364
        int boundaryId = intersection.boundaryId();
365
366
367
368
369
370
371
372
373
374
375
376
        if (boundaryId < 1 || boundaryId > 4)
        {
            std::cout<<"invalid boundaryId: "<<boundaryId<<std::endl;
            DUNE_THROW(Dune::InvalidStateException, "Invalid " << boundaryId);
        }
        if (boundaryId == dirichletBoundary_)
            values.setAllDirichlet();
        else
            values.setAllNeumann();
    }

    /*!
377
378
     * \brief Evaluates the boundary conditions for a Dirichlet
     *        boundary segment
379
     *
380
381
     * \param values Stores the Dirichlet values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} ] \f$
382
     * \param globalPos The global position
383
384
385
386
387
     */
    void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
    {
        initial_(values, globalPos);
    }
388

389
    /*!
390
     * \brief Evaluate the boundary conditions for a Neumann
391
392
     *        boundary segment.
     *
393
394
      * \param values Stores the Neumann values for the conservation equations in
     *               \f$ [ \textnormal{unit of conserved quantity} / (m^(dim-1) \cdot s )] \f$
395
     * \param element The finite element
396
     * \param fvGeometry The finite volume geometry of the element
397
     * \param intersection The intersection between element and boundary
398
     * \param scvIdx The local index of the sub-control volume
399
400
     * \param boundaryFaceIdx The index of the boundary face
     *
401
402
     * The \a values store the mass flux of each phase normal to the boundary.
     * Negative values indicate an inflow.
403
404
     *
     * Depending on whether useMoles is set on true or false, the flux has to be given either in
405
     * kg/(m^2*s) or mole/(m^2*s) in the input file!! Convertion with molar mass obtained from fluid system FluidSystem::molarMass(nCompIdx)
406
407
408
     */
    void neumann(PrimaryVariables &values,
                 const Element &element,
409
                 const FVElementGeometry &fvGeometry,
410
                 const Intersection &intersection,
411
412
413
                 int scvIdx,
                 int boundaryFaceIdx) const
    {
414
        int boundaryId = intersection.boundaryId();
415
416
417
418

        values = 0;
        if (boundaryId == injectionBottom_)
        {
419
            values[contiCO2EqIdx] = -injectionRate_; //see above: either give in kg/(m^2*s) or mole/(m^2*s) depending on useMoles
420
421
422
423
424
425
426
427
428
429
430
        }
    }

    // \}

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

    /*!
431
     * \brief Evaluates the initial values for a control volume
432
     *
433
434
435
     * \param values Stores the initial values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variables} ] \f$
     * \param globalPos The global position
436
     */
437
438
    void initialAtPos(PrimaryVariables &values,
                      const GlobalPosition &globalPos) const
439
440
441
442
443
    {
        initial_(values, globalPos);
    }

    /*!
444
     * \brief Returns the initial phase state for a control volume.
445
     *
446
     * \param vertex The vertex
447
     * \param vIdxGlobal The global index of the vertex
448
449
     * \param globalPos The global position
     */
450
    int initialPhasePresence(const Vertex &vertex,
451
                             int &vIdxGlobal,
452
                             const GlobalPosition &globalPos) const
453
    { return Indices::wPhaseOnly; }
454
455
456
457

    // \}

private:
458
459
460
461
462
463
464
465
466
    /*!
     * \brief Evaluates the initial values for a control volume
     *
     * The internal method for the initial condition
     *
     * \param values Stores the initial values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variables} ] \f$
     * \param globalPos The global position
     */
467
468
469
470
471
472
    void initial_(PrimaryVariables &values,
                  const GlobalPosition &globalPos) const
    {
        Scalar temp = temperature_(globalPos);
        Scalar densityW = FluidSystem::Brine::liquidDensity(temp, 1e7);

473
        Scalar pl =  1e5 - densityW*this->gravity()[dimWorld-1]*(depthBOR_ - globalPos[dimWorld-1]);
474
475
476
477
478
479
        Scalar moleFracLiquidCO2 = 0.00;
        Scalar moleFracLiquidBrine = 1.0 - moleFracLiquidCO2;

        Scalar meanM =
            FluidSystem::molarMass(BrineIdx)*moleFracLiquidBrine +
            FluidSystem::molarMass(CO2Idx)*moleFracLiquidCO2;
480
481
482
483
484
485
486
487
488
        if(!useMoles) //mass-fraction formulation
        {
        	Scalar massFracLiquidCO2 = moleFracLiquidCO2*FluidSystem::molarMass(CO2Idx)/meanM;
        	values[Indices::switchIdx] = massFracLiquidCO2;
        }
        else //mole-fraction formulation
        {
        	values[Indices::switchIdx] = moleFracLiquidCO2;
        }
489
490
491
492
493
        values[Indices::pressureIdx] = pl;
    }

    Scalar temperature_(const GlobalPosition globalPos) const
    {
494
        Scalar T = 283.0 + (depthBOR_ - globalPos[dimWorld-1])*0.03; 
495
        return T;
496
    }
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519

    Scalar depthBOR_;
    Scalar injectionRate_;
    Scalar eps_;

    int nTemperature_;
    int nPressure_;

    std::string name_ ;

    Scalar pressureLow_, pressureHigh_;
    Scalar temperatureLow_, temperatureHigh_;

    int injectionTop_;
    int injectionBottom_;
    int dirichletBoundary_;
    int noFlowBoundary_;

    const IntersectionToVertexBC<TypeTag> intersectionToVertexBC_;
};
} //end namespace

#endif