heterogeneousproblemni.hh 20.8 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/implicit/co2/co2model.hh>
37
#include <dumux/material/fluidsystems/brineco2fluidsystem.hh>
38
#include <dumux/implicit/common/implicitporousmediaproblem.hh>
39
#include <dumux/implicit/box/intersectiontovertexbc.hh>
40
#include <test/implicit/co2/heterogeneousspatialparameters.hh>
41
42
43
44
45
46
47
48
49

#include "heterogeneousco2tables.hh"

#define ISOTHERMAL 0

namespace Dumux
{

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

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


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

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

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

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

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

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


/*!
 * \ingroup CO2NIModel
94
 * \ingroup ImplicitTestProblems
95
 * \brief Definition of a problem, where CO2 is injected in a reservoir.
96
97
98
99
100
101
102
103
104
105
106
107
108
 *
 * 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.
 *
109
110
111
 * 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.
 *
112
113
 * 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>
114
 */
115
template <class TypeTag >
116
class HeterogeneousNIProblem : public ImplicitPorousMediaProblem<TypeTag>
117
{
118
    typedef ImplicitPorousMediaProblem<TypeTag> ParentType;
119
120
121
122
123
124
125
126

    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;

127
128
    static const bool useMoles = GET_PROP_VALUE(TypeTag, UseMoles);

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

    // copy some indices for convenience
136
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
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
168
169
    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;
170
171
    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };
    enum { dofCodim = isBox ? dim : 0 };
172
173
174
175
176
177
178
179

public:
    /*!
     * \brief The constructor
     *
     * \param timeManager The time manager
     * \param gridView The grid view
     */
180
    HeterogeneousNIProblem(TimeManager &timeManager,
181
                     const GridView &gridView)
182
        : ParentType(timeManager, GridCreator::grid().leafGridView()),
183
184
185
186
187
188
          //Boundary Id Setup:
          injectionTop_ (1),
          injectionBottom_(2),
          dirichletBoundary_(3),
          noFlowBoundary_(4),
          intersectionToVertexBC_(*this)
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
    {
        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_);
238
239
240
241
242
243
244
245
246
247

        //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;
        }
248
249
250
    }

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

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

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

            //create required scalar fields
283
284
285
286
287
            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);
288
289
290
            (*boxVolume) = 0;

            //Fill the scalar fields with values
291

292
293
            ScalarField *rank = this->resultWriter().allocateManagedBuffer(numElements);

294
            FVElementGeometry fvGeometry;
295
296
            VolumeVariables volVars;

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


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

            //pass the scalar fields to the vtkwriter
324
325
326
327
328
            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);
329
330
331
332
333
334
335
336
337
338
339

        }

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

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

355
    /*!
356
     * \brief Returns the source term
357
     *
358
359
     * \param values Stores the source values for the conservation equations in
     *               \f$ [ \textnormal{unit of primary variable} / (m^\textrm{dim} \cdot s )] \f$
360
     * \param globalPos The global position
361
362
363
364
365
366
     *
     * 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.
367
     */
368
369
370
371
372
373
374
375
376
    void sourceAtPos(PrimaryVariables &values,
                const GlobalPosition &globalPos) const
    {
        values = 0;
    }

    /*!
     * \name Boundary conditions
     */
377

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

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

426
    /*!
427
     * \brief Evaluate the boundary conditions for a Neumann
428
429
     *        boundary segment.
     *
430
431
      * \param values Stores the Neumann values for the conservation equations in
     *               \f$ [ \textnormal{unit of conserved quantity} / (m^(dim-1) \cdot s )] \f$
432
     * \param element The finite element
433
     * \param fvGeometry The finite volume geometry of the element
434
     * \param intersection The intersection between element and boundary
435
     * \param scvIdx The local index of the sub-control volume
436
437
     * \param boundaryFaceIdx The index of the boundary face
     *
438
439
     * The \a values store the mass flux of each phase normal to the boundary.
     * Negative values indicate an inflow.
440
441
442
     *
     * 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!!
443
444
445
     */
    void neumann(PrimaryVariables &values,
                 const Element &element,
446
                 const FVElementGeometry &fvGeometry,
447
                 const Intersection &intersection,
448
449
450
                 int scvIdx,
                 int boundaryFaceIdx) const
    {
451
        int boundaryId = intersection.boundaryId();
452
453
454
455

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

    // \}

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

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

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

    // \}

private:
499
500
501
502
503
504
505
506
507
    /*!
     * \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
     */
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
533
534
    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
    {
535
        Scalar T = 283.0 + (depthBOR_ - globalPos[dim-1])*0.03; 
536
        return T;
537
    }
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562

    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