Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
spalding_device.F90
Go to the documentation of this file.
1! Copyright (c) 2025, 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, c_rp
36 use, intrinsic :: iso_c_binding, only : c_ptr
37 use utils, only : neko_error
38 implicit none
39 private
40
41#ifdef HAVE_HIP
42 interface
43 subroutine hip_spalding_compute(u_d, v_d, w_d, &
44 ind_r_d, ind_s_d, ind_t_d, ind_e_d, &
45 n_x_d, n_y_d, n_z_d, nu_d, h_d, &
46 tau_x_d, tau_y_d, tau_z_d, n_nodes, lx, &
47 kappa, B, tstep) &
48 bind(c, name = 'hip_spalding_compute')
49 use, intrinsic :: iso_c_binding, only : c_ptr, c_int
50 use num_types, only : c_rp
51 implicit none
52 type(c_ptr), value :: u_d, v_d, w_d
53 type(c_ptr), value :: ind_r_d, ind_s_d, ind_t_d, ind_e_d
54 type(c_ptr), value :: n_x_d, n_y_d, n_z_d, h_d, nu_d
55 real(c_rp) :: kappa, B
56 type(c_ptr), value :: tau_x_d, tau_y_d, tau_z_d
57 integer(c_int) :: n_nodes, lx, tstep
58 end subroutine hip_spalding_compute
59 end interface
60#elif HAVE_CUDA
61 interface
62 subroutine cuda_spalding_compute(u_d, v_d, w_d, &
63 ind_r_d, ind_s_d, ind_t_d, ind_e_d, &
64 n_x_d, n_y_d, n_z_d, nu_d, h_d, &
65 tau_x_d, tau_y_d, tau_z_d, n_nodes, lx, &
66 kappa, B, tstep) &
67 bind(c, name = 'cuda_spalding_compute')
68 use, intrinsic :: iso_c_binding, only : c_ptr, c_int
69 use num_types, only : c_rp
70 implicit none
71 type(c_ptr), value :: u_d, v_d, w_d
72 type(c_ptr), value :: ind_r_d, ind_s_d, ind_t_d, ind_e_d
73 type(c_ptr), value :: n_x_d, n_y_d, n_z_d, h_d, nu_d
74 real(c_rp) :: kappa, B
75 type(c_ptr), value :: tau_x_d, tau_y_d, tau_z_d
76 integer(c_int) :: n_nodes, lx, tstep
77 end subroutine cuda_spalding_compute
78 end interface
79#elif HAVE_OPENCL
80#endif
82
83contains
87 subroutine spalding_compute_device(u_d, v_d, w_d, &
88 ind_r_d, ind_s_d, ind_t_d, ind_e_d, &
89 n_x_d, n_y_d, n_z_d, nu_d, h_d, tau_x_d, tau_y_d, tau_z_d, &
90 n_nodes, lx, kappa, B, tstep)
91 integer, intent(in) :: n_nodes, lx, tstep
92 type(c_ptr), intent(in) :: u_d, v_d, w_d
93 type(c_ptr), intent(in) :: ind_r_d, ind_s_d, ind_t_d, ind_e_d
94 type(c_ptr), intent(in) :: n_x_d, n_y_d, n_z_d, h_d, nu_d
95 type(c_ptr), intent(inout) :: tau_x_d, tau_y_d, tau_z_d
96 real(kind=rp), intent(in) :: kappa, b
97
98#if HAVE_HIP
99 call hip_spalding_compute(u_d, v_d, w_d, &
100 ind_r_d, ind_s_d, ind_t_d, ind_e_d, &
101 n_x_d, n_y_d, n_z_d, nu_d, h_d, &
102 tau_x_d, tau_y_d, tau_z_d, n_nodes, lx, kappa, b, tstep)
103#elif HAVE_CUDA
104 call cuda_spalding_compute(u_d, v_d, w_d, &
105 ind_r_d, ind_s_d, ind_t_d, ind_e_d, &
106 n_x_d, n_y_d, n_z_d, nu_d, h_d, &
107 tau_x_d, tau_y_d, tau_z_d, n_nodes, lx, kappa, b, tstep)
108#elif HAVE_OPENCL
109 call neko_error("OPENCL is not implemented for Spalding's model")
110#else
111 call neko_error('No device backend configured')
112#endif
113
114 end subroutine spalding_compute_device
115end module spalding_device
integer, parameter, public c_rp
Definition num_types.f90:13
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements the device kernel for the spalding_t type.
subroutine, public spalding_compute_device(u_d, v_d, w_d, ind_r_d, ind_s_d, ind_t_d, ind_e_d, n_x_d, n_y_d, n_z_d, nu_d, h_d, tau_x_d, tau_y_d, tau_z_d, n_nodes, lx, kappa, b, tstep)
Compute the wall shear stress on device using Spalding's model.
Utilities.
Definition utils.f90:35
void cuda_spalding_compute(void *u_d, void *v_d, void *w_d, void *ind_r_d, void *ind_s_d, void *ind_t_d, void *ind_e_d, void *n_x_d, void *n_y_d, void *n_z_d, void *nu_d, void *h_d, void *tau_x_d, void *tau_y_d, void *tau_z_d, int *n_nodes, int *lx, real *kappa, real *B, int *tstep)
Definition spalding.cu:43