heterogeneousproblemni.hh 20.7 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
29
 */
#ifndef DUMUX_HETEROGENEOUS_NI_PROBLEM_NI_HH
#define DUMUX_HETEROGENEOUS_NI_PROBLEM_NI_HH

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

#include <dune/grid/io/file/dgfparser/dgfs.hh>
34
35
#include <dumux/implicit/2p2c/2p2cmodel.hh>
#include <dumux/implicit/co2/co2volumevariables.hh>
36
#include <dumux/material/fluidsystems/brineco2fluidsystem.hh>
37
#include <dumux/implicit/common/implicitporousmediaproblem.hh>
38
#include <dumux/implicit/box/intersectiontovertexbc.hh>
39
#include <test/implicit/co2/heterogeneousspatialparameters.hh>
40
41
42
43
44
45
46
47
48

#include "heterogeneousco2tables.hh"

#define ISOTHERMAL 0

namespace Dumux
{

template <class TypeTag>
49
class HeterogeneousNIProblem;
50
51
52

namespace Properties
{
53
54
55
NEW_TYPE_TAG(HeterogeneousNIProblem, INHERITS_FROM(TwoPTwoCNI, HeterogeneousSpatialParams));
NEW_TYPE_TAG(HeterogeneousNIBoxProblem, INHERITS_FROM(BoxModel, HeterogeneousNIProblem));
NEW_TYPE_TAG(HeterogeneousNICCProblem, INHERITS_FROM(CCModel, HeterogeneousNIProblem));
56
57
58
59


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

// Set the problem property
66
SET_PROP(HeterogeneousNIProblem, Problem)
67
{
68
    typedef Dumux::HeterogeneousNIProblem<TypeTag> type;
69
70
71
};

// Set fluid configuration
72
SET_PROP(HeterogeneousNIProblem, FluidSystem)
73
74
75
76
77
{
    typedef Dumux::BrineCO2FluidSystem<TypeTag> type;
};

// Set the CO2 table to be used; in this case not the the default table
78
SET_TYPE_PROP(HeterogeneousNIProblem, CO2Table, Dumux::Heterogeneous::CO2Tables);
79
// Set the salinity mass fraction of the brine in the reservoir
80
SET_SCALAR_PROP(HeterogeneousNIProblem, ProblemSalinity, 1e-1);
81
82

//! the CO2 Model and VolumeVariables properties
83
SET_TYPE_PROP(HeterogeneousNIProblem, IsothermalVolumeVariables, CO2VolumeVariables<TypeTag>);
84

85
// Use Moles
86
SET_BOOL_PROP(HeterogeneousNIProblem, UseMoles, false);
87
88
89
90
91
}


/*!
 * \ingroup CO2NIModel
92
 * \ingroup ImplicitTestProblems
93
 * \brief Definition of a problem, where CO2 is injected in a reservoir.
94
95
96
97
98
99
100
101
102
103
104
105
106
 *
 * 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 followed by 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.
 *
107
108
109
 * 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.
 *
110
111
 * To run the simulation execute the following line in shell (works with the box and cell centered spatial discretization method):
 * <tt>./test_ccco2ni </tt> or <tt>./test_boxco2ni </tt>
112
 */
113
template <class TypeTag >
114
class HeterogeneousNIProblem : public ImplicitPorousMediaProblem<TypeTag>
115
{
116
    typedef ImplicitPorousMediaProblem<TypeTag> ParentType;
117
118
119
120
121
122
123
124

    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;

125
126
    static const bool useMoles = GET_PROP_VALUE(TypeTag, UseMoles);

127
128
129
130
131
132
133
    enum {
        // Grid and world dimension
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld
    };

    // copy some indices for convenience
134
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
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
161
162
163
164
165
166
167
    enum {
        lPhaseIdx = Indices::wPhaseIdx,
        gPhaseIdx = Indices::nPhaseIdx,


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

        conti0EqIdx = Indices::conti0EqIdx,
        contiCO2EqIdx = conti0EqIdx + CO2Idx,
#if !ISOTHERMAL
        temperatureIdx = CO2Idx +1,
        energyEqIdx = contiCO2EqIdx +1,
#endif

    };


    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;
    typedef typename GET_PROP_TYPE(TypeTag, PTAG(CO2Table)) CO2Table;
    typedef Dumux::CO2<Scalar, CO2Table> CO2;
168
169
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };
    enum { dofCodim = isBox ? dim : 0 };
170
171
172
173
174
175
176
177

public:
    /*!
     * \brief The constructor
     *
     * \param timeManager The time manager
     * \param gridView The grid view
     */
178
    HeterogeneousNIProblem(TimeManager &timeManager,
179
                     const GridView &gridView)
180
        : ParentType(timeManager, GridCreator::grid().leafGridView()),
181
182
183
184
185
186
          //Boundary Id Setup:
          injectionTop_ (1),
          injectionBottom_(2),
          dirichletBoundary_(3),
          noFlowBoundary_(4),
          intersectionToVertexBC_(*this)
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
    {
        try
        {
            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);
            injectionPressure_ = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.InjectionPressure);
            injectionTemperature_ = GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.InjectionTemperature);
        }
        catch (Dumux::ParameterException &e) {
            std::cerr << e << ". Abort!\n";
            exit(1) ;
        }
        catch (...) {
            std::cerr << "Unknown exception thrown!\n";
            exit(1);
        }

        /* 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();
        FluidSystem::init(/*Tmin=*/temperatureLow_,
                          /*Tmax=*/temperatureHigh_,
                          /*nT=*/nTemperature_,
                          /*pmin=*/pressureLow_,
                          /*pmax=*/pressureHigh_,
                          /*np=*/nPressure_);
236
237
238
239
240
241
242
243
244
245

        //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;
        }
246
247
248
    }

    /*!
249
250
251
     * \brief User defined output after the time integration
     *
     * Will be called diretly after the time integration.
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
     */
    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";
        }
    }

267
    /*!
268
269
270
     * \brief Append all quantities of interest which can be derived
     *        from the solution of the current time step to the VTK
     *        writer.
271
     */
272
273
274
    void addOutputVtkFields()
        {
            typedef Dune::BlockVector<Dune::FieldVector<double, 1> > ScalarField;
275
276
277
278

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

            //create required scalar fields
281
282
283
284
285
            ScalarField *Kxx = this->resultWriter().allocateManagedBuffer(numElements);
            ScalarField *cellPorosity = this->resultWriter().allocateManagedBuffer(numElements);
            ScalarField *boxVolume = this->resultWriter().allocateManagedBuffer(numDofs);
            ScalarField *enthalpyW = this->resultWriter().allocateManagedBuffer(numDofs);
            ScalarField *enthalpyN = this->resultWriter().allocateManagedBuffer(numDofs);
286
287
288
            (*boxVolume) = 0;

            //Fill the scalar fields with values
289

290
291
            ScalarField *rank = this->resultWriter().allocateManagedBuffer(numElements);

292
            FVElementGeometry fvGeometry;
293
294
            VolumeVariables volVars;

295
296
297
            ElementIterator eIt = this->gridView().template begin<0>();
            ElementIterator eEndIt = this->gridView().template end<0>();
            for (; eIt != eEndIt; ++eIt)
298
            {
299
300
                int eIdx = this->elementMapper().map(*eIt);
                (*rank)[eIdx] = this->gridView().comm().rank();
301
                fvGeometry.update(this->gridView(), *eIt);
302
303


304
                for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx)
305
                {
306
                    int globalIdx = this->model().dofMapper().map(*eIt, scvIdx, dofCodim);
307
308
                    volVars.update(this->model().curSol()[globalIdx],
                                   *this,
309
                                   *eIt,
310
311
                                   fvGeometry,
                                   scvIdx,
312
                                   false);
313
                    (*boxVolume)[globalIdx] += fvGeometry.subContVol[scvIdx].volume;
314
315
316
                    (*enthalpyW)[globalIdx] = volVars.enthalpy(lPhaseIdx);
                    (*enthalpyN)[globalIdx] = volVars.enthalpy(gPhaseIdx);
                }
317
318
                (*Kxx)[eIdx] = this->spatialParams().intrinsicPermeability(*eIt, fvGeometry, /*element data*/ 0);
                (*cellPorosity)[eIdx] = this->spatialParams().porosity(*eIt, fvGeometry, /*element data*/ 0);
319
320
321
            }

            //pass the scalar fields to the vtkwriter
322
323
324
325
326
            this->resultWriter().attachDofData(*boxVolume, "boxVolume", isBox);
            this->resultWriter().attachDofData(*Kxx, "Kxx", false); //element data
            this->resultWriter().attachDofData(*cellPorosity, "cellwisePorosity", false); //element data
            this->resultWriter().attachDofData(*enthalpyW, "enthalpyW", isBox);
            this->resultWriter().attachDofData(*enthalpyN, "enthalpyN", isBox);
327
328
329
330
331
332
333
334
335
336
337

        }

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

338
#if ISOTHERMAL
339
340
341
    /*!
     * \brief Returns the temperature within the domain.
     *
342
343
344
345
     * \param globalPos The position
     *
     * This problem assumes a geothermal gradient with 
     * a surface temperature of 10 degrees Celsius.
346
     */
347
    Scalar temperatureAtPos(const GlobalPosition &globalPos) const
348
349
    {
        return temperature_(globalPos);
350
    }
351
#endif
352

353
    /*!
354
     * \brief Returns the source term
355
     *
356
357
     * \param values Stores the source values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} / (m^\textrm{dim} \cdot s )] \f$
358
     * \param globalPos The global position
359
360
361
362
363
364
     *
     * Depending on whether useMoles is set on true or false, the flux has to be given either in
     * kg/(m^3*s) or mole/(m^3*s) in the input file!!
     *
     * Note that the energy balance is always calculated in terms of specific enthalpies [J/kg]
     * and that the Neumann fluxes have to be specified accordingly.
365
     */
366
367
368
369
370
371
372
373
374
    void sourceAtPos(PrimaryVariables &values,
                const GlobalPosition &globalPos) const
    {
        values = 0;
    }

    /*!
     * \name Boundary conditions
     */
375

376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
    /*!
     * \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
394
     * \param intersection specifies the intersection at which boundary
395
     *           condition is to set
396
     */
397
    void boundaryTypes(BoundaryTypes &values, const Intersection &intersection) const
398
    {
399
        int boundaryId = intersection.boundaryId();
400
401
402
403
404
405
406
407
408
409
410
        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();
    }

411
    /*!
412
413
     * \brief Evaluates the boundary conditions for a Dirichlet
     *        boundary segment
414
     *
415
416
     * \param values Stores the Dirichlet values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} ] \f$
417
418
419
420
421
422
     * \param globalPos The global position
     */
    void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
    {
        initial_(values, globalPos);
    }
423

424
    /*!
425
     * \brief Evaluate the boundary conditions for a Neumann
426
427
     *        boundary segment.
     *
428
429
      * \param values Stores the Neumann values for the conservation equations in
     *               \f$ [ \textnormal{unit of conserved quantity} / (m^(dim-1) \cdot s )] \f$
430
     * \param element The finite element
431
     * \param fvGeometry The finite volume geometry of the element
432
     * \param intersection The intersection between element and boundary
433
     * \param scvIdx The local index of the sub-control volume
434
435
     * \param boundaryFaceIdx The index of the boundary face
     *
436
437
     * The \a values store the mass flux of each phase normal to the boundary.
     * Negative values indicate an inflow.
438
439
440
     *
     * Depending on whether useMoles is set on true or false, the flux has to be given either in
     * kg/(m^2*s) or mole/(m^2*s) in the input file!!
441
442
443
     */
    void neumann(PrimaryVariables &values,
                 const Element &element,
444
                 const FVElementGeometry &fvGeometry,
445
                 const Intersection &intersection,
446
447
448
                 int scvIdx,
                 int boundaryFaceIdx) const
    {
449
        int boundaryId = intersection.boundaryId();
450
451
452
453

        values = 0;
        if (boundaryId == injectionBottom_)
        {
454
            values[contiCO2EqIdx] = -injectionRate_; ///FluidSystem::molarMass(CO2Idx); // kg/(s*m^2) or mole/(m^2*s) !!
455
#if !ISOTHERMAL
456
457
            values[energyEqIdx] = -injectionRate_/*kg/(m^2 s)*/*CO2::gasEnthalpy(
                                    injectionTemperature_, injectionPressure_)/*J/kg*/; // W/(m^2)
458
459
460
461
462
463
464
465
466
467
468
469
#endif
        }
    }

    // \}

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

    /*!
470
     * \brief Evaluates the initial values for a control volume
471
     *
472
473
474
     * \param values Stores the initial values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variables} ] \f$
     * \param globalPos The global position
475
     */
476
477
    void initialAtPos(PrimaryVariables &values,
                      const GlobalPosition &globalPos) const
478
479
480
481
482
    {
        initial_(values, globalPos);
    }

    /*!
483
     * \brief Returns the initial phase state for a control volume.
484
     *
485
     * \param vertex The vertex
486
     * \param globalIdx The global index of the vertex
487
488
     * \param globalPos The global position
     */
489
    int initialPhasePresence(const Vertex &vertex,
490
491
                             int &globalIdx,
                             const GlobalPosition &globalPos) const
492
    { return Indices::wPhaseOnly; }
493
494
495
496

    // \}

private:
497
498
499
500
501
502
503
504
505
    /*!
     * \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
     */
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
    void initial_(PrimaryVariables &values,
                  const GlobalPosition &globalPos) const
    {
        Scalar temp = temperature_(globalPos);
        Scalar densityW = FluidSystem::Brine::liquidDensity(temp, 1e7);

        Scalar pl =  1e5 - densityW*this->gravity()[dim-1]*(depthBOR_ - globalPos[dim-1]);
        Scalar moleFracLiquidCO2 = 0.00;
        Scalar moleFracLiquidBrine = 1.0 - moleFracLiquidCO2;

        Scalar meanM =
            FluidSystem::molarMass(BrineIdx)*moleFracLiquidBrine +
            FluidSystem::molarMass(CO2Idx)*moleFracLiquidCO2;

        Scalar massFracLiquidCO2 = moleFracLiquidCO2*FluidSystem::molarMass(CO2Idx)/meanM;

        values[Indices::pressureIdx] = pl;
        values[Indices::switchIdx] = massFracLiquidCO2;
#if !ISOTHERMAL
            values[temperatureIdx] = temperature_(globalPos); //K
#endif


    }

    Scalar temperature_(const GlobalPosition globalPos) const
    {
533
        Scalar T = 283.0 + (depthBOR_ - globalPos[dim-1])*0.03; 
534
        return T;
535
    }
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560

    Scalar depthBOR_;
    Scalar injectionRate_;
    Scalar injectionPressure_;
    Scalar injectionTemperature_;
    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