Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
richardson_cpu.f90
Go to the documentation of this file.
1! Copyright (c) 2026, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use num_types, only : rp
36 use utils, only : neko_error, neko_warning
37 use logger, only : log_size, neko_log
38 use math, only : glsum, glmin, glmax
39 implicit none
40 private
41
43
44 abstract interface
45 function tau_interface(magu, Ri_b, h, z0, l, kappa) result(tau)
46 import rp
47 real(kind=rp), intent(in) :: magu, ri_b, h, z0, l, kappa
48 real(kind=rp) :: tau
49 end function tau_interface
50
51 function heat_flux_interface(ti, ts, Ri_b, h, magu, z0h, Pr,&
52 l, utau, kappa) result(heat_flux)
53 import rp
54 real(kind=rp), intent(in) :: ts, ti, ri_b, h, magu
55 real(kind=rp), intent(in) :: z0h, pr, l, utau, kappa
56 real(kind=rp) :: heat_flux
57 end function heat_flux_interface
58
59 end interface
60
61 ! These will point to the correct functions
62 ! depending on stability regime and bc_type.
63 procedure(tau_interface), pointer :: tau_ptr => null()
64 procedure(heat_flux_interface), pointer :: heat_flux_ptr => null()
65
66contains
67
69 subroutine compute_ri_b(bc_type, g_dot_n, hi, ti, ts, magu, kappa, q, Ri_b)
70 character(len=*), intent(in) :: bc_type
71 real(kind=rp), intent(in) :: hi, ti, ts
72 real(kind=rp), intent(in) :: magu, kappa, g_dot_n
73 real(kind=rp), intent(inout) :: q, ri_b
74
75 select case (bc_type)
76 case ("neumann")
77 ri_b = - g_dot_n*hi / ti*q / (magu**3*kappa**2)
78 case ("dirichlet")
79 ri_b = g_dot_n*hi/ti*(ti - ts)/magu**2
80 case default
81 call neko_error("Invalid specified temperature b.c. type ('neumann' or 'dirichlet'?)")
82 end select
83 end subroutine compute_ri_b
84
86 subroutine assign_bc_value(bc_type,bc_value,q,ts,ti,kappa,utau,z0h,hi)
87 character(len=*), intent(in) :: bc_type
88 real(kind=rp), intent(in) :: hi, ti, kappa, utau, z0h,bc_value
89 real(kind=rp), intent(inout) :: q,ts
90
91 select case (bc_type)
92 case ("neumann")
93 ! ts not used
94 q = bc_value
95 case ("dirichlet")
96 ts = bc_value
97 q = kappa*utau*(ts - ti)/log(hi/z0h)
98 case default
99 call neko_error("Invalid specified temperature b.c. type ('neumann' or 'dirichlet'?)")
100 end select
101 end subroutine assign_bc_value
102
104 subroutine set_stability_regime(Ri_b,Ri_threshold)
105 real(kind=rp), intent(in) :: ri_b, ri_threshold
106
107 if (ri_b > ri_threshold) then
110 elseif (ri_b < -ri_threshold) then
113 else
116 end if
117 end subroutine set_stability_regime
118
121 subroutine richardson_compute_cpu(u, v, w, temp, ind_r, ind_s, ind_t, ind_e, &
122 n_x, n_y, n_z, h, tau_x, tau_y, tau_z, n_nodes, lx, nelv, &
123 kappa, mu_w, rho_w, g_vec, Pr, z0, z0h_in, bc_type, bc_value, tstep, &
124 Ri_b_diagn, L_ob_diagn, utau_diagn, magu_diagn, ti_diagn, ts_diagn,&
125 q_diagn, h_x_idx, h_y_idx, h_z_idx)
126 integer, intent(in) :: n_nodes, lx, nelv, tstep
127 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, v, w, temp
128 integer, intent(in), dimension(n_nodes) :: ind_r, ind_s, ind_t, ind_e
129 real(kind=rp), dimension(n_nodes), intent(in) :: n_x, n_y, n_z, h
130 real(kind=rp), intent(in) :: kappa, z0, z0h_in, bc_value, pr
131 real(kind=rp), dimension(3), intent(in) :: g_vec
132 real(kind=rp), dimension(n_nodes), intent(in) :: mu_w, rho_w
133 real(kind=rp) :: g_dot_n
134 character(len=*), intent(in) :: bc_type
135 real(kind=rp), dimension(n_nodes), intent(inout) :: tau_x, tau_y, tau_z
136 integer :: i
137 real(kind=rp) :: ui, vi, wi, hi, rho, mu
138 real(kind=rp) :: normu, z0h
139 real(kind=rp) :: l
140 real(kind=rp), parameter :: tol = 0.001_rp
141 real(kind=rp), parameter :: nr_step = 0.001_rp
142 real(kind=rp), parameter :: ri_threshold = 0.0001_rp
143 character(len=LOG_SIZE) :: log_buf
144 real(kind=rp) :: utau, ri_b, l_ob, magu, q, ti, ts
145 real(kind=rp), dimension(n_nodes), intent(inout) :: ri_b_diagn, l_ob_diagn
146 real(kind=rp), dimension(n_nodes), intent(inout) :: utau_diagn, magu_diagn
147 real(kind=rp), dimension(n_nodes), intent(inout) :: ti_diagn, ts_diagn
148 real(kind=rp), dimension(n_nodes), intent(inout) :: q_diagn
149 integer, dimension(n_nodes), intent(in) :: h_x_idx
150 integer, dimension(n_nodes), intent(in) :: h_y_idx
151 integer, dimension(n_nodes), intent(in) :: h_z_idx
152
153 do i=1, n_nodes
154 ! Sample the variables
155 ui = u(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
156 vi = v(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
157 wi = w(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
158 ti = temp(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
159 hi = h(i)
160 rho = rho_w(i)
161 mu = mu_w(i)
162
163 ! Project on horizontal directions
164 normu = ui * n_x(i) + vi * n_y(i) + wi * n_z(i)
165 ui = ui - normu * n_x(i)
166 vi = vi - normu * n_y(i)
167 wi = wi - normu * n_z(i)
168
169 ! Compute velocity magnitude
170 magu = sqrt(ui**2 + vi**2 + wi**2)
171 magu = max(magu, 1.0e-6_rp)
172 utau = magu*kappa / log(hi/z0)
173
174 ! Compute thermal roughness length from Zilitinkevich, 1995
175 if (z0h_in < 0) then
176 ! z0h_in is interpreted as -C_Zil (Zilitinkevich constant) for z0h
177 z0h = z0 * exp(z0h_in*sqrt((utau*z0)/(mu/rho)))
178 else
179 z0h = z0h_in
180 end if
181
182 ! Get q, ts based on bc_type
183 ! Maybe redundant, but needed to initialise Rib
184 call assign_bc_value(bc_type,bc_value,q,ts,ti,kappa,utau,z0h,hi)
185
186 ! Compute g along the normal (generalisation for hills and similar)
187 g_dot_n = abs(g_vec(1)*n_x(i) + g_vec(2)*n_y(i) + g_vec(3)*n_z(i))
188
189 ! Compute Richardson and set stability accordingly
190 call compute_ri_b(bc_type, g_dot_n, hi, ti, ts, magu, kappa, q, ri_b)
191 call set_stability_regime(ri_b, ri_threshold)
192
193 ! Set length scale
194 l = kappa * hi
195 ! Compute u*
196 utau = sqrt(tau_ptr(magu, ri_b, hi, z0, l, kappa))
197 select case (bc_type)
198 case ("neumann")
199 !!! TEMPORARY: neutral log-law approximation
200 ts = ti - (q * pr * log(hi/z0h)) / (max(utau, 1e-6_rp) * kappa)
201 q = q
202 case ("dirichlet")
203 ! Compute q
204 q = heat_flux_ptr(ti, ts, ri_b, hi, magu, z0h, pr, l, utau, kappa)
205 case default
206 call neko_error("Invalid specified temperature b.c. type ('neumann' or 'dirichlet'?)")
207 end select
208
209 ! Distribute according to the velocity vector and bound magu to avoid 0 division
210 magu = max(magu, 1.0e-6_rp)
211 tau_x(i) = -rho*utau**2 * ui / magu
212 tau_y(i) = -rho*utau**2 * vi / magu
213 tau_z(i) = -rho*utau**2 * wi / magu
214 if (abs(ri_b) <= ri_threshold) then
215 ! Neutral (L_ob undefined)
216 l_ob = 1e10_rp
217 else
218 l_ob = -(ts*utau**3)/(kappa*g_dot_n*q)
219 end if
220
221 ri_b_diagn(i) = ri_b
222 l_ob_diagn(i) = l_ob
223 utau_diagn(i) = utau
224 magu_diagn(i) = magu
225 ti_diagn(i) = ti
226 ts_diagn(i) = temp(ind_r(i)-h_x_idx(i), ind_s(i)-h_y_idx(i), ind_t(i)-h_z_idx(i), ind_e(i))
227 q_diagn(i) = q
228 end do
229
230 end subroutine richardson_compute_cpu
231
234 function tau_stable(magu, Ri_b, h, z0, l, kappa) result(tau)
235 real(kind=rp), intent(in) :: magu, ri_b, h, z0, l, kappa
236 real(kind=rp) :: tau
237
238 tau = magu**2/(log(h/z0)**2) * f_tau_stable(ri_b)/ &
239 f_tau_stable(0.0_rp) * (l/h)**2
240 end function tau_stable
241
242 function heat_flux_stable(ti, ts, Ri_b, h, magu, z0h, Pr,&
243 l, utau, kappa) result(heat_flux)
244 real(kind=rp), intent(in) :: ts, ti, ri_b, h, magu
245 real(kind=rp), intent(in) :: z0h, pr, l, utau, kappa
246 real(kind=rp) :: heat_flux
247
248 heat_flux = (ti - ts)/(log(h/z0h)) * &
249 f_theta_stable(ri_b)/abs(f_theta_stable(0.0_rp)) * &
250 (l/h) * utau/pr
251 end function heat_flux_stable
252
253 function f_tau_stable(Ri_b) result(f_tau)
254 real(kind=rp), intent(in) :: ri_b
255 real(kind=rp) :: f_tau
256
257 f_tau = 0.17 * (0.25 + 0.75 / (1.0 + 4.0*ri_b))
258 end function f_tau_stable
259
260 function f_theta_stable(Ri_b) result(f_theta)
261 real(kind=rp), intent(in) :: ri_b
262 real(kind=rp) :: f_theta
263
264 f_theta = -0.145 / (1.0 + 4.0 * ri_b)
265 end function f_theta_stable
266
269 function tau_convective(magu, Ri_b, h, z0, l, kappa) result(tau)
270 real(kind=rp), intent(in) :: magu, ri_b, h, z0, l, kappa
271 real(kind=rp) :: tau
272 real(kind=rp) :: a, b, c
273
274 a = kappa / log(h/z0)
275 b = 2.0
276 c = 7.4 * a**2 * b * (h/z0)**0.5
277
278 tau = a**2 * magu**2 * f_tau_convective(ri_b, c)
279 end function tau_convective
280
281 function heat_flux_convective(ti, ts, Ri_b, h, magu, z0h, Pr,&
282 l, utau, kappa) result(heat_flux)
283 real(kind=rp), intent(in) :: ts, ti, ri_b, h, magu
284 real(kind=rp), intent(in) :: z0h, pr, l, utau, kappa
285 real(kind=rp) :: heat_flux
286 real(kind=rp) :: a, b, c
287
288 a = kappa / log(h/z0h)
289 b = 2.0
290 c = 5.3 * a**2 * b * (h/z0h)**0.5
291
292 heat_flux = - a**2 / 0.74 * magu * &
293 (ti - ts) * f_theta_convective(ri_b, c)
294
295 end function heat_flux_convective
296
297 function f_tau_convective(Ri_b, c) result(f_tau)
298 real(kind=rp), intent(in) :: ri_b, c
299 real(kind=rp) :: f_tau
300
301 f_tau = 1.0 - 2*ri_b / (1.0 + c * abs(ri_b)**0.5)
302 end function f_tau_convective
303
304 function f_theta_convective(Ri_b, c) result(f_theta)
305 real(kind=rp), intent(in) :: ri_b, c
306 real(kind=rp) :: f_theta
307
308 f_theta = 1.0 - 2*ri_b / (1.0 + c * abs(ri_b)**0.5)
309 end function f_theta_convective
310
312 function tau_neutral(magu, Ri_b, h, z0, l, kappa) result(tau)
313 real(kind=rp), intent(in) :: magu, ri_b, h, z0, l, kappa
314 real(kind=rp) :: tau
315
316 tau = (kappa*magu/log(h/z0))**2
317 end function tau_neutral
318
319 function heat_flux_neutral(ti, ts, Ri_b, h, magu, z0h, Pr,&
320 l, utau, kappa) result(heat_flux)
321 real(kind=rp), intent(in) :: ts, ti, ri_b, h, magu
322 real(kind=rp), intent(in) :: z0h, pr, l, utau, kappa
323 real(kind=rp) :: heat_flux
324
325 heat_flux = kappa*utau * (ti - ts)/log(h/z0h)
326 end function heat_flux_neutral
327
328end module richardson_cpu
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:549
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:566
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:596
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements the CPU kernel for the richardson_t type.
subroutine assign_bc_value(bc_type, bc_value, q, ts, ti, kappa, utau, z0h, hi)
Initialises q when the temperature surface bc is dirichlet.
procedure(heat_flux_interface), pointer heat_flux_ptr
real(kind=rp) function f_theta_convective(ri_b, c)
real(kind=rp) function f_tau_stable(ri_b)
real(kind=rp) function heat_flux_neutral(ti, ts, ri_b, h, magu, z0h, pr, l, utau, kappa)
subroutine, public richardson_compute_cpu(u, v, w, temp, ind_r, ind_s, ind_t, ind_e, n_x, n_y, n_z, h, tau_x, tau_y, tau_z, n_nodes, lx, nelv, kappa, mu_w, rho_w, g_vec, pr, z0, z0h_in, bc_type, bc_value, tstep, ri_b_diagn, l_ob_diagn, utau_diagn, magu_diagn, ti_diagn, ts_diagn, q_diagn, h_x_idx, h_y_idx, h_z_idx)
Main routine to compute the surface stresses based on richardson.
real(kind=rp) function f_theta_stable(ri_b)
real(kind=rp) function tau_convective(magu, ri_b, h, z0, l, kappa)
Similarity laws and corrections for the UNSTABLE (convective) regime: Based on Louis 1979.
subroutine set_stability_regime(ri_b, ri_threshold)
Sets the stability regime based on the Richardson number value (quite arbitrary).
real(kind=rp) function f_tau_convective(ri_b, c)
procedure(tau_interface), pointer tau_ptr
real(kind=rp) function heat_flux_stable(ti, ts, ri_b, h, magu, z0h, pr, l, utau, kappa)
real(kind=rp) function tau_stable(magu, ri_b, h, z0, l, kappa)
Similarity laws and corrections for the STABLE regime: Based on Mauritsen et al. 2007.
real(kind=rp) function heat_flux_convective(ti, ts, ri_b, h, magu, z0h, pr, l, utau, kappa)
real(kind=rp) function tau_neutral(magu, ri_b, h, z0, l, kappa)
Similarity laws and corrections for the NEUTRAL regime:
subroutine compute_ri_b(bc_type, g_dot_n, hi, ti, ts, magu, kappa, q, ri_b)
Computes the Richardson number.
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:392
#define max(a, b)
Definition tensor.cu:40