Line | Branch | Exec | Source |
---|---|---|---|
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 | // SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup FrictionLaws | ||
10 | * \brief Implementation of the friction law after Manning. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_MATERIAL_FLUIDMATRIX_FRICTIONLAW_MANNING_HH | ||
14 | #define DUMUX_MATERIAL_FLUIDMATRIX_FRICTIONLAW_MANNING_HH | ||
15 | |||
16 | #include <algorithm> | ||
17 | #include <cmath> | ||
18 | #include <dune/common/math.hh> | ||
19 | #include "frictionlaw.hh" | ||
20 | |||
21 | namespace Dumux { | ||
22 | /*! | ||
23 | * \addtogroup FrictionLaws | ||
24 | * \copydetails Dumux::FrictionLawManning | ||
25 | */ | ||
26 | |||
27 | /*! | ||
28 | * \ingroup FrictionLaws | ||
29 | * \brief Implementation of the friction law after Manning. | ||
30 | * | ||
31 | * ### Manning | ||
32 | * | ||
33 | * This friction law calculates the stress between the flowing fluid and the bottom, | ||
34 | * which is called bottom shear stress, using the Manning friction law: | ||
35 | * | ||
36 | * \f$\tau_{x} = \frac{g}{(\frac{h^{1/6}}{n})^2} u \sqrt{u^2 + v^2}\f$ and | ||
37 | * \f$\tau_{y} = \frac{g}{(\frac{h^{1/6}}{n})^2} v \sqrt{u^2 + v^2}\f$ | ||
38 | * | ||
39 | * with the gravity constant \f$\mathrm{g}\f$ in \f$\mathrm{[m/s^2]}\f$, the water depth | ||
40 | * \f$\mathrm{h}\f$ in \f$\mathrm{[m]}\f$ and the Manning friction coefficient | ||
41 | * \f$\mathrm{n}\f$ in \f$\mathrm{[s/m^{1/3}]}\f$. | ||
42 | * | ||
43 | * The bottom shear stress is needed to calculate on the one hand the loss of | ||
44 | * momentum due to bottom friction and on the other hand the bedload transport rate. | ||
45 | * | ||
46 | * The LET mobility model is used to limit the friction for small water | ||
47 | * depths if a roughness height > 0.0 is provided (default roughnessHeight = 0.0). | ||
48 | */ | ||
49 | |||
50 | template <typename VolumeVariables> | ||
51 | class FrictionLawManning : public FrictionLaw<VolumeVariables> | ||
52 | { | ||
53 | using Scalar = typename VolumeVariables::PrimaryVariables::value_type; | ||
54 | public: | ||
55 | /*! | ||
56 | * \brief Constructor | ||
57 | * | ||
58 | * \param gravity Gravity constant (in m/s^2) | ||
59 | * \param manningN Manning friction coefficient (in s/m^(1/3) | ||
60 | * \param roughnessHeight roughness height for limiting (in m) default = 0.0 | ||
61 | */ | ||
62 | 2 | FrictionLawManning(const Scalar gravity, const Scalar manningN, const Scalar roughnessHeight=0.0) | |
63 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2 | : gravity_(gravity), manningN_(manningN), roughnessHeight_(roughnessHeight) {} |
64 | |||
65 | /*! | ||
66 | * \brief Compute the bottom shear stress. | ||
67 | * | ||
68 | * \param volVars Volume variables | ||
69 | * | ||
70 | * Compute the bottom shear stress due to bottom friction. | ||
71 | * The bottom shear stress is a projection of the shear stress tensor onto the river bed. | ||
72 | * It can therefore be represented by a (tangent) vector with two entries. | ||
73 | * | ||
74 | * \return shear stress in N/m^2. First entry is the x-component, the second the y-component. | ||
75 | */ | ||
76 | 1378282 | Dune::FieldVector<Scalar, 2> bottomShearStress(const VolumeVariables& volVars) const final | |
77 | { | ||
78 | using std::pow; | ||
79 | using Dune::power; | ||
80 | using std::hypot; | ||
81 | |||
82 | 1378282 | Dune::FieldVector<Scalar, 2> shearStress(0.0); | |
83 | |||
84 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1378282 times.
|
1378282 | const Scalar artificialWaterDepth = this->limitRoughH(roughnessHeight_, volVars.waterDepth()); |
85 | // c has units of m^(1/2)/s so c^2 has units of m/s^2 | ||
86 | 1378282 | const Scalar c = pow(volVars.waterDepth() + artificialWaterDepth, 1.0/6.0) * 1.0/(manningN_); | |
87 | 1378282 | const Scalar uv = hypot(volVars.velocity(0), volVars.velocity(1)); | |
88 | 1378282 | const Scalar dimensionlessFactor = gravity_/(c*c); | |
89 | |||
90 | 1378282 | shearStress[0] = dimensionlessFactor * volVars.velocity(0) * uv * volVars.density(); | |
91 | 1378282 | shearStress[1] = dimensionlessFactor * volVars.velocity(1) * uv * volVars.density(); | |
92 | |||
93 | 1378282 | return shearStress; | |
94 | } | ||
95 | |||
96 | private: | ||
97 | Scalar gravity_; | ||
98 | Scalar manningN_; | ||
99 | Scalar roughnessHeight_; | ||
100 | }; | ||
101 | |||
102 | } // end namespace Dumux | ||
103 | |||
104 | #endif | ||
105 |