/home/runner/work/amr-wind/amr-wind/amr-wind/boundary_conditions/wall_models/ShearStressSimple.H Source File

AMR-Wind API: /home/runner/work/amr-wind/amr-wind/amr-wind/boundary_conditions/wall_models/ShearStressSimple.H Source File
AMR-Wind API v0.1.0
CFD solver for wind plant simulations
Loading...
Searching...
No Matches
ShearStressSimple.H
Go to the documentation of this file.
1#ifndef SHEARSTRESSSIMPLE_H
2#define SHEARSTRESSSIMPLE_H
3
7
8namespace amr_wind {
9
11{
13 : utau2(ll.utau_mean * ll.utau_mean), wspd_mean(ll.wspd_mean), m_ll(ll)
14 {}
15
16 AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real get_shear(
17 amrex::Real u,
18 amrex::Real /* wspd */,
19 amrex::Real /*u_dx*/,
20 amrex::Real /*v_dx*/,
21 amrex::Real /*x_c*/,
22 amrex::Real /*unit_nor*/) const
23 {
24 return u / wspd_mean * utau2;
25 };
26
27 amrex::Real utau2;
28 amrex::Real wspd_mean;
30};
32{
33 explicit SimpleShearLogLaw(const amr_wind::LogLaw& ll) : m_ll(ll) {}
34
35 AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real get_shear(
36 amrex::Real u,
37 amrex::Real wspd,
38 amrex::Real /*u_dx*/,
39 amrex::Real /*v_dx*/,
40 amrex::Real /*x_c*/,
41 amrex::Real /*unit_nor*/) const
42 {
43 amrex::Real utau = m_ll.get_utau(wspd);
44 return utau * utau * u / wspd;
45 };
46
48};
49
51{
53 const amr_wind::LogLaw& ll, const amr_wind::MOSD& md)
54 : m_ll(ll), m_md(md)
55 {}
56
57 AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real get_shear(
58 amrex::Real u,
59 amrex::Real wspd,
60 amrex::Real u_dx,
61 amrex::Real v_dx,
62 amrex::Real x_c,
63 amrex::Real unit_nor) const
64 {
65 amrex::Real utau = m_ll.get_utau(wspd);
66 amrex::Real tau_vis = utau * utau * u / wspd;
67
68 amrex::Real tau_wave = m_md.get_dyn_tau(u_dx, v_dx, x_c, unit_nor);
69
70 return tau_vis + tau_wave;
71 };
72
75};
76
77} // namespace amr_wind
78
79#endif /* SHEARSTRESSSIMPLE_H */
Definition BCInterface.cpp:7
Definition LogLaw.H:8
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real get_utau(amrex::Real wspd) const
Definition LogLaw.H:33
Definition MOSD.H:8
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE amrex::Real get_dyn_tau(const amrex::Real u_dx, const amrex::Real v_dx, const amrex::Real xc, const amrex::Real unit_nor) const
Definition MOSD.H:19
Definition ShearStressSimple.H:32
SimpleShearLogLaw(const amr_wind::LogLaw &ll)
Definition ShearStressSimple.H:33
const amr_wind::LogLaw m_ll
Definition ShearStressSimple.H:47
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real get_shear(amrex::Real u, amrex::Real wspd, amrex::Real, amrex::Real, amrex::Real, amrex::Real) const
Definition ShearStressSimple.H:35
Definition ShearStressSimple.H:51
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real get_shear(amrex::Real u, amrex::Real wspd, amrex::Real u_dx, amrex::Real v_dx, amrex::Real x_c, amrex::Real unit_nor) const
Definition ShearStressSimple.H:57
const amr_wind::LogLaw m_ll
Definition ShearStressSimple.H:73
const amr_wind::MOSD m_md
Definition ShearStressSimple.H:74
SimpleShearMOSD(const amr_wind::LogLaw &ll, const amr_wind::MOSD &md)
Definition ShearStressSimple.H:52
Definition ShearStressSimple.H:11
amrex::Real wspd_mean
Definition ShearStressSimple.H:28
const amr_wind::LogLaw m_ll
Definition ShearStressSimple.H:29
SimpleShearSchumann(const amr_wind::LogLaw &ll)
Definition ShearStressSimple.H:12
AMREX_GPU_DEVICE AMREX_FORCE_INLINE amrex::Real get_shear(amrex::Real u, amrex::Real, amrex::Real, amrex::Real, amrex::Real, amrex::Real) const
Definition ShearStressSimple.H:16
amrex::Real utau2
Definition ShearStressSimple.H:27