Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
most.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!
35module most
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
46 use most_cpu, only : most_compute_cpu
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
63 type, public, extends(wall_model_t) :: most_t
65 real(kind=rp) :: kappa
67 real(kind=rp) :: pr
69 real(kind=rp) :: z0
71 real(kind=rp) :: z0h_in
72 ! The dynamic viscosity at the boundary
73 type(vector_t) :: mu_w
74 ! The fluid density at the boundary
75 type(vector_t) :: rho_w
77 real(kind=rp) :: g(3)
79 character(len=:), allocatable :: bc_type
81 real(kind=rp) :: bc_value
83 character(len=:), allocatable :: scalar_name
85 type(vector_t) :: ri_b, l_ob, utau, magu, ti, ts, q
86 integer, allocatable :: h_x_idx(:), h_y_idx(:), h_z_idx(:)
87 type(c_ptr) :: h_x_idx_d = c_null_ptr
88 type(c_ptr) :: h_y_idx_d = c_null_ptr
89 type(c_ptr) :: h_z_idx_d = c_null_ptr
90 contains
92 procedure, pass(this) :: init => most_init
95 procedure, pass(this) :: partial_init => most_partial_init
97 procedure, pass(this) :: finalize => most_finalize
99 procedure, pass(this) :: init_from_components => &
102 procedure, pass(this) :: free => most_free
104 procedure, pass(this) :: compute => most_compute
105 ! Extract fluid properties at the wall (mu and rho)
106 procedure, pass(this) :: extract_properties => most_extract_properties
107 end type most_t
108
109contains
117 subroutine most_init(this, scheme_name, coef, msk, facet, &
118 h_index, json)
119 class(most_t), intent(inout) :: this
120 character(len=*), intent(in) :: scheme_name
121 type(coef_t), intent(in) :: coef
122 integer, intent(in) :: msk(:)
123 integer, intent(in) :: facet(:)
124 integer, intent(in) :: h_index
125 type(json_file), intent(inout) :: json
126 real(kind=rp) :: kappa, z0, z0h_in, pr
127 character(len=:), allocatable :: bc_type
128 character(len=:), allocatable :: scalar_name
129 real(kind=rp) :: bc_value
130 real(kind=rp), allocatable :: g_tmp(:)
131 real(kind=rp) :: g(3)
132 logical :: if_time_dependent
133
134 call json_get_or_lookup_or_default(json, "kappa", kappa, 0.4_rp)
135 call json_get_or_lookup_or_default(json, "Pr", pr, 1.0_rp)
136 call json_get_or_lookup(json, "z0", z0)
137 ! If z0h is specified and positive, z0h will be constant and equal to
138 ! what's specified in the case file.
139 ! If z0h is specified and negative, the Zilitinkevich 1995 formulation
140 ! is used, with the specified value acting as -C_Zil.
141 ! If z0h is not specified, assign to have the same value as z0.
142 call json_get_or_lookup_or_default(json, "z0h", z0h_in, z0)
143 call json_get(json, "type_of_temp_bc", bc_type)
144 call json_get(json, "scalar_field", scalar_name)
145 call json_get_or_lookup(json, "bottom_bc_flux_or_temp", bc_value)
146 call json_get_or_default(json, "time_dependent_temp_bc", &
147 if_time_dependent, .false.)
148
149 if (if_time_dependent) then
150 call neko_const_registry%add_real_scalar(bc_value, "bc_value")
151 end if
152
153 call json_get_or_lookup(json, "g", g_tmp)
154 if (size(g_tmp) == 3) then
155 g = g_tmp
156 else
157 call neko_error("MOST WM: The gravity vector should have exactly 3 components")
158 end if
159 deallocate(g_tmp)
160
161 call this%init_from_components(scheme_name, scalar_name, coef, msk, facet, h_index, &
162 kappa, g, pr, z0, z0h_in, bc_type, bc_value)
163 deallocate(bc_type)
164 deallocate(scalar_name)
165 end subroutine most_init
166
170 subroutine most_partial_init(this, coef, json)
171 class(most_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("MOST WM: The gravity vector should have exactly 3 components")
199 end if
200 deallocate(g_tmp)
201
202 call neko_log%section('Wall model')
203 write(log_buf, '(A, A)') 'Model : MOST'
204 call neko_log%message(log_buf)
205 write(log_buf, '(A, A)') 'scalar_name : ', trim(this%scalar_name)
206 call neko_log%message(log_buf)
207 write(log_buf, '(A, A)') 'bc_type : ', trim(this%bc_type)
208 call neko_log%message(log_buf)
209 write(log_buf, '(A, E15.7)') 'bc_value : ', this%bc_value
210 call neko_log%message(log_buf)
211 write(log_buf, '(A, E15.7)') 'kappa : ', this%kappa
212 call neko_log%message(log_buf)
213 write(log_buf, '(A, E15.7)') 'z0 : ', this%z0
214 call neko_log%message(log_buf)
215 write(log_buf, '(A, E15.7)') 'z0h : ', this%z0h_in
216 call neko_log%message(log_buf)
217 write(log_buf, '(A, E15.7)') 'Pr : ', this%Pr
218 call neko_log%message(log_buf)
219 write(log_buf, '(A, 3(E15.7,1X))') 'g : ', this%g
220 call neko_log%message(log_buf)
221 call neko_log%end_section()
222
223 end subroutine most_partial_init
224
228 subroutine most_finalize(this, msk, facet)
229 class(most_t), intent(inout) :: this
230 integer, intent(in) :: msk(:)
231 integer, intent(in) :: facet(:)
232 integer :: i, fid
233
234
235 call this%finalize_base(msk, facet)
236
237 call this%Ri_b%init(this%n_nodes)
238 call this%L_ob%init(this%n_nodes)
239 call this%utau%init(this%n_nodes)
240 call this%magu%init(this%n_nodes)
241 call this%ti%init(this%n_nodes)
242 call this%ts%init(this%n_nodes)
243 call this%q%init(this%n_nodes)
244
245 call this%mu_w%init(this%n_nodes)
246 call this%rho_w%init(this%n_nodes)
247
248 allocate(this%h_x_idx(this%n_nodes))
249 allocate(this%h_y_idx(this%n_nodes))
250 allocate(this%h_z_idx(this%n_nodes))
251 if (neko_bcknd_device .eq. 1) then
252 call device_map(this%h_x_idx, this%h_x_idx_d, this%n_nodes)
253 call device_map(this%h_y_idx, this%h_y_idx_d, this%n_nodes)
254 call device_map(this%h_z_idx, this%h_z_idx_d, this%n_nodes)
255 end if
256
257 do i = 1, this%n_nodes
258 fid = this%facet(i)
259 select case (fid)
260 case (1)
261 this%h_x_idx(i) = this%h_index
262 this%h_y_idx(i) = 0
263 this%h_z_idx(i) = 0
264 case (2)
265 this%h_x_idx(i) = - this%h_index
266 this%h_y_idx(i) = 0
267 this%h_z_idx(i) = 0
268 case (3)
269 this%h_x_idx(i) = 0
270 this%h_y_idx(i) = this%h_index
271 this%h_z_idx(i) = 0
272 case (4)
273 this%h_x_idx(i) = 0
274 this%h_y_idx(i) = - this%h_index
275 this%h_z_idx(i) = 0
276 case (5)
277 this%h_x_idx(i) = 0
278 this%h_y_idx(i) = 0
279 this%h_z_idx(i) = this%h_index
280 case (6)
281 this%h_x_idx(i) = 0
282 this%h_y_idx(i) = 0
283 this%h_z_idx(i) = - this%h_index
284 case default
285 call neko_error("The face index is not correct (most.f90)")
286 end select
287
288 end do
289 if (neko_bcknd_device .eq. 1) then
290 call device_memcpy(this%h_x_idx, this%h_x_idx_d, this%n_nodes, &
291 host_to_device, sync = .false.)
292 call device_memcpy(this%h_y_idx, this%h_y_idx_d, this%n_nodes, &
293 host_to_device, sync = .false.)
294 call device_memcpy(this%h_z_idx, this%h_z_idx_d, this%n_nodes, &
295 host_to_device, sync = .false.)
296 end if
297 end subroutine most_finalize
298
300 subroutine most_extract_properties(this)
301 class(most_t), intent(inout) :: this
302
303 if (neko_bcknd_device .eq. 1) then
304 call device_masked_gather_copy_0(this%mu_w%x_d, this%mu%x_d, this%msk_d, &
305 this%mu%size(), this%mu_w%size())
306 call device_masked_gather_copy_0(this%rho_w%x_d, this%rho%x_d, this%msk_d, &
307 this%rho%size(), this%rho_w%size())
308 else
309 call masked_gather_copy_0(this%mu_w%x, this%mu%x, this%msk, &
310 this%mu%size(), this%mu_w%size())
311 call masked_gather_copy_0(this%rho_w%x, this%rho%x, this%msk, &
312 this%rho%size(), this%rho_w%size())
313 end if
314 end subroutine most_extract_properties
315
329 subroutine most_init_from_components(this, scheme_name, scalar_name, &
330 coef, msk, facet, h_index, kappa, g, Pr, z0, &
331 z0h_in, bc_type, bc_value)
332 class(most_t), intent(inout) :: this
333 character(len=*), intent(in) :: scheme_name
334 character(len=*), intent(in) :: bc_type
335 character(len=*), intent(in) :: scalar_name
336 type(coef_t), intent(in) :: coef
337 integer, intent(in) :: msk(:)
338 integer, intent(in) :: facet(:)
339 integer, intent(in) :: h_index
340 real(kind=rp), intent(in) :: g(3)
341 real(kind=rp) :: g_mag, g_dot_n, cos_alpha, max_ang
342 integer :: i
343 real(kind=rp), intent(in) :: kappa
344 real(kind=rp), intent(in) :: z0, z0h_in, bc_value, pr
345 character(len=LOG_SIZE) :: log_buf
346
347 call this%free()
348 call this%init_base(scheme_name, coef, msk, facet, h_index)
349
350 this%kappa = kappa
351 this%g = g
352 this%Pr = pr
353 this%z0 = z0
354 this%z0h_in = z0h_in
355 this%bc_type = bc_type
356 this%bc_value = bc_value
357 this%scalar_name = scalar_name
358
359 call this%mu_w%init(this%n_nodes)
360 call this%rho_w%init(this%n_nodes)
361
363 g_mag = sqrt(sum(g**2))
364 if (g_mag < 1.0e-6_rp) then
365 call neko_error("MOST WM: Gravity magnitude is zero. Check your input configuration.")
366 end if
367
369 max_ang = 0.0_rp
370 do i = 1, this%n_nodes
371 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))
372 cos_alpha = g_dot_n / g_mag
373 max_ang = max(max_ang, acos(min(1.0_rp, cos_alpha)))
374 end do
375 max_ang = max_ang * 180.0_rp / (4.0_rp * atan(1.0_rp))
376 if (max_ang > 8.0_rp) then
377 write(log_buf, '(A, F6.2, A)') "MOST WM: Significant gravity-normal misalignment (max ", &
378 max_ang, " deg). Stability corrections will use projected gravity."
379 call neko_warning(trim(log_buf))
380 end if
381
383 if (any(this%h%x(1:this%n_nodes) .le. this%z0)) then
384 call neko_error("MOST WM: Sampling height h must be greater than roughness z0.")
385 else if ( (this%z0h_in .gt. 0.0_rp) .and. (any(this%h%x(1:this%n_nodes) .le. this%z0h_in)) ) then
386 call neko_error("MOST WM: Sampling height h must be greater than thermal roughness z0h.")
387 else if (this%z0 .eq. 0.0_rp) then
388 call neko_error("MOST WM: Roughness z0 must be greater than 0.")
389 else if (this%z0h_in .eq. 0.0_rp) then
390 call neko_error("MOST WM: Thermal roughness z0h must be greater than 0.")
391 end if
392
393 end subroutine most_init_from_components
394
396 subroutine most_free(this)
397 class(most_t), intent(inout) :: this
398
399 if (allocated(this%bc_type)) then
400 deallocate(this%bc_type)
401 end if
402
403 if (allocated(this%scalar_name)) then
404 deallocate(this%scalar_name)
405 end if
406
407 call this%mu_w%free()
408 call this%rho_w%free()
409 call this%free_base()
410
411 call this%Ri_b%free()
412 call this%L_ob%free()
413 call this%utau%free()
414 call this%magu%free()
415 call this%ti%free()
416 call this%ts%free()
417 call this%q%free()
418
419 if (allocated(this%h_x_idx)) then
420 deallocate(this%h_x_idx)
421 end if
422 if (c_associated(this%h_x_idx_d)) then
423 call device_free(this%h_x_idx_d)
424 end if
425
426 if (allocated(this%h_y_idx)) then
427 deallocate(this%h_y_idx)
428 end if
429 if (c_associated(this%h_y_idx_d)) then
430 call device_free(this%h_y_idx_d)
431 end if
432
433 if (allocated(this%h_z_idx)) then
434 deallocate(this%h_z_idx)
435 end if
436 if (c_associated(this%h_z_idx_d)) then
437 call device_free(this%h_z_idx_d)
438 end if
439
440 end subroutine most_free
441
445 subroutine most_compute(this, t, tstep)
446 class(most_t), intent(inout) :: this
447 real(kind=rp), intent(in) :: t
448 integer, intent(in) :: tstep
449 type(field_t), pointer :: u
450 type(field_t), pointer :: v
451 type(field_t), pointer :: w
452 type(field_t), pointer :: temp
453 real(kind=rp), pointer :: updated_bc_value
454
455 ! Extract boundary values for mu and rho
456 call this%extract_properties()
457
458 u => neko_registry%get_field("u")
459 v => neko_registry%get_field("v")
460 w => neko_registry%get_field("w")
461 temp => neko_registry%get_field(this%scalar_name)
462
463 if (neko_const_registry%real_scalar_exists("bc_value")) then
464 updated_bc_value => neko_const_registry%get_real_scalar("bc_value")
465 this%bc_value = updated_bc_value
466 end if
467
468 if (neko_bcknd_device .eq. 1) then
469 call most_compute_device(u%x_d, v%x_d, w%x_d, temp%x_d, &
470 this%ind_r_d, this%ind_s_d, this%ind_t_d, &
471 this%ind_e_d, this%n_x%x_d, this%n_y%x_d, this%n_z%x_d, &
472 this%h%x_d, this%tau_x%x_d, this%tau_y%x_d, &
473 this%tau_z%x_d, this%n_nodes, u%Xh%lx, this%kappa, &
474 this%mu_w%x_d, this%rho_w%x_d, this%g, this%Pr, this%z0, this%z0h_in, &
475 this%bc_type, this%bc_value, tstep, this%Ri_b%x_d, &
476 this%L_ob%x_d, this%utau%x_d, this%magu%x_d, this%ti%x_d, this%ts%x_d,&
477 this%q%x_d, this%h_x_idx_d, this%h_y_idx_d, this%h_z_idx_d)
478 else
479 call most_compute_cpu(u%x, v%x, w%x, temp%x, this%ind_r, &
480 this%ind_s, this%ind_t, this%ind_e, this%n_x%x, &
481 this%n_y%x, this%n_z%x, this%h%x, this%tau_x%x, &
482 this%tau_y%x, this%tau_z%x, this%n_nodes, u%Xh%lx, &
483 u%msh%nelv, this%kappa, this%mu_w%x, this%rho_w%x, &
484 this%g, this%Pr, this%z0, this%z0h_in, this%bc_type, &
485 this%bc_value, tstep, this%Ri_b%x, this%L_ob%x, &
486 this%utau%x, this%magu%x, this%ti%x, this%ts%x, &
487 this%q%x, this%h_x_idx, this%h_y_idx, this%h_z_idx)
488 end if
489
490 call most_log_diagnostics(this%Ri_b, this%L_ob, &
491 this%utau, this%magu, this%ti, this%ts, this%q, &
492 this%n_nodes, this%bc_value)
493 end subroutine most_compute
494
495 subroutine most_log_diagnostics(Ri_b, L_ob, utau, magu, ti, ts, q, &
496 n_nodes, bc_value)
497 character(len=LOG_SIZE) :: log_buf
498 integer, intent(in) :: n_nodes
499 real(kind=rp), intent(in) :: bc_value
500 type(vector_t), intent(in) :: ri_b, l_ob, utau
501 type(vector_t), intent(in) :: magu, ti, ts, q
502
503 call neko_log%section("Wall model diagnostics")
504 write(log_buf, '(A)') '--- sum --- ' !min max ---'
505 call neko_log%message(trim(log_buf))
506 write(log_buf,'(A,3E15.7)') "Ri_b: ",&
507 vector_glsum(ri_b, n_nodes) !, &
508 ! vector_glmin(Ri_b, n_nodes), vector_glmax(Ri_b, n_nodes)
509 call neko_log%message(trim(log_buf))
510
511 write(log_buf,'(A,3E15.7)') "L_ob: ", &
512 vector_glsum(l_ob, n_nodes)!, &
513 ! vector_glmin(L_ob, n_nodes), vector_glmax(L_ob, n_nodes)
514 call neko_log%message(trim(log_buf))
515
516 write(log_buf,'(A,3E15.7)') "utau: ", &
517 vector_glsum(utau, n_nodes)!, &
518 ! vector_glmin(utau, n_nodes), vector_glmax(utau, n_nodes)
519 call neko_log%message(trim(log_buf))
520
521 write(log_buf,'(A,3E15.7)') "magu: ", &
522 vector_glsum(magu, n_nodes)!, &
523 ! vector_glmin(magu, n_nodes), vector_glmax(magu, n_nodes)
524 call neko_log%message(trim(log_buf))
525
526 write(log_buf,'(A,3E15.7)') "ti: ", &
527 vector_glsum(ti, n_nodes)!, &
528 ! vector_glmin(ti, n_nodes), vector_glmax(ti, n_nodes)
529 call neko_log%message(trim(log_buf))
530
531 write(log_buf,'(A,3E15.7)') "ts: ", &
532 vector_glsum(ts, n_nodes)!, &
533 ! vector_glmin(ts, n_nodes), vector_glmax(ts, n_nodes)
534 call neko_log%message(trim(log_buf))
535
536 write(log_buf,'(A,3E15.7)') "q: ", &
537 vector_glsum(q, n_nodes)!, &
538 ! vector_glmin(q, n_nodes), vector_glmax(q, n_nodes)
539 call neko_log%message(trim(log_buf))
540
541 write(log_buf,'(A,E15.7)') "bc_value: ", bc_value
542 call neko_log%message(trim(log_buf))
543
544 call neko_log%end_section()
545
546 end subroutine most_log_diagnostics
547
548
549end module most
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
__global__ void most_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
Implements the CPU kernel for the most_t type.
Definition most_cpu.f90:34
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
Implements the device kernel for the most_t type.
subroutine, public most_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 most_t.
Definition most.f90:35
subroutine most_init(this, scheme_name, coef, msk, facet, h_index, json)
Constructor from JSON.
Definition most.f90:119
subroutine most_finalize(this, msk, facet)
Finalize the construction using the mask and facet arrays of the bc.
Definition most.f90:229
subroutine most_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.
Definition most.f90:332
subroutine most_free(this)
Destructor for the most_t (base) class.
Definition most.f90:397
subroutine most_log_diagnostics(ri_b, l_ob, utau, magu, ti, ts, q, n_nodes, bc_value)
Definition most.f90:497
subroutine most_extract_properties(this)
Extract the values of rho and mu at the boundary.
Definition most.f90:301
subroutine most_partial_init(this, coef, json)
Constructor from JSON.
Definition most.f90:171
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
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 based on the Monin-Obukhov Similarity Theory for atmospheric boundary layer flows....
Definition most.f90:63
Base abstract type for wall-stress models for wall-modelled LES.
#define max(a, b)
Definition tensor.cu:40