38 use json_module,
only : json_file
54 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr, c_associated
64 real(kind=
rp) :: kappa
70 real(kind=
rp) :: z0h_in
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
98 procedure, pass(this) :: init_from_components => &
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
142 call json_get(json,
"type_of_temp_bc", bc_type)
143 call json_get(json,
"scalar_field", scalar_name)
146 if_time_dependent, .false.)
148 if (if_time_dependent)
then
153 if (
size(g_tmp) == 3)
then
156 call neko_error(
"Richardson WM: The gravity vector should " &
157 //
"have exactly 3 components")
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)
164 deallocate(scalar_name)
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
178 call this%partial_init_base(coef, json)
183 call json_get(json,
"type_of_temp_bc", this%bc_type)
184 call json_get(json,
"scalar_field", this%scalar_name)
187 if_time_dependent, .false.)
189 if (if_time_dependent)
then
195 if (
size(g_tmp) == 3)
then
198 call neko_error(
"Richardson WM: The gravity vector should " &
199 //
"have exactly 3 components")
204 write(log_buf,
'(A, A)')
'Model : Richardson'
206 write(log_buf,
'(A, A)')
'scalar_name : ', trim(this%scalar_name)
208 write(log_buf,
'(A, A)')
'bc_type : ', trim(this%bc_type)
210 write(log_buf,
'(A, E15.7)')
'bc_value : ', this%bc_value
212 write(log_buf,
'(A, E15.7)')
'kappa : ', this%kappa
214 write(log_buf,
'(A, E15.7)')
'z0 : ', this%z0
216 write(log_buf,
'(A, E15.7)')
'z0h : ', this%z0h_in
218 write(log_buf,
'(A, E15.7)')
'Pr : ', this%Pr
220 write(log_buf,
'(A, 3(E15.7,1X))')
'g : ', this%g
231 integer,
intent(in) :: msk(:)
232 integer,
intent(in) :: facet(:)
236 call this%finalize_base(msk, facet)
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)
246 call this%mu_w%init(this%n_nodes)
247 call this%rho_w%init(this%n_nodes)
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))
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)
258 do i = 1, this%n_nodes
262 this%h_x_idx(i) = this%h_index
266 this%h_x_idx(i) = - this%h_index
271 this%h_y_idx(i) = this%h_index
275 this%h_y_idx(i) = - this%h_index
280 this%h_z_idx(i) = this%h_index
284 this%h_z_idx(i) = - this%h_index
286 call neko_error(
"The face index is not correct (richardson.f90)")
291 call device_memcpy(this%h_x_idx, this%h_x_idx_d, this%n_nodes, &
293 call device_memcpy(this%h_y_idx, this%h_y_idx_d, this%n_nodes, &
295 call device_memcpy(this%h_z_idx, this%h_z_idx_d, this%n_nodes, &
306 this%mu%size(), this%mu_w%size())
308 this%rho%size(), this%rho_w%size())
311 this%mu%size(), this%mu_w%size())
313 this%rho%size(), this%rho_w%size())
331 scalar_name, coef, msk, facet, h_index, kappa, g, Pr, z0, &
332 z0h_in, bc_type, bc_value)
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
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
349 call this%init_base(scheme_name, coef, msk, facet, h_index)
356 this%bc_type = bc_type
357 this%bc_value = bc_value
358 this%scalar_name = scalar_name
360 call this%mu_w%init(this%n_nodes)
361 call this%rho_w%init(this%n_nodes)
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.")
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)))
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."
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.")
405 if (
allocated(this%bc_type))
then
406 deallocate(this%bc_type)
409 if (
allocated(this%scalar_name))
then
410 deallocate(this%scalar_name)
413 call this%mu_w%free()
414 call this%rho_w%free()
415 call this%free_base()
417 call this%Ri_b%free()
418 call this%L_ob%free()
419 call this%utau%free()
420 call this%magu%free()
425 if (
allocated(this%h_x_idx))
then
426 deallocate(this%h_x_idx)
428 if (c_associated(this%h_x_idx_d))
then
432 if (
allocated(this%h_y_idx))
then
433 deallocate(this%h_y_idx)
435 if (c_associated(this%h_y_idx_d))
then
439 if (
allocated(this%h_z_idx))
then
440 deallocate(this%h_z_idx)
442 if (c_associated(this%h_z_idx_d))
then
453 real(kind=
rp),
intent(in) :: t
454 integer,
intent(in) :: tstep
459 real(kind=
rp),
pointer :: updated_bc_value
462 call this%extract_properties()
471 this%bc_value = updated_bc_value
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)
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)
497 this%utau, this%magu, this%ti, this%ts, this%q, &
498 this%n_nodes, this%bc_value)
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
509 call neko_log%section(
"Wall model diagnostics")
510 write(log_buf,
'(A)')
'--- sum --- '
511 call neko_log%message(trim(log_buf))
512 write(log_buf,
'(A, E15.7)')
"Ri_b: ", &
515 call neko_log%message(trim(log_buf))
517 write(log_buf,
'(A, E15.7)')
"L_ob: ", &
520 call neko_log%message(trim(log_buf))
522 write(log_buf,
'(A, E15.7)')
"utau: ", &
525 call neko_log%message(trim(log_buf))
527 write(log_buf,
'(A, E15.7)')
"magu: ", &
530 call neko_log%message(trim(log_buf))
532 write(log_buf,
'(A, E15.7)')
"ti: ", &
535 call neko_log%message(trim(log_buf))
537 write(log_buf,
'(A, E15.7)')
"ts: ", &
540 call neko_log%message(trim(log_buf))
542 write(log_buf,
'(A, E15.7)')
"q: ", &
545 call neko_log%message(trim(log_buf))
547 write(log_buf,
'(A, E15.7)')
"bc_value: ", bc_value
548 call neko_log%message(trim(log_buf))
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
__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)
Copy data between host and device (or device and device)
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.
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.
integer, parameter, public host_to_device
subroutine, public device_free(x_d)
Deallocate memory on the device.
Utilities for retrieving parameters from the case files.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
type(registry_t), target, public neko_const_registry
This registry is used to store user-defined scalars and vectors, provided under the constants section...
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.
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.
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
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 .
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
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.