fvspatialparams1p.hh 6.82 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
22
23
24
25
 *   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 problems using the
 *        fv method.
 */
26
27
#ifndef DUMUX_FV_SPATIAL_PARAMS_ONE_P_HH
#define DUMUX_FV_SPATIAL_PARAMS_ONE_P_HH
28
29
30
31

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

32
#include <dumux/common/basicproperties.hh>
33

34
#include <dune/common/deprecated.hh>
35
36
37
38
39
40
41
#include <dune/common/fmatrix.hh>

namespace Dumux
{
// forward declation of property tags
namespace Properties
{
42
NEW_PROP_TAG(SpatialParams);
43
44
45
46
47
48
49
50
51
52
53
54
}


/*!
 * \ingroup SpatialParameters
 */

/**
 * \brief The base class for spatial parameters of problems using the
 *        fv method.
 */
template<class TypeTag>
55
class FVSpatialParamsOneP
56
57
58
{
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
59
    typedef typename GET_PROP_TYPE(TypeTag, SpatialParams) Implementation;
60
61
62

    enum
    {
Christoph Grueninger's avatar
Christoph Grueninger committed
63
        dim = GridView::dimension,
64
65
66
67
68
        dimWorld = GridView::dimensionworld
    };

    typedef typename GridView::template Codim<0>::Entity Element;
    typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition;
Natalie Schroeder's avatar
Natalie Schroeder committed
69
    typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimWorldMatrix;
70
71

public:
72
    FVSpatialParamsOneP(const GridView &gridView)
73
74
75
    {
    }

76
    ~FVSpatialParamsOneP()
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
    {
    }

    /*!
     * \brief Averages the intrinsic permeability (Scalar).
     * \param K1 intrinsic permeability of the first element
     * \param K2 intrinsic permeability of the second element
     */
    Scalar meanK(Scalar K1, Scalar K2) const
    {
        const Scalar K = Dumux::harmonicMean(K1, K2);
        return K;
    }

    /*!
     * \brief Averages the intrinsic permeability (Scalar).
     * \param result averaged intrinsic permeability
     * \param K1 intrinsic permeability of the first element
     * \param K2 intrinsic permeability of the second element
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
97
    void meanK(DimWorldMatrix &result, Scalar K1, Scalar K2) const
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
    {
        const Scalar K = Dumux::harmonicMean(K1, K2);
        for (int i = 0; i < dimWorld; ++i)
        {
            for (int j = 0; j < dimWorld; ++j)
                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 element
     * \param K2 intrinsic permeability of the second element
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
114
    void meanK(DimWorldMatrix &result, const DimWorldMatrix &K1, const DimWorldMatrix &K2) const
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
    {
        // entry-wise harmonic mean at the main diagonal and arithmetic mean at the off-diagonal
        for (int i = 0; i < dimWorld; ++i)
        {
            result[i][i] = harmonicMean(K1[i][i], K2[i][i]);
            for (int j = 0; j < dimWorld; ++j)
            {
                if (i != j)
                {
                    result[i][j] = 0.5 * (K1[i][j] + K2[i][j]);
                }
            }
        }
    }

    /*!
     * \brief Dummy function that can be used if only one value exist (boundaries).
     * \param result intrinsic permeability
     * \param K intrinsic permeability of the element
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
135
    void meanK(DimWorldMatrix &result, Scalar K) const
136
137
138
139
140
141
142
143
144
145
146
147
148
149
    {
        for (int i = 0; i < dimWorld; ++i)
        {
            for (int j = 0; j < dimWorld; ++j)
                result[i][j] = 0;
            result[i][i] = K;
        }
    }

    /*!
     * \brief Dummy function that can be used if only one value exist (boundaries).
     * \param result intrinsic permeability
     * \param K intrinsic permeability of the element
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
150
    void meanK(DimWorldMatrix &result, const DimWorldMatrix &K) const
151
152
153
154
155
156
157
158
159
160
    {
        result = K;
    }

    /*!
     * \brief Function for defining the intrinsic (absolute) permeability.
     *
     * \return intrinsic (absolute) permeability
     * \param element The element
     */
Natalie Schroeder's avatar
Natalie Schroeder committed
161
    const DimWorldMatrix& intrinsicPermeability (const Element& element) const
162
163
164
165
166
167
168
169
170
171
    {
        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
172
    const DimWorldMatrix& intrinsicPermeabilityAtPos (const GlobalPosition& globalPos) const
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
    {
        DUNE_THROW(Dune::InvalidStateException,
                   "The spatial parameters do not provide "
                   "a intrinsicPermeabilityAtPos() method.");
    }

    /*!
     * \brief Function for defining the porosity.
     *
     * \return porosity
     * \param element The element
     */
    Scalar porosity(const Element& element) 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.");
    }

203
    double DUNE_DEPRECATED tortuosity(const Element &element) const
204
205
206
207
    {
            return -1;
    }

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

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

} // namespace Dumux

#endif