3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
sequentialfv1p.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_SEQUENTIAL_FV_SPATIAL_PARAMS_ONE_P_HH
26#define DUMUX_SEQUENTIAL_FV_SPATIAL_PARAMS_ONE_P_HH
27
29#include <dumux/common/math.hh>
30
31#include <dune/common/fmatrix.hh>
32
33namespace Dumux {
34
40template<class TypeTag>
42{
47
48 enum
49 {
50 dim = GridView::dimension,
51 dimWorld = GridView::dimensionworld
52 };
53
54 using Element = typename GridView::template Codim<0>::Entity;
55 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
56 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
57
58public:
59 SequentialFVSpatialParamsOneP(const Problem& problem)
60 {
61 }
67 Scalar meanK(Scalar K1, Scalar K2) const
68 {
69 const Scalar K = harmonicMean(K1, K2);
70 return K;
71 }
72
79 void meanK(DimWorldMatrix &result, Scalar K1, Scalar K2) const
80 {
81 const Scalar K = harmonicMean(K1, K2);
82 for (int i = 0; i < dimWorld; ++i)
83 {
84 for (int j = 0; j < dimWorld; ++j)
85 result[i][j] = 0;
86 result[i][i] = K;
87 }
88 }
89
96 void meanK(DimWorldMatrix &result, const DimWorldMatrix &K1, const DimWorldMatrix &K2) const
97 {
98 // entry-wise harmonic mean at the main diagonal and arithmetic mean at the off-diagonal
99 for (int i = 0; i < dimWorld; ++i)
100 {
101 result[i][i] = harmonicMean(K1[i][i], K2[i][i]);
102 for (int j = 0; j < dimWorld; ++j)
103 {
104 if (i != j)
105 {
106 result[i][j] = 0.5 * (K1[i][j] + K2[i][j]);
107 }
108 }
109 }
110 }
111
117 void meanK(DimWorldMatrix &result, Scalar K) const
118 {
119 for (int i = 0; i < dimWorld; ++i)
120 {
121 for (int j = 0; j < dimWorld; ++j)
122 result[i][j] = 0;
123 result[i][i] = K;
124 }
125 }
126
132 void meanK(DimWorldMatrix &result, const DimWorldMatrix &K) const
133 {
134 result = K;
135 }
136
143 const DimWorldMatrix& intrinsicPermeability (const Element& element) const
144 {
145 return asImp_().intrinsicPermeabilityAtPos(element.geometry().center());
146 }
147
154 const DimWorldMatrix& intrinsicPermeabilityAtPos (const GlobalPosition& globalPos) const
155 {
156 DUNE_THROW(Dune::InvalidStateException,
157 "The spatial parameters do not provide "
158 "a intrinsicPermeabilityAtPos() method.");
159 }
160
167 Scalar porosity(const Element& element) const
168 {
169 return asImp_().porosityAtPos(element.geometry().center());
170 }
171
178 Scalar porosityAtPos(const GlobalPosition& globalPos) const
179 {
180 DUNE_THROW(Dune::InvalidStateException,
181 "The spatial parameters do not provide "
182 "a porosityAtPos() method.");
183 }
184
185private:
186 Implementation &asImp_()
187 {
188 return *static_cast<Implementation*> (this);
189 }
190
191 const Implementation &asImp_() const
192 {
193 return *static_cast<const Implementation*> (this);
194 }
195};
196
197} // namespace Dumux
198
199#endif
Define some often used mathematical functions.
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:69
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
The base class for spatial parameters of problems using the fv method.
Definition: sequentialfv1p.hh:42
Scalar porosity(const Element &element) const
Function for defining the porosity.
Definition: sequentialfv1p.hh:167
Scalar porosityAtPos(const GlobalPosition &globalPos) const
Function for defining the porosity.
Definition: sequentialfv1p.hh:178
SequentialFVSpatialParamsOneP(const Problem &problem)
Definition: sequentialfv1p.hh:59
const DimWorldMatrix & intrinsicPermeability(const Element &element) const
Function for defining the intrinsic (absolute) permeability.
Definition: sequentialfv1p.hh:143
void meanK(DimWorldMatrix &result, const DimWorldMatrix &K1, const DimWorldMatrix &K2) const
Averages the intrinsic permeability (Tensor).
Definition: sequentialfv1p.hh:96
const DimWorldMatrix & intrinsicPermeabilityAtPos(const GlobalPosition &globalPos) const
Function for defining the intrinsic (absolute) permeability.
Definition: sequentialfv1p.hh:154
void meanK(DimWorldMatrix &result, const DimWorldMatrix &K) const
Dummy function that can be used if only one value exist (boundaries).
Definition: sequentialfv1p.hh:132
Scalar meanK(Scalar K1, Scalar K2) const
Averages the intrinsic permeability (Scalar).
Definition: sequentialfv1p.hh:67
void meanK(DimWorldMatrix &result, Scalar K1, Scalar K2) const
Averages the intrinsic permeability (Scalar).
Definition: sequentialfv1p.hh:79
void meanK(DimWorldMatrix &result, Scalar K) const
Dummy function that can be used if only one value exist (boundaries).
Definition: sequentialfv1p.hh:117
Declares all properties used in Dumux.