Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
wall_model.f90
Go to the documentation of this file.
1! Copyright (c) 2024, 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 num_types, only : rp
37 use field, only : field_t
38 use json_module, only : json_file
39 use registry, only : neko_registry
40 use dofmap, only : dofmap_t
41 use coefs, only : coef_t
44 use vector, only : vector_t
46 use math, only : glmin, glmax
47 use comm, only : pe_rank
49 use file, only : file_t
50 use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
54 implicit none
55 private
56
58 type, abstract, public :: wall_model_t
60 type(coef_t), pointer :: coef => null()
62 type(dofmap_t), pointer :: dof => null()
64 type(field_t), pointer :: mu => null()
66 type(field_t), pointer :: rho => null()
70 character(len=:), allocatable :: scheme_name
72 integer, pointer :: msk(:) => null()
73 type(c_ptr) :: msk_d = c_null_ptr
75 integer, pointer :: facet(:) => null()
77 type(vector_t) :: tau_x
79 type(vector_t) :: tau_y
81 type(vector_t) :: tau_z
83 type(vector_t) :: n_x
85 type(vector_t) :: n_y
87 type(vector_t) :: n_z
89 integer, allocatable :: ind_r(:)
90 type(c_ptr) :: ind_r_d = c_null_ptr
92 integer, allocatable :: ind_s(:)
93 type(c_ptr) :: ind_s_d = c_null_ptr
95 integer, allocatable :: ind_t(:)
96 type(c_ptr) :: ind_t_d = c_null_ptr
98 integer, allocatable :: ind_e(:)
99 type(c_ptr) :: ind_e_d = c_null_ptr
101 type(vector_t) :: h
103 integer :: h_index = 0
105 integer :: n_nodes = 0
107 type(field_t), pointer :: tau_field => null()
108 contains
110 procedure, pass(this) :: init_base => wall_model_init_base
112 procedure, pass(this) :: partial_init_base => wall_model_partial_init_base
114 procedure, pass(this) :: finalize_base => wall_model_finalize_base
116 procedure, pass(this) :: free_base => wall_model_free_base
118 procedure, pass(this) :: compute_mag_field => wall_model_compute_mag_field
120 procedure(wall_model_init), pass(this), deferred :: init
128 procedure(wall_model_partial_init), pass(this), deferred :: partial_init
130 procedure(wall_model_finalize), pass(this), deferred :: finalize
132 procedure(wall_model_free), pass(this), deferred :: free
134 procedure(wall_model_compute), pass(this), deferred :: compute
136 procedure, pass(this) :: find_points => wall_model_find_points
137 end type wall_model_t
138
139 abstract interface
140
143 subroutine wall_model_compute(this, t, tstep)
144 import wall_model_t, rp
145 class(wall_model_t), intent(inout) :: this
146 real(kind=rp), intent(in) :: t
147 integer, intent(in) :: tstep
148 end subroutine wall_model_compute
149 end interface
150
151 abstract interface
152
160 subroutine wall_model_init(this, scheme_name, coef, msk, facet, &
161 h_index, json)
162 import wall_model_t, json_file, dofmap_t, coef_t, rp
163 class(wall_model_t), intent(inout) :: this
164 character(len=*), intent(in) :: scheme_name
165 type(coef_t), intent(in) :: coef
166 integer, intent(in) :: msk(:)
167 integer, intent(in) :: facet(:)
168 integer, intent(in) :: h_index
169 type(json_file), intent(inout) :: json
170 end subroutine wall_model_init
171 end interface
172
173 abstract interface
174
178 subroutine wall_model_partial_init(this, coef, json)
179 import wall_model_t, json_file, dofmap_t, coef_t, rp
180 class(wall_model_t), intent(inout) :: this
181 type(coef_t), intent(in) :: coef
182 type(json_file), intent(inout) :: json
183 end subroutine wall_model_partial_init
184 end interface
185
186 abstract interface
187
190 subroutine wall_model_finalize(this, msk, facet)
191 import wall_model_t
192 class(wall_model_t), intent(inout) :: this
193 integer, intent(in) :: msk(:)
194 integer, intent(in) :: facet(:)
195 end subroutine wall_model_finalize
196 end interface
197
198 abstract interface
199
200 subroutine wall_model_free(this)
201 import wall_model_t
202 class(wall_model_t), intent(inout) :: this
203 end subroutine wall_model_free
204 end interface
205
206 interface
207
216 module subroutine wall_model_factory(object, scheme_name, coef, msk, &
217 facet, json)
218 class(wall_model_t), allocatable, intent(inout) :: object
219 character(len=*), intent(in) :: scheme_name
220 type(coef_t), intent(in) :: coef
221 integer, intent(in) :: msk(:)
222 integer, intent(in) :: facet(:)
223 type(json_file), intent(inout) :: json
224 end subroutine wall_model_factory
225 end interface
226
227 interface
228
231 module subroutine wall_model_allocator(object, type_name)
232 class(wall_model_t), allocatable, intent(inout) :: object
233 character(len=:), allocatable, intent(in) :: type_name
234 end subroutine wall_model_allocator
235 end interface
236
237 !
238 ! Machinery for injecting user-defined types
239 !
240
244 abstract interface
245 subroutine wall_model_allocate(obj)
246 import wall_model_t
247 class(wall_model_t), allocatable, intent(inout) :: obj
248 end subroutine wall_model_allocate
249 end interface
250
251 interface
252
253 module subroutine register_wall_model(type_name, allocator)
254 character(len=*), intent(in) :: type_name
255 procedure(wall_model_allocate), pointer, intent(in) :: allocator
256 end subroutine register_wall_model
257 end interface
258
259 ! A name-allocator pair for user-defined types. A helper type to define a
260 ! registry of custom allocators.
261 type allocator_entry
262 character(len=20) :: type_name
263 procedure(wall_model_allocate), pointer, nopass :: allocator
264 end type allocator_entry
265
267 type(allocator_entry), allocatable :: wall_model_registry(:)
268
270 integer :: wall_model_registry_size = 0
271
272 public :: wall_model_factory, wall_model_allocator, register_wall_model, &
273 wall_model_allocate
274
275contains
282 subroutine wall_model_init_base(this, scheme_name, coef, msk, facet, index)
283 class(wall_model_t), intent(inout) :: this
284 type(coef_t), target, intent(in) :: coef
285 integer, target, intent(in) :: msk(0:)
286 integer, target, intent(in) :: facet(0:)
287 character(len=*) :: scheme_name
288 integer, intent(in) :: index
289
290 call this%free_base
291
292 this%coef => coef
293 this%dof => coef%dof
294 this%h_index = index
295 this%scheme_name = trim(scheme_name)
296 this%mu => neko_registry%get_field_by_name(this%scheme_name // "_mu")
297 this%rho => neko_registry%get_field_by_name(this%scheme_name // &
298 "_rho")
299
300 call neko_registry%add_field(this%dof, "tau", &
301 ignore_existing = .true.)
302 this%tau_field => neko_registry%get_field("tau")
303
304 call this%finalize_base(msk, facet)
305 end subroutine wall_model_init_base
306
311 subroutine wall_model_partial_init_base(this, coef, json)
312 class(wall_model_t), intent(inout) :: this
313 type(coef_t), target, intent(in) :: coef
314 type(json_file), intent(inout) :: json
315
316 call this%free_base()
317
318 this%coef => coef
319 this%dof => coef%dof
320 call json_get_or_lookup(json, "h_index", this%h_index)
321 call json_get(json, "scheme_name", this%scheme_name)
322
323 this%mu => neko_registry%get_field_by_name(this%scheme_name // "_mu")
324 this%rho => neko_registry%get_field_by_name(this%scheme_name // &
325 "_rho")
326
327 call neko_registry%add_field(this%dof, "tau", &
328 ignore_existing = .true.)
329 this%tau_field => neko_registry%get_field("tau")
330 end subroutine wall_model_partial_init_base
331
332 subroutine wall_model_finalize_base(this, msk, facet)
333 class(wall_model_t), intent(inout) :: this
334 integer, target, intent(in) :: msk(0:)
335 integer, target, intent(in) :: facet(:)
336
337 this%msk(0:msk(0)) => msk
338 if (neko_bcknd_device .eq. 1) this%msk_d = device_get_ptr(msk)
339 this%facet(0:msk(0)) => facet
340
341 call this%tau_x%init(this%msk(0))
342 call this%tau_y%init(this%msk(0))
343 call this%tau_z%init(this%msk(0))
344
345 allocate(this%ind_r(this%msk(0)))
346 allocate(this%ind_s(this%msk(0)))
347 allocate(this%ind_t(this%msk(0)))
348 allocate(this%ind_e(this%msk(0)))
349
350 call this%h%init(this%msk(0))
351 call this%n_x%init(this%msk(0))
352 call this%n_y%init(this%msk(0))
353 call this%n_z%init(this%msk(0))
354
355 call this%find_points()
356
357 ! Initialize pointers for device
358 if (neko_bcknd_device .eq. 1) then
359 call device_map(this%ind_r, this%ind_r_d, this%n_nodes)
360 call device_map(this%ind_s, this%ind_s_d, this%n_nodes)
361 call device_map(this%ind_t, this%ind_t_d, this%n_nodes)
362 call device_map(this%ind_e, this%ind_e_d, this%n_nodes)
363 call device_memcpy(this%ind_r, this%ind_r_d, this%n_nodes, &
364 host_to_device, sync = .false.)
365 call device_memcpy(this%ind_s, this%ind_s_d, this%n_nodes, &
366 host_to_device, sync = .false.)
367 call device_memcpy(this%ind_t, this%ind_t_d, this%n_nodes, &
368 host_to_device, sync = .false.)
369 call device_memcpy(this%ind_e, this%ind_e_d, this%n_nodes, &
370 host_to_device, sync = .false.)
371 end if
372
373 end subroutine wall_model_finalize_base
374
376 subroutine wall_model_free_base(this)
377 class(wall_model_t), intent(inout) :: this
378
379 nullify(this%coef)
380 nullify(this%msk)
381 nullify(this%facet)
382 nullify(this%tau_field)
383 nullify(this%mu)
384 nullify(this%rho)
385 this%msk_d = c_null_ptr
386
387 call this%tau_x%free()
388 call this%tau_y%free()
389 call this%tau_z%free()
390
391
392 if (allocated(this%ind_r)) then
393 if (neko_bcknd_device .eq. 1) then
394 call device_unmap(this%ind_r, this%ind_r_d)
395 end if
396 deallocate(this%ind_r)
397 end if
398 if (allocated(this%ind_s)) then
399 if (neko_bcknd_device .eq. 1) then
400 call device_unmap(this%ind_s, this%ind_s_d)
401 end if
402 deallocate(this%ind_s)
403 end if
404 if (allocated(this%ind_t)) then
405 if (neko_bcknd_device .eq. 1) then
406 call device_unmap(this%ind_t, this%ind_t_d)
407 end if
408 deallocate(this%ind_t)
409 end if
410 if (allocated(this%ind_e)) then
411 if (neko_bcknd_device .eq. 1) then
412 call device_unmap(this%ind_e, this%ind_e_d)
413 end if
414 deallocate(this%ind_e)
415 end if
416
417
418 if (allocated(this%scheme_name)) then
419 deallocate(this%scheme_name)
420 end if
421
422
423 call this%h%free()
424 call this%n_x%free()
425 call this%n_y%free()
426 call this%n_z%free()
427
428 nullify(this%dof)
429 end subroutine wall_model_free_base
430
432 subroutine wall_model_find_points(this)
433 class(wall_model_t), intent(inout) :: this
434 integer :: n_nodes, fid, idx(4), i, linear
435 real(kind=rp) :: normal(3), p(3), x, y, z, xw, yw, zw, magp
436 real(kind=rp) :: hmin, hmax
437 type(field_t), pointer :: h_field
438 type(file_t) :: h_file
439 character(len=LOG_SIZE), allocatable :: log_msg
440
441 n_nodes = this%msk(0)
442 this%n_nodes = n_nodes
443
444 call neko_registry%add_field(this%coef%dof, "sampling_height", &
445 ignore_existing = .true.)
446
447 h_field => neko_registry%get_field_by_name("sampling_height")
448
449 do i = 1, n_nodes
450 linear = this%msk(i)
451 fid = this%facet(i)
452 idx = nonlinear_index(linear, this%coef%Xh%lx, this%coef%Xh%ly,&
453 this%coef%Xh%lz)
454 normal = this%coef%get_normal(idx(1), idx(2), idx(3), idx(4), fid)
455
456 this%n_x%x(i) = normal(1)
457 this%n_y%x(i) = normal(2)
458 this%n_z%x(i) = normal(3)
459
460 ! inward normal
461 normal = -normal
462
463 select case (fid)
464 case (1)
465 this%ind_r(i) = idx(1) + this%h_index
466 this%ind_s(i) = idx(2)
467 this%ind_t(i) = idx(3)
468 case (2)
469 this%ind_r(i) = idx(1) - this%h_index
470 this%ind_s(i) = idx(2)
471 this%ind_t(i) = idx(3)
472 case (3)
473 this%ind_r(i) = idx(1)
474 this%ind_s(i) = idx(2) + this%h_index
475 this%ind_t(i) = idx(3)
476 case (4)
477 this%ind_r(i) = idx(1)
478 this%ind_s(i) = idx(2) - this%h_index
479 this%ind_t(i) = idx(3)
480 case (5)
481 this%ind_r(i) = idx(1)
482 this%ind_s(i) = idx(2)
483 this%ind_t(i) = idx(3) + this%h_index
484 case (6)
485 this%ind_r(i) = idx(1)
486 this%ind_s(i) = idx(2)
487 this%ind_t(i) = idx(3) - this%h_index
488 case default
489 call neko_error("The face index is not correct ")
490 end select
491 this%ind_e(i) = idx(4)
492
493 ! Location of the wall node
494 xw = this%dof%x(idx(1), idx(2), idx(3), idx(4))
495 yw = this%dof%y(idx(1), idx(2), idx(3), idx(4))
496 zw = this%dof%z(idx(1), idx(2), idx(3), idx(4))
497
498 ! Location of the sampling point
499 x = this%dof%x(this%ind_r(i), this%ind_s(i), this%ind_t(i), &
500 this%ind_e(i))
501 y = this%dof%y(this%ind_r(i), this%ind_s(i), this%ind_t(i), &
502 this%ind_e(i))
503 z = this%dof%z(this%ind_r(i), this%ind_s(i), this%ind_t(i), &
504 this%ind_e(i))
505
506 ! Vector from the sampling point to the wall
507 p(1) = x - xw
508 p(2) = y - yw
509 p(3) = z - zw
510
511 ! Total distance to the sampling point
512 magp = sqrt(p(1)**2 + p(2)**2 + p(3)**2)
513
514 ! Project on the normal direction to get h
515 this%h%x(i) = p(1)*normal(1) + p(2)*normal(2) + p(3)*normal(3)
516
517 h_field%x(linear,1,1,1) = this%h%x(i)
518
519 ! Look at how much the total distance distance from the normal and warn
520 ! if significant
521 if ((this%h%x(i) - magp) / magp > 0.1) then
522 write(log_msg,*) "Significant misalignment between wall normal and"
523 call neko_log%message(log_msg, neko_log_debug)
524 write(log_msg,*) "sampling point direction at wall node", xw, yw, zw
525 call neko_log%message(log_msg, neko_log_debug)
526 end if
527 end do
528
529! hmin = glmin(this%h%x, n_nodes)
530! hmax = glmax(this%h%x, n_nodes)
531! if (pe_rank .eq. 0) then
532! write(*, "(A, F10.4, F10.4)") " h min / max:", hmin, hmax
533! end if
534
535 if (neko_bcknd_device .eq. 1) then
536 call device_memcpy(this%h%x, this%h%x_d, n_nodes, host_to_device,&
537 sync = .false.)
538 call device_memcpy(this%n_x%x, this%n_x%x_d, n_nodes, host_to_device, &
539 sync = .false.)
540 call device_memcpy(this%n_y%x, this%n_y%x_d, n_nodes, host_to_device, &
541 sync = .false.)
542 call device_memcpy(this%n_z%x, this%n_z%x_d, n_nodes, host_to_device, &
543 sync = .true.)
544 end if
545
546 ! Each wall_model bc will do a write unfortunately... But very helpful
547 ! for setup debugging.
548 call h_file%init("sampling_height.fld")
549 call h_file%write(h_field)
550 end subroutine wall_model_find_points
551
553 class(wall_model_t), intent(inout) :: this
554 integer :: i, m
555 real(kind=rp) :: magtau
556
557 m = this%msk(0)
558 if (m > 0) then
559 if (neko_bcknd_device .eq. 1) then
560 call wall_model_compute_mag_field_device(this%tau_x%x_d, &
561 this%tau_y%x_d, &
562 this%tau_z%x_d, &
563 this%tau_field%x_d, &
564 this%msk_d, m)
565 else
566 do i = 1, m
567 magtau = sqrt(this%tau_x%x(i)**2 + &
568 this%tau_y%x(i)**2 + &
569 this%tau_z%x(i)**2)
570 this%tau_field%x(this%msk(i),1,1,1) = magtau
571 end do
572 end if
573 end if
574
575 end subroutine wall_model_compute_mag_field
576
577end module wall_model
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
__global__ void wall_model_compute_mag_field(const T *__restrict__ tau_x_d, const T *__restrict__ tau_y_d, const T *__restrict__ tau_z_d, T *__restrict__ tau_field_d, const int *__restrict__ msk_d, const int m)
Return the device pointer for an associated Fortran array.
Definition device.F90:108
Map a Fortran array to a device (allocate and associate)
Definition device.F90:78
Copy data between host and device (or device and device)
Definition device.F90:72
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:84
Retrieves a parameter by name or throws an error.
Compute wall shear stress.
Finilzation of partial construction, similar to bc_t
Common constructor.
Partial constructor from JSON, meant to work as the first stage of initialization before the finalize...
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
integer, public pe_rank
MPI rank.
Definition comm.F90:58
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:48
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:90
type(log_t), public neko_log
Global log stream.
Definition log.f90:80
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:648
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:686
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
Utilities.
Definition utils.f90:35
Defines a vector.
Definition vector.f90:34
Implements the device kernel for the wall_model_t type.
subroutine, public wall_model_compute_mag_field_device(tau_x_d, tau_y_d, tau_z_d, tau_field_d, msk_d, m)
Compute the wall shear stress's magnitude on device.
Implements wall_model_t.
subroutine wall_model_finalize_base(this, msk, facet)
subroutine wall_model_init_base(this, scheme_name, coef, msk, facet, index)
Constructor for the wall_model_t (base) class.
subroutine wall_model_partial_init_base(this, coef, json)
Partial initialization based on JSON, prior to knowing the mask and facets.
subroutine wall_model_find_points(this)
Find sampling points based on the requested index.
subroutine wall_model_free_base(this)
Destructor for the wall_model_t (base) class.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:56
Base abstract type for wall-stress models for wall-modelled LES.