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 abstract base class for friction laws. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_MATERIAL_FLUIDMATRIX_FRICTIONLAW_NIKURADSE_HH | ||
14 | #define DUMUX_MATERIAL_FLUIDMATRIX_FRICTIONLAW_NIKURADSE_HH | ||
15 | |||
16 | /*! | ||
17 | * \file | ||
18 | * \ingroup FrictionLaws | ||
19 | * \brief Implementation of the friction law after Nikuradse. | ||
20 | */ | ||
21 | #include <algorithm> | ||
22 | #include <cmath> | ||
23 | #include <dune/common/math.hh> | ||
24 | #include "frictionlaw.hh" | ||
25 | |||
26 | namespace Dumux { | ||
27 | /*! | ||
28 | * \addtogroup FrictionLaws | ||
29 | * \copydetails Dumux::FrictionLawNikuradse | ||
30 | */ | ||
31 | |||
32 | /*! | ||
33 | * \ingroup FrictionLaws | ||
34 | * \brief Implementation of the friction law after Nikuradse. | ||
35 | * | ||
36 | * ### Nikuradse | ||
37 | * | ||
38 | * This friction law calculates the stress between the flowing fluid and the bottom, | ||
39 | * which is called bottom shear stress, using the Nikuradse \cite Nikuradse1950 friction law | ||
40 | * | ||
41 | *\f$\tau_{x} = \frac{\kappa^2}{(ln\frac{12h}{ks})^2} u \sqrt{u^2 + v^2}\f$ and | ||
42 | *\f$\tau_{y} = \frac{\kappa^2}{(ln\frac{12h}{ks})^2} v \sqrt{u^2 + v^2}\f$ | ||
43 | * | ||
44 | * with the dimensionless Karman's constant \f$\mathrm{\kappa}\f$, the quivalent sand roughness | ||
45 | * \f$\mathrm{ks}\f$ in \f$\mathrm{[m]}\f$ and the water depth \f$\mathrm{h}\f$ | ||
46 | * in \f$\mathrm{[m]}\f$. | ||
47 | * | ||
48 | * The bottom shear stress is needed to calculate on the one hand the loss of | ||
49 | * momentum due to bottom friction and on the other hand the bedload transport rate. | ||
50 | * | ||
51 | * The LET mobility model is used to limit the friction for small water | ||
52 | * depths if a roughness height > 0.0 is provided (default roughnessHeight = 0.0). | ||
53 | * | ||
54 | */ | ||
55 | |||
56 | template <typename VolumeVariables> | ||
57 | class FrictionLawNikuradse : public FrictionLaw<VolumeVariables> | ||
58 | { | ||
59 | using Scalar = typename VolumeVariables::PrimaryVariables::value_type; | ||
60 | public: | ||
61 | /*! | ||
62 | * \brief Constructor | ||
63 | * | ||
64 | * \param ks Equivalent sand roughness (in m) | ||
65 | * \param roughnessHeight roughness height (in m) default = 0.0 | ||
66 | */ | ||
67 | 1 | FrictionLawNikuradse(const Scalar ks, const Scalar roughnessHeight=0.0) | |
68 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
|
1 | : ks_(ks), roughnessHeight_(roughnessHeight) {} |
69 | |||
70 | /*! | ||
71 | * \brief Compute the bottom shear stress. | ||
72 | * | ||
73 | * \param volVars Volume variables | ||
74 | * | ||
75 | * Compute the bottom shear stress due to bottom friction. | ||
76 | * The bottom shear stress is a projection of the shear stress tensor onto the river bed. | ||
77 | * It can therefore be represented by a (tangent) vector with two entries. | ||
78 | * | ||
79 | * \return shear stress in N/m^2. First entry is the x-component, the second the y-component. | ||
80 | */ | ||
81 | 1071465 | Dune::FieldVector<Scalar, 2> bottomShearStress(const VolumeVariables& volVars) const final | |
82 | { | ||
83 | using Dune::power; | ||
84 | using std::log; | ||
85 | using std::hypot; | ||
86 | |||
87 | 1071465 | Dune::FieldVector<Scalar, 2> shearStress(0.0); | |
88 | |||
89 |
1/2✓ Branch 0 taken 1071465 times.
✗ Branch 1 not taken.
|
1071465 | const Scalar artificialWaterDepth = this->limitRoughH(roughnessHeight_, volVars.waterDepth()); |
90 | 1071465 | const Scalar karmanConstant = 0.41; // Karman's constant is dimensionless | |
91 | 1071465 | const Scalar dimensionlessFactor = power(karmanConstant, 2)/power(log((12*(volVars.waterDepth() + artificialWaterDepth))/ks_), 2); | |
92 | 1071465 | const Scalar uv = hypot(volVars.velocity(0),volVars.velocity(1)); | |
93 | |||
94 | 1071465 | shearStress[0] = dimensionlessFactor * volVars.velocity(0) * uv * volVars.density(); | |
95 | 1071465 | shearStress[1] = dimensionlessFactor * volVars.velocity(1) * uv * volVars.density(); | |
96 | |||
97 | 1071465 | return shearStress; | |
98 | } | ||
99 | |||
100 | private: | ||
101 | Scalar ks_; | ||
102 | Scalar roughnessHeight_; | ||
103 | }; | ||
104 | |||
105 | } // end namespace Dumux | ||
106 | |||
107 | #endif | ||
108 |