Neko  0.8.1
A portable framework for high-order spectral element flow simulations
scalar_residual_device.F90
Go to the documentation of this file.
1 ! Copyright (c) 2022, 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 !
34  use scalar_residual
35  use gather_scatter
36  use operators
37  use device_math
38  use device_mathops
39  use, intrinsic :: iso_c_binding
40  implicit none
41  private
42 
43  type, public, extends(scalar_residual_t) :: scalar_residual_device_t
44  contains
45  procedure, nopass :: compute => scalar_residual_device_compute
47 
48 #ifdef HAVE_HIP
49  interface
50  subroutine scalar_residual_update_hip(s_res_d,f_s_d, n) &
51  bind(c, name='scalar_residual_update_hip')
52  use, intrinsic :: iso_c_binding
53  implicit none
54  type(c_ptr), value :: s_res_d
55  type(c_ptr), value :: f_s_d
56  integer(c_int) :: n
57  end subroutine scalar_residual_update_hip
58  end interface
59 #elif HAVE_CUDA
60 
61  interface
62  subroutine scalar_residual_update_cuda(s_res_d,f_s_d, n) &
63  bind(c, name='scalar_residual_update_cuda')
64  use, intrinsic :: iso_c_binding
65  implicit none
66  type(c_ptr), value :: s_res_d
67  type(c_ptr), value :: f_s_d
68  integer(c_int) :: n
69  end subroutine scalar_residual_update_cuda
70  end interface
71 #elif HAVE_OPENCL
72 
73  interface
74  subroutine scalar_residual_update_opencl(s_res_d,f_s_d, n) &
75  bind(c, name='scalar_residual_update_opencl')
76  use, intrinsic :: iso_c_binding
77  implicit none
78  type(c_ptr), value :: s_res_d
79  type(c_ptr), value :: f_s_d
80  integer(c_int) :: n
81  end subroutine scalar_residual_update_opencl
82  end interface
83 #endif
84 
85 
86 contains
87 
88 
89  subroutine scalar_residual_device_compute(Ax, s, s_res, f_Xh, c_Xh, msh, Xh, &
90  lambda, rhocp, bd, dt, n)
91  class(ax_t), intent(in) :: Ax
92  type(mesh_t), intent(inout) :: msh
93  type(space_t), intent(inout) :: Xh
94  type(field_t), intent(inout) :: s
95  type(field_t), intent(inout) :: s_res
96  type(field_t), intent(inout) :: f_Xh
97  type(coef_t), intent(inout) :: c_Xh
98  real(kind=rp), intent(in) :: lambda
99  real(kind=rp), intent(in) :: rhocp
100  real(kind=rp), intent(in) :: bd
101  real(kind=rp), intent(in) :: dt
102  integer, intent(in) :: n
103 
104  call device_cfill(c_xh%h1_d, lambda, n)
105  call device_cfill(c_xh%h2_d, rhocp * (bd / dt), n)
106  c_xh%ifh2 = .true.
107 
108  call ax%compute(s_res%x, s%x, c_xh, msh, xh)
109 
110 #ifdef HAVE_HIP
111  call scalar_residual_update_hip(s_res%x_d, f_xh%x_d, n)
112 #elif HAVE_CUDA
113  call scalar_residual_update_cuda(s_res%x_d, f_xh%x_d, n)
114 #elif HAVE_OPENCL
115  call scalar_residual_update_opencl(s_res%x_d, f_xh%x_d, n)
116 #endif
117 
118  end subroutine scalar_residual_device_compute
119 
120 end module scalar_residual_device
subroutine, public device_cfill(a_d, c, n)
Gather-scatter.
Operators.
Definition: operators.f90:34
subroutine scalar_residual_device_compute(Ax, s, s_res, f_Xh, c_Xh, msh, Xh, lambda, rhocp, bd, dt, n)
Defines the residual for the scalar transport equation.
void scalar_residual_update_opencl(void *s_res, void *f_s, int *n)
void scalar_residual_update_cuda(void *s_res, void *f_s, int *n)
Abstract type to compute scalar residual.