Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
richardson.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!
33!
36 use field, only : field_t
37 use num_types, only : rp
38 use json_module, only : json_file
39 use coefs, only : coef_t
41 use wall_model, only : wall_model_t
49 use logger, only : log_size, neko_log
50 use vector, only : vector_t
52 use math, only: masked_gather_copy_0
54 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
56 implicit none
57 private
58
62 type, public, extends(wall_model_t) :: richardson_t
64 real(kind=rp) :: kappa
66 real(kind=rp) :: pr
68 real(kind=rp) :: z0
70 real(kind=rp) :: z0h_in
72 real(kind=rp) :: g(3)
73 ! The dynamic viscosity at the boundary
74 type(vector_t) :: mu_w
75 ! The fluid density at the boundary
76 type(vector_t) :: rho_w
78 character(len=:), allocatable :: bc_type
80 real(kind=rp) :: bc_value
82 character(len=:), allocatable :: scalar_name
84 type(vector_t) :: ri_b, l_ob, utau, magu, ti, ts, q
85 integer, allocatable :: h_x_idx(:), h_y_idx(:), h_z_idx(:)
86 type(c_ptr) :: h_x_idx_d = c_null_ptr
87 type(c_ptr) :: h_y_idx_d = c_null_ptr
88 type(c_ptr) :: h_z_idx_d = c_null_ptr
89 contains
91 procedure, pass(this) :: init => richardson_init
94 procedure, pass(this) :: partial_init => richardson_partial_init
96 procedure, pass(this) :: finalize => richardson_finalize
98 procedure, pass(this) :: init_from_components => &
101 procedure, pass(this) :: free => richardson_free
103 procedure, pass(this) :: compute => richardson_compute
104 ! Extract fluid properties at the wall (mu and rho)
105 procedure, pass(this) :: extract_properties => richardson_extract_properties
106 end type richardson_t
107
108contains
116 subroutine richardson_init(this, scheme_name, coef, msk, facet, &
117 h_index, json)
118 class(richardson_t), intent(inout) :: this
119 character(len=*), intent(in) :: scheme_name
120 type(coef_t), intent(in) :: coef
121 integer, intent(in) :: msk(:)
122 integer, intent(in) :: facet(:)
123 integer, intent(in) :: h_index
124 type(json_file), intent(inout) :: json
125 real(kind=rp) :: kappa, z0, z0h_in, pr
126 character(len=:), allocatable :: bc_type
127 character(len=:), allocatable :: scalar_name
128 real(kind=rp) :: bc_value
129 real(kind=rp), allocatable :: g_tmp(:)
130 real(kind=rp) :: g(3)
131 logical :: if_time_dependent
132
133 call json_get_or_lookup_or_default(json, "kappa", kappa, 0.4_rp)
134 call json_get_or_lookup_or_default(json, "Pr", pr, 1.0_rp)
135 call json_get_or_lookup(json, "z0", z0)
136 ! If z0h is specified and positive, z0h will be constant and equal to
137 ! what's specified in the case file.
138 ! If z0h is specified and negative, the Zilitinkevich 1995 formulation
139 ! is used, with the specified value acting as -C_Zil.
140 ! If z0h is not specified, assign to have the same value as z0.
141 call json_get_or_lookup_or_default(json, "z0h", z0h_in, z0)
142 call json_get(json, "type_of_temp_bc", bc_type)
143 call json_get(json, "scalar_field", scalar_name)
144 call json_get_or_lookup(json, "bottom_bc_flux_or_temp", bc_value)
145 call json_get_or_default(json, "time_dependent_temp_bc", &
146 if_time_dependent, .false.)
147
148 if (if_time_dependent) then
149 call neko_const_registry%add_real_scalar(bc_value, "bc_value")
150 end if
151
152 call json_get_or_lookup(json, "g", g_tmp)
153 if (size(g_tmp) == 3) then
154 g = g_tmp
155 else
156 call neko_error("Richardson WM: The gravity vector should " &
157 // "have exactly 3 components")
158 end if
159 deallocate(g_tmp)
160
161 call this%init_from_components(scheme_name, scalar_name, coef, &
162 msk, facet, h_index, kappa, g, pr, z0, z0h_in, bc_type, bc_value)
163 deallocate(bc_type)
164 deallocate(scalar_name)
165 end subroutine richardson_init
166
170 subroutine richardson_partial_init(this, coef, json)
171 class(richardson_t), intent(inout) :: this
172 type(coef_t), intent(in) :: coef
173 type(json_file), intent(inout) :: json
174 real(kind=rp), allocatable :: g_tmp(:)
175 character(len=LOG_SIZE) :: log_buf
176 logical :: if_time_dependent
177
178 call this%partial_init_base(coef, json)
179 call json_get_or_lookup_or_default(json, "kappa", this%kappa, 0.4_rp)
180 call json_get_or_lookup_or_default(json, "Pr", this%Pr, 1.0_rp)
181 call json_get_or_lookup(json, "z0", this%z0)
182 call json_get_or_lookup_or_default(json, "z0h", this%z0h_in, this%z0)
183 call json_get(json, "type_of_temp_bc", this%bc_type)
184 call json_get(json, "scalar_field", this%scalar_name)
185 call json_get_or_lookup(json, "bottom_bc_flux_or_temp", this%bc_value)
186 call json_get_or_default(json, "time_dependent_temp_bc", &
187 if_time_dependent, .false.)
188
189 if (if_time_dependent) then
190 call neko_const_registry%add_real_scalar(this%bc_value, "bc_value")
191 end if
192
193
194 call json_get_or_lookup(json, "g", g_tmp)
195 if (size(g_tmp) == 3) then
196 this%g = g_tmp
197 else
198 call neko_error("Richardson WM: The gravity vector should " &
199 // "have exactly 3 components")
200 end if
201 deallocate(g_tmp)
202
203 call neko_log%section('Wall model')
204 write(log_buf, '(A, A)') 'Model : Richardson'
205 call neko_log%message(log_buf)
206 write(log_buf, '(A, A)') 'scalar_name : ', trim(this%scalar_name)
207 call neko_log%message(log_buf)
208 write(log_buf, '(A, A)') 'bc_type : ', trim(this%bc_type)
209 call neko_log%message(log_buf)
210 write(log_buf, '(A, E15.7)') 'bc_value : ', this%bc_value
211 call neko_log%message(log_buf)
212 write(log_buf, '(A, E15.7)') 'kappa : ', this%kappa
213 call neko_log%message(log_buf)
214 write(log_buf, '(A, E15.7)') 'z0 : ', this%z0
215 call neko_log%message(log_buf)
216 write(log_buf, '(A, E15.7)') 'z0h : ', this%z0h_in
217 call neko_log%message(log_buf)
218 write(log_buf, '(A, E15.7)') 'Pr : ', this%Pr
219 call neko_log%message(log_buf)
220 write(log_buf, '(A, 3(E15.7,1X))') 'g : ', this%g
221 call neko_log%message(log_buf)
222 call neko_log%end_section()
223
224 end subroutine richardson_partial_init
225
229 subroutine richardson_finalize(this, msk, facet)
230 class(richardson_t), intent(inout) :: this
231 integer, intent(in) :: msk(:)
232 integer, intent(in) :: facet(:)
233 integer :: i, fid
234
235
236 call this%finalize_base(msk, facet)
237
238 call this%Ri_b%init(this%n_nodes)
239 call this%L_ob%init(this%n_nodes)
240 call this%utau%init(this%n_nodes)
241 call this%magu%init(this%n_nodes)
242 call this%ti%init(this%n_nodes)
243 call this%ts%init(this%n_nodes)
244 call this%q%init(this%n_nodes)
245
246 call this%mu_w%init(this%n_nodes)
247 call this%rho_w%init(this%n_nodes)
248
249 allocate(this%h_x_idx(this%n_nodes))
250 allocate(this%h_y_idx(this%n_nodes))
251 allocate(this%h_z_idx(this%n_nodes))
252 if (neko_bcknd_device .eq. 1) then
253 call device_map(this%h_x_idx, this%h_x_idx_d, this%n_nodes)
254 call device_map(this%h_y_idx, this%h_y_idx_d, this%n_nodes)
255 call device_map(this%h_z_idx, this%h_z_idx_d, this%n_nodes)
256 end if
257
258 do i = 1, this%n_nodes
259 fid = this%facet(i)
260 select case (fid)
261 case (1)
262 this%h_x_idx(i) = this%h_index
263 this%h_y_idx(i) = 0
264 this%h_z_idx(i) = 0
265 case (2)
266 this%h_x_idx(i) = - this%h_index
267 this%h_y_idx(i) = 0
268 this%h_z_idx(i) = 0
269 case (3)
270 this%h_x_idx(i) = 0
271 this%h_y_idx(i) = this%h_index
272 this%h_z_idx(i) = 0
273 case (4)
274 this%h_x_idx(i) = 0
275 this%h_y_idx(i) = - this%h_index
276 this%h_z_idx(i) = 0
277 case (5)
278 this%h_x_idx(i) = 0
279 this%h_y_idx(i) = 0
280 this%h_z_idx(i) = this%h_index
281 case (6)
282 this%h_x_idx(i) = 0
283 this%h_y_idx(i) = 0
284 this%h_z_idx(i) = - this%h_index
285 case default
286 call neko_error("The face index is not correct (richardson.f90)")
287 end select
288 end do
289
290 if (neko_bcknd_device .eq. 1) then
291 call device_memcpy(this%h_x_idx, this%h_x_idx_d, this%n_nodes, &
292 host_to_device, sync = .false.)
293 call device_memcpy(this%h_y_idx, this%h_y_idx_d, this%n_nodes, &
294 host_to_device, sync = .false.)
295 call device_memcpy(this%h_z_idx, this%h_z_idx_d, this%n_nodes, &
296 host_to_device, sync = .true.)
297 end if
298 end subroutine richardson_finalize
299
302 class(richardson_t), intent(inout) :: this
303
304 if (neko_bcknd_device .eq. 1) then
305 call device_masked_gather_copy_0(this%mu_w%x_d, this%mu%x_d, this%msk_d, &
306 this%mu%size(), this%mu_w%size())
307 call device_masked_gather_copy_0(this%rho_w%x_d, this%rho%x_d, this%msk_d, &
308 this%rho%size(), this%rho_w%size())
309 else
310 call masked_gather_copy_0(this%mu_w%x, this%mu%x, this%msk, &
311 this%mu%size(), this%mu_w%size())
312 call masked_gather_copy_0(this%rho_w%x, this%rho%x, this%msk, &
313 this%rho%size(), this%rho_w%size())
314 end if
315 end subroutine richardson_extract_properties
316
330 subroutine richardson_init_from_components(this, scheme_name, &
331 scalar_name, coef, msk, facet, h_index, kappa, g, Pr, z0, &
332 z0h_in, bc_type, bc_value)
333 class(richardson_t), intent(inout) :: this
334 character(len=*), intent(in) :: scheme_name
335 character(len=*), intent(in) :: bc_type
336 character(len=*), intent(in) :: scalar_name
337 type(coef_t), intent(in) :: coef
338 integer, intent(in) :: msk(:)
339 integer, intent(in) :: facet(:)
340 integer, intent(in) :: h_index
341 real(kind=rp), intent(in) :: g(3)
342 real(kind=rp) :: g_mag, g_dot_n, cos_alpha, max_ang
343 integer :: i
344 real(kind=rp), intent(in) :: kappa
345 real(kind=rp), intent(in) :: z0, z0h_in, bc_value, pr
346 character(len=LOG_SIZE) :: log_buf
347
348 call this%free()
349 call this%init_base(scheme_name, coef, msk, facet, h_index)
350
351 this%kappa = kappa
352 this%g = g
353 this%Pr = pr
354 this%z0 = z0
355 this%z0h_in = z0h_in
356 this%bc_type = bc_type
357 this%bc_value = bc_value
358 this%scalar_name = scalar_name
359
360 call this%mu_w%init(this%n_nodes)
361 call this%rho_w%init(this%n_nodes)
362
364 g_mag = sqrt(sum(g**2))
365 if (g_mag < 1.0e-6_rp) then
366 call neko_error("Richardson WM: Gravity magnitude is zero. " &
367 // "Check your input configuration.")
368 end if
369
371 max_ang = 0.0_rp
372 do i = 1, this%n_nodes
373 g_dot_n = abs(g(1)*this%n_x%x(i) + g(2)*this%n_y%x(i) + g(3)*this%n_z%x(i))
374 cos_alpha = g_dot_n / g_mag
375 max_ang = max(max_ang, acos(min(1.0_rp, cos_alpha)))
376 end do
377 max_ang = max_ang * 180.0_rp / (4.0_rp * atan(1.0_rp))
378 if (max_ang > 8.0_rp) then
379 write(log_buf, '(A, F6.2, A)') "Richardson WM: Significant " &
380 // "gravity-normal misalignment (max ", &
381 max_ang, " deg). Stability corrections will use projected gravity."
382 call neko_warning(trim(log_buf))
383 end if
384
386 if (any(this%h%x(1:this%n_nodes) .le. this%z0)) then
387 call neko_error("Richardson WM: Sampling height h must be greater " &
388 // "than roughness z0.")
389 else if ( (this%z0h_in .gt. 0.0_rp) .and. &
390 (any(this%h%x(1:this%n_nodes) .le. this%z0h_in)) ) then
391 call neko_error("Richardson WM: Sampling height h must be greater " &
392 // "than thermal roughness z0h.")
393 else if (this%z0 .eq. 0.0_rp) then
394 call neko_error("Richardson WM: Roughness z0 must be greater than 0.")
395 else if (this%z0h_in .eq. 0.0_rp) then
396 call neko_error("Richardson WM: Thermal roughness z0h must be greater than 0.")
397 end if
398
400
402 subroutine richardson_free(this)
403 class(richardson_t), intent(inout) :: this
404
405 if (allocated(this%bc_type)) then
406 deallocate(this%bc_type)
407 end if
408
409 if (allocated(this%scalar_name)) then
410 deallocate(this%scalar_name)
411 end if
412
413 call this%mu_w%free()
414 call this%rho_w%free()
415 call this%free_base()
416
417 call this%Ri_b%free()
418 call this%L_ob%free()
419 call this%utau%free()
420 call this%magu%free()
421 call this%ti%free()
422 call this%ts%free()
423 call this%q%free()
424
425 if (allocated(this%h_x_idx)) then
426 deallocate(this%h_x_idx)
427 end if
428 if (c_associated(this%h_x_idx_d)) then
429 call device_free(this%h_x_idx_d)
430 end if
431
432 if (allocated(this%h_y_idx)) then
433 deallocate(this%h_y_idx)
434 end if
435 if (c_associated(this%h_y_idx_d)) then
436 call device_free(this%h_y_idx_d)
437 end if
438
439 if (allocated(this%h_z_idx)) then
440 deallocate(this%h_z_idx)
441 end if
442 if (c_associated(this%h_z_idx_d)) then
443 call device_free(this%h_z_idx_d)
444 end if
445
446 end subroutine richardson_free
447
451 subroutine richardson_compute(this, t, tstep)
452 class(richardson_t), intent(inout) :: this
453 real(kind=rp), intent(in) :: t
454 integer, intent(in) :: tstep
455 type(field_t), pointer :: u
456 type(field_t), pointer :: v
457 type(field_t), pointer :: w
458 type(field_t), pointer :: temp
459 real(kind=rp), pointer :: updated_bc_value
460
461 ! Extract boundary values for mu and rho
462 call this%extract_properties()
463
464 u => neko_registry%get_field("u")
465 v => neko_registry%get_field("v")
466 w => neko_registry%get_field("w")
467 temp => neko_registry%get_field(this%scalar_name)
468
469 if (neko_const_registry%real_scalar_exists("bc_value")) then
470 updated_bc_value => neko_const_registry%get_real_scalar("bc_value")
471 this%bc_value = updated_bc_value
472 end if
473
474 if (neko_bcknd_device .eq. 1) then
475 call richardson_compute_device(u%x_d, v%x_d, w%x_d, temp%x_d, &
476 this%ind_r_d, this%ind_s_d, this%ind_t_d, &
477 this%ind_e_d, this%n_x%x_d, this%n_y%x_d, this%n_z%x_d, &
478 this%h%x_d, this%tau_x%x_d, this%tau_y%x_d, &
479 this%tau_z%x_d, this%n_nodes, u%Xh%lx, this%kappa, &
480 this%mu_w%x_d, this%rho_w%x_d, this%g, this%Pr, this%z0, this%z0h_in, &
481 this%bc_type, this%bc_value, tstep, this%Ri_b%x_d, &
482 this%L_ob%x_d, this%utau%x_d, this%magu%x_d, this%ti%x_d, this%ts%x_d,&
483 this%q%x_d, this%h_x_idx_d, this%h_y_idx_d, this%h_z_idx_d)
484 else
485 call richardson_compute_cpu(u%x, v%x, w%x, temp%x, this%ind_r, &
486 this%ind_s, this%ind_t, this%ind_e, this%n_x%x, &
487 this%n_y%x, this%n_z%x, this%h%x, this%tau_x%x, &
488 this%tau_y%x, this%tau_z%x, this%n_nodes, u%Xh%lx, &
489 u%msh%nelv, this%kappa, this%mu_w%x, this%rho_w%x, &
490 this%g, this%Pr, this%z0, this%z0h_in, this%bc_type, &
491 this%bc_value, tstep, this%Ri_b%x, this%L_ob%x, &
492 this%utau%x, this%magu%x, this%ti%x, this%ts%x, &
493 this%q%x, this%h_x_idx, this%h_y_idx, this%h_z_idx)
494 end if
495
496 call richardson_log_diagnostics(this%Ri_b, this%L_ob, &
497 this%utau, this%magu, this%ti, this%ts, this%q, &
498 this%n_nodes, this%bc_value)
499 end subroutine richardson_compute
500
501 subroutine richardson_log_diagnostics(Ri_b, L_ob, utau, magu, &
502 ti, ts, q, n_nodes, bc_value)
503 character(len=LOG_SIZE) :: log_buf
504 integer, intent(in) :: n_nodes
505 real(kind=rp), intent(in) :: bc_value
506 type(vector_t), intent(in) :: ri_b, l_ob, utau
507 type(vector_t), intent(in) :: magu, ti, ts, q
508
509 call neko_log%section("Wall model diagnostics")
510 write(log_buf, '(A)') '--- sum --- ' ! min max ---'
511 call neko_log%message(trim(log_buf))
512 write(log_buf, '(A, E15.7)') "Ri_b: ", &
513 vector_glsum(ri_b, n_nodes) !, &
514 ! vector_glmin(Ri_b, n_nodes), vector_glmax(Ri_b, n_nodes)
515 call neko_log%message(trim(log_buf))
516
517 write(log_buf, '(A, E15.7)') "L_ob: ", &
518 vector_glsum(l_ob, n_nodes) !, &
519 ! vector_glmin(L_ob, n_nodes), vector_glmax(L_ob, n_nodes)
520 call neko_log%message(trim(log_buf))
521
522 write(log_buf, '(A, E15.7)') "utau: ", &
523 vector_glsum(utau, n_nodes) !, &
524 ! vector_glmin(utau, n_nodes), vector_glmax(utau, n_nodes)
525 call neko_log%message(trim(log_buf))
526
527 write(log_buf, '(A, E15.7)') "magu: ", &
528 vector_glsum(magu, n_nodes) !, &
529 ! vector_glmin(magu, n_nodes), vector_glmax(magu, n_nodes)
530 call neko_log%message(trim(log_buf))
531
532 write(log_buf, '(A, E15.7)') "ti: ", &
533 vector_glsum(ti, n_nodes) !, &
534 ! vector_glmin(ti, n_nodes), vector_glmax(ti, n_nodes)
535 call neko_log%message(trim(log_buf))
536
537 write(log_buf, '(A, E15.7)') "ts: ", &
538 vector_glsum(ts, n_nodes) !, &
539 ! vector_glmin(ts, n_nodes), vector_glmax(ts, n_nodes)
540 call neko_log%message(trim(log_buf))
541
542 write(log_buf, '(A, E15.7)') "q: ", &
543 vector_glsum(q, n_nodes) !, &
544 ! vector_glmin(q, n_nodes), vector_glmax(q, n_nodes)
545 call neko_log%message(trim(log_buf))
546
547 write(log_buf, '(A, E15.7)') "bc_value: ", bc_value
548 call neko_log%message(trim(log_buf))
549
550 call neko_log%end_section()
551
552 end subroutine richardson_log_diagnostics
553
554
555end module richardson
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
__global__ void richardson_compute(const T *__restrict__ u_d, const T *__restrict__ v_d, const T *__restrict__ w_d, const T *__restrict__ temp_d, const T *__restrict__ h_d, const T *__restrict__ n_x_d, const T *__restrict__ n_y_d, const T *__restrict__ n_z_d, const int *__restrict__ ind_r_d, const int *__restrict__ ind_s_d, const int *__restrict__ ind_t_d, const int *__restrict__ ind_e_d, T *__restrict__ tau_x_d, T *__restrict__ tau_y_d, T *__restrict__ tau_z_d, int n_nodes, int lx, T kappa, const T *__restrict__ mu_w_d, const T *__restrict__ rho_w_d, T g1, T g2, T g3, T Pr, T z0, T z0h_in, T bc_value, T *__restrict__ Ri_b_diagn, T *__restrict__ L_ob_diagn, T *__restrict__ utau_diagn, T *__restrict__ magu_diagn, T *__restrict__ ti_diagn, T *__restrict__ ts_diagn, T *__restrict__ q_diagn, const int *__restrict__ h_x_idx, const int *__restrict__ h_y_idx, const int *__restrict__ h_z_idx)
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Coefficients.
Definition coef.f90:34
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:225
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
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
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:315
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
type(registry_t), target, public neko_const_registry
This registry is used to store user-defined scalars and vectors, provided under the constants section...
Definition registry.f90:150
Implements the CPU kernel for the richardson_t type.
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.
Implements the device kernel for the richardson_t type.
subroutine, public richardson_compute_device(u_d, v_d, w_d, temp_d, ind_r_d, ind_s_d, ind_t_d, ind_e_d, n_x_d, n_y_d, n_z_d, h_d, tau_x_d, tau_y_d, tau_z_d, n_nodes, lx, kappa, mu_w_d, rho_w_d, g, 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)
Compute the wall shear stress on device using the rough log-law model.
Implements richardson_t.
subroutine richardson_init_from_components(this, scheme_name, scalar_name, coef, msk, facet, h_index, kappa, g, pr, z0, z0h_in, bc_type, bc_value)
Constructor from components.
subroutine richardson_extract_properties(this)
Extract the values of rho and mu at the boundary.
subroutine richardson_log_diagnostics(ri_b, l_ob, utau, magu, ti, ts, q, n_nodes, bc_value)
subroutine richardson_partial_init(this, coef, json)
Constructor from JSON.
subroutine richardson_finalize(this, msk, facet)
Finalize the construction using the mask and facet arrays of the bc.
subroutine richardson_init(this, scheme_name, coef, msk, facet, h_index, json)
Constructor from JSON.
subroutine richardson_free(this)
Destructor for the richardson_t (base) class.
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:392
real(kind=rp) function, public vector_glmin(a, n)
Global minimum of all elements in a vector .
real(kind=rp) function, public vector_glsum(a, n)
real(kind=rp) function, public vector_glmax(a, n)
Global maximum of all elements in a vector .
Defines a vector.
Definition vector.f90:34
Implements wall_model_t.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Wall model similar to the Monin-Obukhov Similarity Theory for atmospheric boundary layer flows,...
Base abstract type for wall-stress models for wall-modelled LES.
#define max(a, b)
Definition tensor.cu:40