Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
most_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
42 public :: most_compute_cpu
43
44 abstract interface
45 function slaw_m_interface(z, L_ob, z0) result(slaw)
46 import rp
47 real(kind=rp), intent(in) :: z, l_ob, z0
48 real(kind=rp) :: slaw
49 end function slaw_m_interface
50
51 function slaw_h_interface(z, L_ob, z0h) result(slaw)
52 import rp
53 real(kind=rp), intent(in) :: z, l_ob, z0h
54 real(kind=rp) :: slaw
55 end function slaw_h_interface
56
57 function corr_m_interface(z, L_ob) result(corr)
58 import rp
59 real(kind=rp), intent(in) :: z, l_ob
60 real(kind=rp) :: corr
61 end function corr_m_interface
62
63 function corr_h_interface(z, L_ob) result(corr)
64 import rp
65 real(kind=rp), intent(in) :: z, l_ob
66 real(kind=rp) :: corr
67 end function corr_h_interface
68
69 function f_interface(Ri_b, z, z0, z0h, Pr, L_ob, slaw_m, slaw_h) result(f)
71 real(kind=rp), intent(in) :: ri_b, z, z0, z0h, l_ob, pr
72 real(kind=rp) :: f
73 procedure(slaw_m_interface) :: slaw_m
74 procedure(slaw_h_interface) :: slaw_h
75 end function f_interface
76
77 function dfdl_interface(l_upper, l_lower, z, z0, z0h, Pr, L_ob,&
78 slaw_m, slaw_h, fd_h) result(dfdl)
80 real(kind=rp), intent(in) :: l_upper, l_lower, z, z0, z0h, l_ob, fd_h, pr
81 real(kind=rp) :: dfdl
82 procedure(slaw_m_interface) :: slaw_m
83 procedure(slaw_h_interface) :: slaw_h
84 end function dfdl_interface
85 end interface
86
87
88 ! These will point to the correct functions
89 ! depending on stability regime and bc_type.
90 procedure(slaw_m_interface), pointer :: slaw_m_ptr => null()
91 procedure(slaw_h_interface), pointer :: slaw_h_ptr => null()
92 procedure(corr_m_interface), pointer :: corr_m_ptr => null()
93 procedure(corr_h_interface), pointer :: corr_h_ptr => null()
94 procedure(f_interface), pointer :: f_ptr => null()
95 procedure(dfdl_interface), pointer :: dfdl_ptr => null()
96
97contains
98
101 subroutine select_bc_operators(bc_type, bc_value, q, ts, ti, kappa, &
102 utau, z0h, hi, Pr)
103 character(len=*), intent(in) :: bc_type
104 real(kind=rp), intent(in) :: hi, ti, kappa, utau, z0h, bc_value, pr
105 real(kind=rp), intent(inout) :: q,ts
106 select case (bc_type)
107 case ("neumann")
108 ! ts not used
109 q = bc_value
112 case ("dirichlet")
113 ts = bc_value
114 q = (kappa/pr)*utau*(ts - ti)/log(hi/z0h)
117 case default
118 call neko_error("Invalid specified temperature b.c. type ('neumann' or 'dirichlet'?)")
119 end select
120 end subroutine select_bc_operators
121
123 subroutine compute_ri_b(bc_type, g_dot_n, hi, ti, ts, magu, kappa, q, Pr, Ri_b)
124 character(len=*), intent(in) :: bc_type
125 real(kind=rp), intent(in) :: hi, ti, ts, pr
126 real(kind=rp), intent(in) :: magu, kappa, g_dot_n
127 real(kind=rp), intent(inout) :: q, ri_b
128
129 select case (bc_type)
130 case ("neumann")
131 ri_b = - g_dot_n*hi / ti*q*pr / (magu**3*kappa**2)
132 case ("dirichlet")
133 ri_b = g_dot_n*hi/ti*(ti - ts)/magu**2
134 case default
135 call neko_error("Invalid specified temperature b.c. type ('neumann' or 'dirichlet'?)")
136 end select
137 end subroutine compute_ri_b
138
140 subroutine set_stability_regime(Ri_b, Ri_threshold)
141 real(kind=rp), intent(in) :: ri_b, ri_threshold
142
143 if (ri_b > ri_threshold) then
148 elseif (ri_b < -ri_threshold) then
153 else
156 end if
157 end subroutine set_stability_regime
158
161 subroutine most_compute_cpu(u, v, w, temp, ind_r, ind_s, ind_t, ind_e, &
162 n_x, n_y, n_z, h, tau_x, tau_y, tau_z, n_nodes, lx, nelv, &
163 kappa, mu_w, rho_w, g_vec, Pr, z0, z0h_in, bc_type, bc_value, tstep, &
164 Ri_b_diagn, L_ob_diagn, utau_diagn, magu_diagn, ti_diagn, ts_diagn,&
165 q_diagn, h_x_idx, h_y_idx, h_z_idx)
166 integer, intent(in) :: n_nodes, lx, nelv, tstep
167 real(kind=rp), dimension(lx, lx, lx, nelv), intent(in) :: u, v, w, temp
168 integer, intent(in), dimension(n_nodes) :: ind_r, ind_s, ind_t, ind_e
169 real(kind=rp), dimension(n_nodes), intent(in) :: n_x, n_y, n_z, h
170 real(kind=rp), intent(in) :: kappa, z0, z0h_in, bc_value, pr
171 real(kind=rp), dimension(3), intent(in) :: g_vec
172 real(kind=rp), dimension(n_nodes), intent(in) :: mu_w, rho_w
173 real(kind=rp) :: g_dot_n
174 character(len=*), intent(in) :: bc_type
175 real(kind=rp), dimension(n_nodes), intent(inout) :: tau_x, tau_y, tau_z
176 real(kind=rp) :: ui, vi, wi, hi, rho, mu
177 real(kind=rp) :: normu, z0h
178 real(kind=rp) :: l_upper, l_lower, l_old
179 real(kind=rp) :: f, dfdl, fd_h, l_new, l_sign
180 integer :: i, count
181 integer, parameter :: max_count = 50
182 real(kind=rp), parameter :: tol = 0.001_rp
183 real(kind=rp), parameter :: nr_step = 0.001_rp
184 real(kind=rp), parameter :: ri_threshold = 0.0001_rp
185 character(len=LOG_SIZE) :: log_buf
186 real(kind=rp) :: utau, ri_b, l_ob, magu, q, ti, ts
187 real(kind=rp), dimension(n_nodes), intent(inout) :: ri_b_diagn, l_ob_diagn
188 real(kind=rp), dimension(n_nodes), intent(inout) :: utau_diagn, magu_diagn
189 real(kind=rp), dimension(n_nodes), intent(inout) :: ti_diagn, ts_diagn
190 real(kind=rp), dimension(n_nodes), intent(inout) :: q_diagn
191 integer, dimension(n_nodes), intent(in) :: h_x_idx
192 integer, dimension(n_nodes), intent(in) :: h_y_idx
193 integer, dimension(n_nodes), intent(in) :: h_z_idx
194
195 do i=1, n_nodes
196 ! Sample the variables
197 ui = u(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
198 vi = v(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
199 wi = w(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
200 ti = temp(ind_r(i), ind_s(i), ind_t(i), ind_e(i))
201 hi = h(i)
202 rho = rho_w(i)
203 mu = mu_w(i)
204
205 ! Project on horizontal directions
206 normu = ui * n_x(i) + vi * n_y(i) + wi * n_z(i)
207 ui = ui - normu * n_x(i)
208 vi = vi - normu * n_y(i)
209 wi = wi - normu * n_z(i)
210
211 ! Compute velocity magnitude
212 magu = sqrt(ui**2 + vi**2 + wi**2)
213 magu = max(magu, 1.0e-6_rp)
214 utau = magu*kappa / log(hi/z0)
215
216 ! Compute thermal roughness length from Zilitinkevich, 1995
217 if (z0h_in < 0) then
218 ! z0h_in is interpreted as -C_Zil (Zilitinkevich constant) for z0h
219 z0h = z0 * exp(z0h_in*sqrt((utau*z0)/(mu/rho)))
220 else
221 z0h = z0h_in
222 end if
223
224 ! Get q, Ri_b, f_ptr, dfdl_ptr based on bc_type
225 ! Maybe redundant, but needed to initialise Rib
226 call select_bc_operators(bc_type, bc_value, q, ts, ti, kappa, utau, z0h, hi, pr)
227
228 ! Compute g along the normal (generalisation for hills and similar)
229 g_dot_n = abs(g_vec(1)*n_x(i) + g_vec(2)*n_y(i) + g_vec(3)*n_z(i))
230
231 ! Compute Richardson and set stability accordingly
232 call compute_ri_b(bc_type, g_dot_n, hi, ti, ts, magu, kappa, q, pr, ri_b)
233 call set_stability_regime(ri_b, ri_threshold)
234
235 if (abs((ri_b)) <= ri_threshold) then
236 ! Neutral (L_ob undefined)
237 l_ob = 0.0_rp
238 else
239 ! Determine target regime sign
240 if ((ri_b) > 0.0_rp) then
241 l_ob = hi / max(ri_b, ri_threshold) ! Stable guess
242 l_sign = 1.0_rp
243 else
244 l_ob = hi / min(ri_b, -ri_threshold) ! Convective guess
245 l_sign = -1.0_rp
246 end if
247
248 l_old = 1.0e10_rp
249 count = 0
250
251 ! Find Obukhov length
252 do while ((abs(l_old - l_ob)/abs(l_ob) > tol) .and. (count < max_count))
253 ! Switch between stable and convective based on bulk Richardson (Ri_b)
254 l_old = l_ob
255 count = count + 1
256 fd_h = nr_step*l_ob
257 l_upper = l_ob + fd_h
258 l_lower = l_ob - fd_h
259
260 ! Compute L_ob based on stability and bc_type
261 f = f_ptr(ri_b, hi, z0, z0h, pr, l_ob, slaw_m_ptr, slaw_h_ptr)
262 dfdl = dfdl_ptr(l_upper, l_lower, hi, z0, z0h, pr, l_ob, &
263 slaw_m_ptr, slaw_h_ptr, fd_h)
264 if (abs(dfdl) < 1.0e-12_rp) call neko_error("Division by zero in dfdl")
265 l_new = l_ob - f/dfdl
266 ! Avoid regime crossing during Newton iter (otherwise crash)
267 if (l_new*l_sign <= 0.0_rp) then
268 ! "damp update" (stay on same side)
269 l_new = 0.5_rp * l_ob
270 end if
271 ! Bound L_ob
272 l_ob = sign(max(abs(l_new), 1.0e-8_rp), l_sign)
273 l_ob = sign(min(abs(l_ob), 1.0e8_rp), l_sign)
274 end do
275
276 if (abs(l_ob) > 5e5_rp .or. abs(l_ob) < 1e-6_rp) then
277 count = max_count
278 call neko_warning("Obukhov length did not converge (MOST wall model)")
279 end if
280
281 if (.not. associated(f_ptr) .or. .not. associated(dfdl_ptr)) then
282 call neko_error("Unassociated pointer for f or dfdl")
283 end if
284 end if
285
286 ! Based on stability and bc_type, compute utau/q
287 select case (bc_type)
288 case ("neumann")
289 ! Compute u* with the new Obukhov length
290 utau = kappa*magu/slaw_m_ptr(hi, l_ob, z0)
291 case ("dirichlet")
292 ! Compute u* with the new Obukhov length
293 utau = kappa*magu/slaw_m_ptr(hi, l_ob, z0)
294 ! and compute q from here
295 q = kappa/pr*utau*(ts - ti)/slaw_h_ptr(hi, l_ob, z0h)
296 case default
297 call neko_error("Invalid specified temperature b.c. type ('neumann' or 'dirichlet'?)")
298 end select
299
300 ! Distribute according to the velocity vector and bound magu to avoid 0 division
301 magu = max(magu, 1.0e-6_rp)
302 tau_x(i) = -rho * utau**2 * ui / magu
303 tau_y(i) = -rho * utau**2 * vi / magu
304 tau_z(i) = -rho * utau**2 * wi / magu
305
306 ri_b_diagn(i) = ri_b
307 l_ob_diagn(i) = l_ob
308 utau_diagn(i) = utau
309 magu_diagn(i) = magu
310 ti_diagn(i) = ti
311 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))
312 q_diagn(i) = q
313 end do
314 end subroutine most_compute_cpu
315
322 function slaw_m_stable(z, L_ob, z0) result(slaw)
323 real(kind=rp), intent(in) :: z, l_ob, z0
324 real(kind=rp) :: slaw
325
326 slaw = log(z/z0)-corr_m_stable(z,l_ob)+corr_m_stable(z0,l_ob)
327 end function slaw_m_stable
328
329 function slaw_h_stable(z,L_ob,z0h) result(slaw)
330 real(kind=rp), intent(in) :: z,l_ob,z0h
331 real(kind=rp) :: slaw
332
333 slaw = log(z/z0h)-corr_h_stable(z,l_ob)+corr_h_stable(z0h,l_ob)
334 end function slaw_h_stable
335
336 function corr_m_stable(z, L_ob) result(corr)
337 real(kind=rp), intent(in) :: z, l_ob
338 real(kind=rp) :: corr
339 real(kind=rp) :: a, b, c, d
340 real(kind=rp) :: zeta
341 zeta = z/l_ob
342 ! Coefficients specific to Cheng & Brutsaert (2005)
343 a = 1.0_rp
344 b = 2.0_rp/3.0_rp
345 c = 5.0_rp
346 d = 0.35_rp
347 corr = - a*zeta - b*(zeta-c/d)*exp(-d*zeta) - b*c/d
348 end function corr_m_stable
349
350 function corr_h_stable(z, L_ob) result(corr)
351 real(kind=rp), intent(in) :: z, l_ob
352 real(kind=rp) :: corr
353 real(kind=rp) :: a, b, c, d
354 real(kind=rp) :: zeta
355
356 zeta = z/l_ob
357 ! Coefficients specific to Cheng & Brutsaert (2005)
358 a = 1.0_rp
359 b = 2.0_rp/3.0_rp
360 c = 5.0_rp
361 d = 0.35_rp
362 corr = -b * (zeta-c/d)*exp(-d*zeta)-(1.0_rp+ 2.0_rp/3.0_rp * a * zeta)**1.5_rp-b*c/d + 1.0_rp
363 end function corr_h_stable
364
370 function slaw_m_convective(z, L_ob, z0) result(slaw)
371 real(kind=rp), intent(in) :: z, l_ob, z0
372 real(kind=rp) :: slaw
373
374 slaw = log(z/z0) - corr_m_convective(z, l_ob) + corr_m_convective(z0, l_ob)
375 end function slaw_m_convective
376
377 function slaw_h_convective(z, L_ob, z0h) result(slaw)
378 real(kind=rp), intent(in) :: z, l_ob, z0h
379 real(kind=rp) :: slaw
380
381 slaw = log(z/z0h) - corr_h_convective(z, l_ob) + corr_h_convective(z0h, l_ob)
382 end function slaw_h_convective
383
384 function corr_m_convective(z, L_ob) result(corr)
385 real(kind=rp), intent(in) :: z, l_ob
386 real(kind=rp) :: xi, pi, zeta
387 real(kind=rp) :: corr
388
389 zeta = z/l_ob
390 pi = 4*atan(1.0_rp)
391 ! Standard Dyer-Businger coefficient gamma = 16.0
392 xi = (1.0_rp - 16.0_rp*zeta)**0.25_rp
393 corr = 2*log(0.5_rp*(1 + xi)) + log(0.5_rp*(1 + xi**2)) - 2*atan(xi) + pi/2
394 end function corr_m_convective
395
396 function corr_h_convective(z, L_ob) result(corr)
397 real(kind=rp), intent(in) :: z, l_ob
398 real(kind=rp) :: zeta, pi, xi
399 real(kind=rp) :: corr
400
401 zeta = z/l_ob
402 pi = 4*atan(1.0_rp)
403 ! Standard Dyer-Businger coefficient gamma = 16.0
404 xi = (1.0_rp - 16.0_rp*zeta)**0.25_rp
405 corr = 2*log(0.5_rp*(1 + xi**2))
406 end function corr_h_convective
407
409 function slaw_m_neutral(z, L_ob, z0) result(slaw)
410 real(kind=rp), intent(in) :: z, l_ob, z0
411 real(kind=rp) :: slaw
412
413 slaw = log(z/z0)
414 end function slaw_m_neutral
415
416 function slaw_h_neutral(z, L_ob, z0h) result(slaw)
417 real(kind=rp), intent(in) :: z, l_ob, z0h
418 real(kind=rp) :: slaw
419
420 slaw = log(z/z0h)
421 end function slaw_h_neutral
422
424 function f_neumann(Ri_b, z, z0, z0h, Pr, L_ob, slaw_m, slaw_h) result(f)
425 real(kind=rp), intent(in) :: ri_b, z, z0, z0h, l_ob, pr
426 procedure(slaw_m_interface) :: slaw_m
427 procedure(slaw_h_interface) :: slaw_h
428 real(kind=rp) :: f
429
430 f = (ri_b - pr*z/l_ob/slaw_m(z, l_ob, z0)**3)
431 end function f_neumann
432
433 function dfdl_neumann(l_upper, l_lower, z, z0, z0h, Pr, L_ob, &
434 slaw_m, slaw_h, fd_h) result(dfdl)
435 real(kind=rp), intent(in) :: l_upper, l_lower, z, z0, z0h, l_ob, fd_h, pr
436 procedure(slaw_m_interface) :: slaw_m
437 procedure(slaw_h_interface) :: slaw_h
438 real(kind=rp) :: dfdl
439
440 dfdl = (-z/l_upper/slaw_m(z, l_upper, z0)**3) ! conv
441 dfdl = dfdl + (z/l_lower/slaw_m(z, l_lower, z0)**3) ! conv
442 dfdl = pr*dfdl/(2*fd_h)
443 end function dfdl_neumann
444
445 function f_dirichlet(Ri_b, z, z0, z0h, Pr, L_ob, slaw_m, slaw_h) result(f)
446 real(kind=rp), intent(in) :: ri_b, z, z0, z0h, l_ob, pr
447 procedure(slaw_m_interface) :: slaw_m
448 procedure(slaw_h_interface) :: slaw_h
449 real(kind=rp) :: f
450
451 f = (ri_b - pr*z/l_ob*slaw_h(z, l_ob, z0h)/slaw_m(z, l_ob, z0)**2) ! conv
452 end function f_dirichlet
453
454 function dfdl_dirichlet(l_upper, l_lower, z, z0, z0h, Pr, L_ob, &
455 slaw_m, slaw_h, fd_h) result(dfdl)
456 real(kind=rp), intent(in) :: l_upper, l_lower, z, z0, z0h, l_ob, fd_h, pr
457 procedure(slaw_m_interface) :: slaw_m
458 procedure(slaw_h_interface) :: slaw_h
459 real(kind=rp) :: dfdl
460
461 dfdl = (-z/l_upper*slaw_h(z, l_upper, z0h)/slaw_m(z, l_upper, z0)**2) ! conv
462 dfdl = dfdl + (z/l_lower*slaw_h(z, l_lower, z0h)/slaw_m(z, l_lower, z0)**2) ! conv
463 dfdl = pr*dfdl/(2*fd_h)
464 end function dfdl_dirichlet
465
466
467end module most_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), parameter, public pi
Definition math.f90:76
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
Implements the CPU kernel for the most_t type.
Definition most_cpu.f90:34
real(kind=rp) function slaw_m_stable(z, l_ob, z0)
Similarity laws and corrections for the STABLE regime: REFERENCE: Holtslag, A. A. M....
Definition most_cpu.f90:323
subroutine, public most_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 MOST.
Definition most_cpu.f90:166
real(kind=rp) function corr_m_convective(z, l_ob)
Definition most_cpu.f90:385
real(kind=rp) function corr_h_convective(z, l_ob)
Definition most_cpu.f90:397
real(kind=rp) function f_neumann(ri_b, z, z0, z0h, pr, l_ob, slaw_m, slaw_h)
Simialrity laws (different for neumann and dirichlet bc's)
Definition most_cpu.f90:425
real(kind=rp) function slaw_m_neutral(z, l_ob, z0)
Similarity laws and corrections for the NEUTRAL regime:
Definition most_cpu.f90:410
procedure(slaw_m_interface), pointer slaw_m_ptr
Definition most_cpu.f90:90
subroutine compute_ri_b(bc_type, g_dot_n, hi, ti, ts, magu, kappa, q, pr, ri_b)
Computes the Richardson number.
Definition most_cpu.f90:124
real(kind=rp) function corr_m_stable(z, l_ob)
Definition most_cpu.f90:337
subroutine set_stability_regime(ri_b, ri_threshold)
Sets the stability regime based on the Richardson number value (quite arbitrary).
Definition most_cpu.f90:141
procedure(corr_h_interface), pointer corr_h_ptr
Definition most_cpu.f90:93
procedure(dfdl_interface), pointer dfdl_ptr
Definition most_cpu.f90:95
real(kind=rp) function slaw_h_neutral(z, l_ob, z0h)
Definition most_cpu.f90:417
procedure(slaw_h_interface), pointer slaw_h_ptr
Definition most_cpu.f90:91
real(kind=rp) function dfdl_dirichlet(l_upper, l_lower, z, z0, z0h, pr, l_ob, slaw_m, slaw_h, fd_h)
Definition most_cpu.f90:456
real(kind=rp) function slaw_h_stable(z, l_ob, z0h)
Definition most_cpu.f90:330
procedure(corr_m_interface), pointer corr_m_ptr
Definition most_cpu.f90:92
real(kind=rp) function slaw_h_convective(z, l_ob, z0h)
Definition most_cpu.f90:378
real(kind=rp) function corr_h_stable(z, l_ob)
Definition most_cpu.f90:351
real(kind=rp) function f_dirichlet(ri_b, z, z0, z0h, pr, l_ob, slaw_m, slaw_h)
Definition most_cpu.f90:446
procedure(f_interface), pointer f_ptr
Definition most_cpu.f90:94
real(kind=rp) function slaw_m_convective(z, l_ob, z0)
Similarity laws and corrections for the UNSTABLE (convective) regime: REFERENCE: Dyer,...
Definition most_cpu.f90:371
subroutine select_bc_operators(bc_type, bc_value, q, ts, ti, kappa, utau, z0h, hi, pr)
Selects different expressions for the similarity functions in MOST based on the type of bottom bounda...
Definition most_cpu.f90:103
real(kind=rp) function dfdl_neumann(l_upper, l_lower, z, z0, z0h, pr, l_ob, slaw_m, slaw_h, fd_h)
Definition most_cpu.f90:435
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
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