injectionspatialparams.hh 10.3 KB
Newer Older
1
2
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
Bernd Flemisch's avatar
Bernd Flemisch committed
3
4
5
/*****************************************************************************
 *   Copyright (C) 2008-2009 by Klaus Mosthaf                                *
 *   Copyright (C) 2008-2009 by Andreas Lauser                               *
Andreas Lauser's avatar
Andreas Lauser committed
6
 *   Institute for Modelling Hydraulic and Environmental Systems             *
Bernd Flemisch's avatar
Bernd Flemisch committed
7
8
9
 *   University of Stuttgart, Germany                                        *
 *   email: <givenname>.<name>@iws.uni-stuttgart.de                          *
 *                                                                           *
10
 *   This program is free software: you can redistribute it and/or modify    *
Bernd Flemisch's avatar
Bernd Flemisch committed
11
 *   it under the terms of the GNU General Public License as published by    *
12
13
 *   the Free Software Foundation, either version 2 of the License, or       *
 *   (at your option) any later version.                                     *
Bernd Flemisch's avatar
Bernd Flemisch committed
14
 *                                                                           *
15
16
17
18
19
20
21
 *   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/>.   *
Bernd Flemisch's avatar
Bernd Flemisch committed
22
 *****************************************************************************/
Melanie Darcis's avatar
doc    
Melanie Darcis committed
23
24
25
/*!
 * \file
 *
26
27
 * \brief Definition of the spatial parameters for the injection
 *        problem which uses the isothermal 2p2c box model
Melanie Darcis's avatar
doc    
Melanie Darcis committed
28
29
 */

30
31
#ifndef DUMUX_INJECTION_SPATIAL_PARAMS_HH
#define DUMUX_INJECTION_SPATIAL_PARAMS_HH
Bernd Flemisch's avatar
Bernd Flemisch committed
32
33
34
35
36
37

#include <dumux/material/spatialparameters/boxspatialparameters.hh>
#include <dumux/material/fluidmatrixinteractions/2p/linearmaterial.hh>
#include <dumux/material/fluidmatrixinteractions/2p/regularizedbrookscorey.hh>
#include <dumux/material/fluidmatrixinteractions/2p/efftoabslaw.hh>

Bernd Flemisch's avatar
Bernd Flemisch committed
38
39
#include <dumux/boxmodels/2p2c/2p2cmodel.hh>

Bernd Flemisch's avatar
Bernd Flemisch committed
40
41
42
namespace Dumux
{

43
44
//forward declaration
template<class TypeTag>
45
class InjectionSpatialParams;
46
47
48
49

namespace Properties
{
// The spatial parameters TypeTag
50
NEW_TYPE_TAG(InjectionSpatialParams);
51
52

// Set the spatial parameters
53
SET_TYPE_PROP(InjectionSpatialParams, SpatialParams, Dumux::InjectionSpatialParams<TypeTag>);
54
55

// Set the material Law
56
SET_PROP(InjectionSpatialParams, MaterialLaw)
57
{
58
 private:
59
60
    // define the material law which is parameterized by effective
    // saturations
61
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
62
    typedef RegularizedBrooksCorey<Scalar> EffMaterialLaw;
63
 public:
64
65
66
67
68
    // define the material law parameterized by absolute saturations
    typedef EffToAbsLaw<EffMaterialLaw> type;
};
}

Melanie Darcis's avatar
doc    
Melanie Darcis committed
69
/*!
70
 * \ingroup TwoPTwoCModel
71
 * \ingroup BoxTestProblems
72
73
 * \brief Definition of the spatial parameters for the injection
 *        problem which uses the isothermal 2p2c box model
Bernd Flemisch's avatar
Bernd Flemisch committed
74
75
 */
template<class TypeTag>
76
class InjectionSpatialParams : public BoxSpatialParameters<TypeTag>
Bernd Flemisch's avatar
Bernd Flemisch committed
77
{
78
    typedef BoxSpatialParameters<TypeTag> ParentType;
79
80
81
    typedef typename GET_PROP_TYPE(TypeTag, Grid) Grid;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
82
    typedef typename Grid::ctype CoordScalar;
83
    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
Bernd Flemisch's avatar
Bernd Flemisch committed
84

Bernd Flemisch's avatar
Bernd Flemisch committed
85
86
87
88
    enum {
        dim=GridView::dimension,
        dimWorld=GridView::dimensionworld,

89
        wPhaseIdx = FluidSystem::wPhaseIdx
Bernd Flemisch's avatar
Bernd Flemisch committed
90
91
92
    };

    typedef Dune::FieldVector<CoordScalar,dimWorld> GlobalPosition;
93
    typedef Dune::FieldVector<CoordScalar,dimWorld> DimVector;
Bernd Flemisch's avatar
Bernd Flemisch committed
94

95
96
    typedef typename GET_PROP_TYPE(TypeTag, FluxVariables) FluxVariables;
    typedef typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables) ElementVolumeVariables;
Bernd Flemisch's avatar
Bernd Flemisch committed
97

98
    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
99
    typedef typename GridView::template Codim<0>::Entity Element;
Bernd Flemisch's avatar
Bernd Flemisch committed
100
101

public:
102
    typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
Bernd Flemisch's avatar
Bernd Flemisch committed
103
104
    typedef typename MaterialLaw::Params MaterialLawParams;

Melanie Darcis's avatar
doc    
Melanie Darcis committed
105
106
107
    /*!
     * \brief The constructor
     *
108
     * \param gridView The grid view
Melanie Darcis's avatar
doc    
Melanie Darcis committed
109
     */
110
    InjectionSpatialParams(const GridView &gridView)
111
        : ParentType(gridView)
Bernd Flemisch's avatar
Bernd Flemisch committed
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
    {
        layerBottom_ = 22.0;

        // intrinsic permeabilities
        fineK_ = 1e-13;
        coarseK_ = 1e-12;

        // porosities
        finePorosity_ = 0.3;
        coarsePorosity_ = 0.3;

        // residual saturations
        fineMaterialParams_.setSwr(0.2);
        fineMaterialParams_.setSnr(0.0);
        coarseMaterialParams_.setSwr(0.2);
        coarseMaterialParams_.setSnr(0.0);

        // parameters for the Brooks-Corey law
        fineMaterialParams_.setPe(1e4);
        coarseMaterialParams_.setPe(1e4);
132
133
        fineMaterialParams_.setLambda(2.0);
        coarseMaterialParams_.setLambda(2.0);
Bernd Flemisch's avatar
Bernd Flemisch committed
134
135
    }

136
    ~InjectionSpatialParams()
Bernd Flemisch's avatar
Bernd Flemisch committed
137
138
139
140
141
142
    {}

    /*!
     * \brief Apply the intrinsic permeability tensor to a pressure
     *        potential gradient.
     *
143
     * \param element The current finite element
144
     * \param fvGeometry The current finite volume geometry of the element
Bernd Flemisch's avatar
Bernd Flemisch committed
145
     * \param scvIdx The index of the sub-control volume
Bernd Flemisch's avatar
Bernd Flemisch committed
146
     */
147
    const Scalar intrinsicPermeability(const Element &element,
148
                                       const FVElementGeometry &fvGeometry,
149
                                       const int scvIdx) const
Bernd Flemisch's avatar
Bernd Flemisch committed
150
    {
151
152
        const GlobalPosition &globalPos = fvGeometry.subContVol[scvIdx].global;
        if (isFineMaterial_(globalPos))
Bernd Flemisch's avatar
Bernd Flemisch committed
153
154
155
156
157
            return fineK_;
        return coarseK_;
    }

    /*!
158
     * \brief Define the porosity \f$[-]\f$ of the spatial parameters
Bernd Flemisch's avatar
Bernd Flemisch committed
159
     *
160
     * \param element The finite element
161
     * \param fvGeometry The finite volume geometry
162
     * \param scvIdx The local index of the sub-control volume where
Bernd Flemisch's avatar
Bernd Flemisch committed
163
164
     *                    the porosity needs to be defined
     */
Nicolas Schwenck's avatar
Nicolas Schwenck committed
165
    Scalar porosity(const Element &element,
166
                    const FVElementGeometry &fvGeometry,
167
                    const int scvIdx) const
Bernd Flemisch's avatar
Bernd Flemisch committed
168
    {
169
170
        const GlobalPosition &globalPos = fvGeometry.subContVol[scvIdx].global;
        if (isFineMaterial_(globalPos))
Bernd Flemisch's avatar
Bernd Flemisch committed
171
172
173
174
175
            return finePorosity_;
        return coarsePorosity_;
    }


Melanie Darcis's avatar
doc    
Melanie Darcis committed
176
    /*!
177
     * \brief return the parameter object for the Brooks-Corey material law which depends on the position
Melanie Darcis's avatar
doc    
Melanie Darcis committed
178
     *
179
180
181
182
     * \param element The current finite element
     * \param fvGeometry The current finite volume geometry of the element
     * \param scvIdx The index of the sub-control volume
     */
183
    const MaterialLawParams& materialLawParams(const Element &element,
184
185
                                               const FVElementGeometry &fvGeometry,
                                               const int scvIdx) const
Bernd Flemisch's avatar
Bernd Flemisch committed
186
    {
187
188
        const GlobalPosition &globalPos = fvGeometry.subContVol[scvIdx].global;
        if (isFineMaterial_(globalPos))
Bernd Flemisch's avatar
Bernd Flemisch committed
189
190
191
192
193
194
195
196
197
            return fineMaterialParams_;
        return coarseMaterialParams_;
    }

    /*!
     * \brief Returns the heat capacity \f$[J/m^3 K]\f$ of the rock matrix.
     *
     * This is only required for non-isothermal models.
     *
198
     * \param element The finite element
199
     * \param fvGeometry The finite volume geometry
200
     * \param scvIdx The local index of the sub-control volume where
Bernd Flemisch's avatar
Bernd Flemisch committed
201
202
     *                    the heat capacity needs to be defined
     */
203
    double heatCapacity(const Element &element,
204
                        const FVElementGeometry &fvGeometry,
205
                        const int scvIdx) const
Bernd Flemisch's avatar
Bernd Flemisch committed
206
207
208
209
    {
        return
            790 // specific heat capacity of granite [J / (kg K)]
            * 2700 // density of granite [kg/m^3]
210
            * (1 - porosity(element, fvGeometry, scvIdx));
Bernd Flemisch's avatar
Bernd Flemisch committed
211
212
213
214
215
216
217
218
    }

    /*!
     * \brief Calculate the heat flux \f$[W/m^2]\f$ through the
     *        rock matrix based on the temperature gradient \f$[K / m]\f$
     *
     * This is only required for non-isothermal models.
     *
Melanie Darcis's avatar
doc    
Melanie Darcis committed
219
     * \param heatFlux The resulting heat flux vector
220
221
     * \param fluxVars The flux variables
     * \param elemVolVars The volume variables
222
223
     * \param tempGrad The temperature gradient
     * \param element The current finite element
224
     * \param fvGeometry The finite volume geometry of the current element
225
     * \param scvfIdx The local index of the sub-control volume face where
Bernd Flemisch's avatar
Bernd Flemisch committed
226
227
     *                    the matrix heat flux should be calculated
     */
Klaus Mosthaf's avatar
Klaus Mosthaf committed
228
    void matrixHeatFlux(DimVector &heatFlux,
229
230
                        const FluxVariables &fluxVars,
                        const ElementVolumeVariables &elemVolVars,
Klaus Mosthaf's avatar
Klaus Mosthaf committed
231
                        const DimVector &tempGrad,
232
                        const Element &element,
233
                        const FVElementGeometry &fvGeometry,
234
                        const int faceIdx) const
Bernd Flemisch's avatar
Bernd Flemisch committed
235
236
237
238
239
    {
        static const Scalar lWater = 0.6;
        static const Scalar lGranite = 2.8;

        // arithmetic mean of the liquid saturation and the porosity
240
241
        const int i = fvGeometry.subContVolFace[faceIdx].i;
        const int j = fvGeometry.subContVolFace[faceIdx].j;
242
243
244
245
        Scalar sW = std::max<Scalar>(0.0, (elemVolVars[i].saturation(wPhaseIdx) +
                                           elemVolVars[j].saturation(wPhaseIdx)) / 2);
        Scalar poro = (porosity(element, fvGeometry, i) +
                       porosity(element, fvGeometry, j)) / 2;
Bernd Flemisch's avatar
Bernd Flemisch committed
246
247
248
249
250
251

        Scalar lsat = pow(lGranite, (1-poro)) * pow(lWater, poro);
        Scalar ldry = pow(lGranite, (1-poro));

        // the heat conductivity of the matrix. in general this is a
        // tensorial value, but we assume isotropic heat conductivity.
252
        Scalar heatCond = ldry + sqrt(sW) * (ldry - lsat);
Bernd Flemisch's avatar
Bernd Flemisch committed
253
254
255

        // the matrix heat flux is the negative temperature gradient
        // times the heat conductivity.
256
        heatFlux = tempGrad;
Bernd Flemisch's avatar
Bernd Flemisch committed
257
258
259
260
        heatFlux *= -heatCond;
    }

private:
261
262
    bool isFineMaterial_(const GlobalPosition &globalPos) const
    { return globalPos[dim-1] > layerBottom_; };
Bernd Flemisch's avatar
Bernd Flemisch committed
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277

    Scalar fineK_;
    Scalar coarseK_;
    Scalar layerBottom_;

    Scalar finePorosity_;
    Scalar coarsePorosity_;

    MaterialLawParams fineMaterialParams_;
    MaterialLawParams coarseMaterialParams_;
};

}

#endif