implicitspatialparams1p.hh 7.71 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
24
25
26
27
28
29
30
31
32
33
// -*- 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
 *
 * \ingroup SpatialParameters
 * \brief The base class for spatial parameters of one-phase problems 
 * using a fully implicit discretization method.
 */
#ifndef DUMUX_IMPLICIT_SPATIAL_PARAMS_ONE_P_HH
#define DUMUX_IMPLICIT_SPATIAL_PARAMS_ONE_P_HH

#include <dumux/common/propertysystem.hh>
#include <dumux/common/math.hh>

#include <dumux/implicit/common/implicitproperties.hh>

34
#include <dune/common/deprecated.hh>
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
#include <dune/common/fmatrix.hh>

namespace Dumux {
// forward declaration of property tags
namespace Properties {
NEW_PROP_TAG(SpatialParams);
NEW_PROP_TAG(SpatialParamsForchCoeff);
}

/*!
 * \ingroup SpatialParameters
 */


/**
 * \brief The base class for spatial parameters of one-phase problems 
 * using a fully implicit discretization method.
 */
template<class TypeTag>
class ImplicitSpatialParamsOneP
{
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, SpatialParams) Implementation;

    enum { dimWorld = GridView::dimensionworld };
61
    enum { dim = GridView::dimension};
62
63
64
65
66

    typedef typename GridView::template Codim<0>::Entity Element;
    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;

    typedef typename GridView::ctype CoordScalar;
Natalie Schroeder's avatar
Natalie Schroeder committed
67
    typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimWorldMatrix;
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
    typedef Dune::FieldVector<CoordScalar,dimWorld> GlobalPosition;

public:
    ImplicitSpatialParamsOneP(const GridView &gridView)
    { }

    ~ImplicitSpatialParamsOneP()
    {}

    /*!
     * \brief Averages the intrinsic permeability (Scalar).
     * \param result averaged intrinsic permeability
     * \param K1 intrinsic permeability of the first node
     * \param K2 intrinsic permeability of the second node
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
83
    void meanK(DimWorldMatrix &result,
84
85
86
87
               Scalar K1,
               Scalar K2) const
    {
        const Scalar K = Dumux::harmonicMean(K1, K2);
Natalie Schroeder's avatar
Natalie Schroeder committed
88
89
        for (int i = 0; i < dimWorld; ++i) {
            for (int j = 0; j < dimWorld; ++j)
90
91
92
93
94
95
96
97
98
99
100
                result[i][j] = 0;
            result[i][i] = K;
        }
    }

    /*!
     * \brief Averages the intrinsic permeability (Tensor).
     * \param result averaged intrinsic permeability
     * \param K1 intrinsic permeability of the first node
     * \param K2 intrinsic permeability of the second node
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
101
102
103
    void meanK(DimWorldMatrix &result,
               const DimWorldMatrix &K1,
               const DimWorldMatrix &K2) const
104
105
106
    {
        // entry-wise harmonic mean. this is almost certainly wrong if
        // you have off-main diagonal entries in your permeabilities!
Natalie Schroeder's avatar
Natalie Schroeder committed
107
108
        for (int i = 0; i < dimWorld; ++i)
            for (int j = 0; j < dimWorld; ++j)
109
110
111
112
113
114
115
116
117
118
119
                result[i][j] = harmonicMean(K1[i][j], K2[i][j]);
    }

    /*!
     * \brief Function for defining the intrinsic (absolute) permeability.
     *
     * \param element The current element
     * \param fvGeometry The current finite volume geometry of the element
     * \param scvIdx The index of the sub-control volume.
     * \return the intrinsic permeability
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
120
    const DimWorldMatrix& intrinsicPermeability (const Element &element,
121
122
123
124
125
126
127
128
129
130
131
132
            const FVElementGeometry &fvGeometry,
            int scvIdx) const
    {
        return asImp_().intrinsicPermeabilityAtPos(element.geometry().center());
    }

    /*!
     * \brief Function for defining the intrinsic (absolute) permeability.
     *
     * \return intrinsic (absolute) permeability
     * \param globalPos The position of the center of the element
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
133
    const DimWorldMatrix& intrinsicPermeabilityAtPos (const GlobalPosition& globalPos) const
134
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
168
169
170
171
172
173
174
175
176
177
    {
        DUNE_THROW(Dune::InvalidStateException,
                   "The spatial parameters do not provide "
                   "a intrinsicPermeabilityAtPos() method.");
    }

    /*!
     * \brief Function for defining the porosity.
     *
     * \param element The current element
     * \param fvGeometry The current finite volume geometry of the element
     * \param scvIdx The index of the sub-control volume.
     * \return porosity
     */
    Scalar porosity(const Element &element,
            const FVElementGeometry &fvGeometry,
            int scvIdx) const
    {
        return asImp_().porosityAtPos(element.geometry().center());
    }

    /*!
     * \brief Function for defining the porosity.
     *
     * \return porosity
     * \param globalPos The position of the center of the element
     */
    Scalar porosityAtPos(const GlobalPosition& globalPos) const
    {
        DUNE_THROW(Dune::InvalidStateException,
                   "The spatial parameters do not provide "
                   "a porosityAtPos() method.");
    }

    /*!
     * \brief Apply the Forchheimer coefficient for inertial forces
     *        calculation.
     *
     *        Source: Ward, J.C. 1964 Turbulent flow in porous media. ASCE J. Hydraul. Div 90.
     *        Actually the Forchheimer coefficient is also a function of the dimensions of the
     *        porous medium. Taking it as a constant is only a first approximation
     *        (Nield, Bejan, Convection in porous media, 2006, p. 10)
     *
     * \param element The current finite element
178
     * \param fvGeometry The current finite volume geometry of the element
179
180
181
182
     * \param scvIdx The index sub-control volume face where the
     *                      intrinsic velocity ought to be calculated.
     */
    Scalar forchCoeff(const Element &element,
183
                    const FVElementGeometry &fvGeometry,
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
                    const unsigned int scvIdx) const
    {
        try
        {
            const Scalar forchCoeff = GET_PARAM_FROM_GROUP(TypeTag, Scalar, SpatialParams, ForchCoeff);

            return forchCoeff ;
        }
        catch (Dumux::ParameterException &e) {
            std::cerr << e << ". Aborted in file "<< __FILE__ << "!\n";
            exit(1) ;
        }
        catch (...) {
            std::cerr << "Unknown exception thrown!\n";
            exit(1) ;
        }
    }

202
    double DUNE_DEPRECATED tortuosity(const Element &element,
203
204
205
206
207
208
                    const FVElementGeometry &fvGeometry,
                    const int scvIdx) const
    {
            return -1;
    }

209
210
211
212
213
214
215
216
217
218
219
private:
    Implementation &asImp_()
    { return *static_cast<Implementation*>(this); }

    const Implementation &asImp_() const
    { return *static_cast<const Implementation*>(this); }
};

} // namespace Dumux

#endif