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-FileCopyrightInfo: 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 Fluidmatrixinteractions | ||
10 | * \copydoc Dumux::FrictionLawManning | ||
11 | */ | ||
12 | #ifndef DUMUX_MATERIAL_FLUIDMATRIX_FRICTIONLAW_MANNING_HH | ||
13 | #define DUMUX_MATERIAL_FLUIDMATRIX_FRICTIONLAW_MANNING_HH | ||
14 | |||
15 | #include <algorithm> | ||
16 | #include <cmath> | ||
17 | #include <dune/common/math.hh> | ||
18 | #include "frictionlaw.hh" | ||
19 | |||
20 | namespace Dumux { | ||
21 | /*! | ||
22 | * \ingroup Fluidmatrixinteractions | ||
23 | * \brief Implementation of the friction law after Manning. | ||
24 | * | ||
25 | * The LET mobility model is used to limit the friction for small water | ||
26 | * depths if a roughness height > 0.0 is provided (default roughnessHeight = 0.0). | ||
27 | */ | ||
28 | |||
29 | template <typename VolumeVariables> | ||
30 | class FrictionLawManning : public FrictionLaw<VolumeVariables> | ||
31 | { | ||
32 | using Scalar = typename VolumeVariables::PrimaryVariables::value_type; | ||
33 | public: | ||
34 | /*! | ||
35 | * \brief Constructor | ||
36 | * | ||
37 | * \param gravity Gravity constant (in m/s^2) | ||
38 | * \param manningN Manning friction coefficient (in s/m^(1/3) | ||
39 | * \param roughnessHeight roughness height for limiting (in m) default = 0.0 | ||
40 | */ | ||
41 | 2 | FrictionLawManning(const Scalar gravity, const Scalar manningN, const Scalar roughnessHeight=0.0) | |
42 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
|
2 | : gravity_(gravity), manningN_(manningN), roughnessHeight_(roughnessHeight) {} |
43 | |||
44 | /*! | ||
45 | * \brief Compute the bottom shear stress. | ||
46 | * | ||
47 | * \param volVars Volume variables | ||
48 | * | ||
49 | * Compute the bottom shear stress due to bottom friction. | ||
50 | * The bottom shear stress is a projection of the shear stress tensor onto the river bed. | ||
51 | * It can therefore be represented by a (tangent) vector with two entries. | ||
52 | * | ||
53 | * \return shear stress in N/m^2. First entry is the x-component, the second the y-component. | ||
54 | */ | ||
55 | 1386896 | Dune::FieldVector<Scalar, 2> bottomShearStress(const VolumeVariables& volVars) const final | |
56 | { | ||
57 | using std::pow; | ||
58 | using Dune::power; | ||
59 | using std::hypot; | ||
60 | |||
61 | 1386896 | Dune::FieldVector<Scalar, 2> shearStress(0.0); | |
62 | |||
63 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 1386896 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1386896 times.
|
2773792 | const Scalar artificialWaterDepth = this->limitRoughH(roughnessHeight_, volVars.waterDepth()); |
64 | // c has units of m^(1/2)/s so c^2 has units of m/s^2 | ||
65 | 2773792 | const Scalar c = pow(volVars.waterDepth() + artificialWaterDepth, 1.0/6.0) * 1.0/(manningN_); | |
66 | 4160688 | const Scalar uv = hypot(volVars.velocity(0), volVars.velocity(1)); | |
67 | 1386896 | const Scalar dimensionlessFactor = gravity_/(c*c); | |
68 | |||
69 | 5547584 | shearStress[0] = dimensionlessFactor * volVars.velocity(0) * uv * volVars.density(); | |
70 | 5547584 | shearStress[1] = dimensionlessFactor * volVars.velocity(1) * uv * volVars.density(); | |
71 | |||
72 | 1386896 | return shearStress; | |
73 | } | ||
74 | |||
75 | private: | ||
76 | Scalar gravity_; | ||
77 | Scalar manningN_; | ||
78 | Scalar roughnessHeight_; | ||
79 | }; | ||
80 | |||
81 | } // end namespace Dumux | ||
82 | |||
83 | #endif | ||
84 |