Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
gs_cpu.f90
Go to the documentation of this file.
1! Copyright (c) 2020-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!
34module gs_cpu
35 use num_types, only : rp
36 use gs_bcknd, only : gs_bcknd_t
38 use, intrinsic :: iso_c_binding, only : c_ptr
39 implicit none
40 private
41
43 type, public, extends(gs_bcknd_t) :: gs_cpu_t
44 contains
45 procedure, pass(this) :: init => gs_cpu_init
46 procedure, pass(this) :: free => gs_cpu_free
47 procedure, pass(this) :: gather => gs_gather_cpu
48 procedure, pass(this) :: scatter => gs_scatter_cpu
49 end type gs_cpu_t
50
51contains
52
54 subroutine gs_cpu_init(this, nlocal, nshared, nlcl_blks, nshrd_blks)
55 class(gs_cpu_t), intent(inout) :: this
56 integer, intent(in) :: nlocal
57 integer, intent(in) :: nshared
58 integer, intent(in) :: nlcl_blks
59 integer, intent(in) :: nshrd_blks
60 end subroutine gs_cpu_init
61
63 subroutine gs_cpu_free(this)
64 class(gs_cpu_t), intent(inout) :: this
65 end subroutine gs_cpu_free
66
68 subroutine gs_gather_cpu(this, v, m, o, dg, u, n, gd, nb, b, bo, op, shrd)
69 integer, intent(in) :: m
70 integer, intent(in) :: n
71 integer, intent(in) :: nb
72 class(gs_cpu_t), intent(inout) :: this
73 real(kind=rp), dimension(m), intent(inout) :: v
74 integer, dimension(m), intent(inout) :: dg
75 real(kind=rp), dimension(n), intent(inout) :: u
76 integer, dimension(m), intent(inout) :: gd
77 integer, dimension(nb), intent(inout) :: b
78 integer, dimension(nb), intent(inout) :: bo
79 integer, intent(in) :: o
80 integer, intent(in) :: op
81 logical, intent(in) :: shrd
82
83 select case(op)
84 case (gs_op_add)
85 call gs_gather_kernel_add(v, m, o, dg, u, n, gd, nb, b, bo)
86 case (gs_op_mul)
87 call gs_gather_kernel_mul(v, m, o, dg, u, n, gd, nb, b, bo)
88 case (gs_op_min)
89 call gs_gather_kernel_min(v, m, o, dg, u, n, gd, nb, b, bo)
90 case (gs_op_max)
91 call gs_gather_kernel_max(v, m, o, dg, u, n, gd, nb, b, bo)
92 end select
93
94 end subroutine gs_gather_cpu
95
98 subroutine gs_gather_kernel_add(v, m, o, dg, u, n, gd, nb, b, bo)
99 integer, intent(in) :: m
100 integer, intent(in) :: n
101 integer, intent(in) :: nb
102 real(kind=rp), dimension(m), intent(inout) :: v
103 integer, dimension(m), intent(inout) :: dg
104 real(kind=rp), dimension(n), intent(inout) :: u
105 integer, dimension(m), intent(inout) :: gd
106 integer, dimension(nb), intent(inout) :: b
107 integer, dimension(nb), intent(inout) :: bo
108 integer, intent(in) :: o
109 integer :: i, j, k, blk_len
110 real(kind=rp) :: tmp
111
112 !$omp do
113 do i = 1, nb
114 k = bo(i)
115 blk_len = b(i)
116 tmp = u(gd(k + 1))
117 do j = 2, blk_len
118 tmp = tmp + u(gd(k + j))
119 end do
120 v(dg(k + 1)) = tmp
121 end do
122 !$omp end do nowait
123
124 if (o .lt. 0) then
125 !$omp do
126 do i = abs(o), m
127 v(dg(i)) = u(gd(i))
128 end do
129 !$omp end do
130 else
131 !$omp do
132 do i = o, m, 2
133 v(dg(i)) = u(gd(i)) + u(gd(i+1))
134 end do
135 !$omp end do
136 end if
137
138 end subroutine gs_gather_kernel_add
139
142 subroutine gs_gather_kernel_mul(v, m, o, dg, u, n, gd, nb, b, bo)
143 integer, intent(in) :: m
144 integer, intent(in) :: n
145 integer, intent(in) :: nb
146 real(kind=rp), dimension(m), intent(inout) :: v
147 integer, dimension(m), intent(inout) :: dg
148 real(kind=rp), dimension(n), intent(inout) :: u
149 integer, dimension(m), intent(inout) :: gd
150 integer, dimension(nb), intent(inout) :: b
151 integer, dimension(nb), intent(inout) :: bo
152 integer, intent(in) :: o
153 integer :: i, j, k, blk_len
154 real(kind=rp) :: tmp
155
156 !$omp do
157 do i = 1, nb
158 k = bo(i)
159 blk_len = b(i)
160 tmp = u(gd(k + 1))
161 do j = 2, blk_len
162 tmp = tmp * u(gd(k + j))
163 end do
164 v(dg(k + 1)) = tmp
165 end do
166 !$omp end do nowait
167
168 if (o .lt. 0) then
169 !$omp do
170 do i = abs(o), m
171 v(dg(i)) = u(gd(i))
172 end do
173 !$omp end do nowait
174 else
175 !$omp do
176 do i = o, m, 2
177 v(dg(i)) = u(gd(i)) * u(gd(i+1))
178 end do
179 !$omp end do nowait
180 end if
181
182 end subroutine gs_gather_kernel_mul
183
186 subroutine gs_gather_kernel_min(v, m, o, dg, u, n, gd, nb, b, bo)
187 integer, intent(in) :: m
188 integer, intent(in) :: n
189 integer, intent(in) :: nb
190 real(kind=rp), dimension(m), intent(inout) :: v
191 integer, dimension(m), intent(inout) :: dg
192 real(kind=rp), dimension(n), intent(inout) :: u
193 integer, dimension(m), intent(inout) :: gd
194 integer, dimension(nb), intent(inout) :: b
195 integer, dimension(nb), intent(inout) :: bo
196 integer, intent(in) :: o
197 integer :: i, j, k, blk_len
198 real(kind=rp) :: tmp
199
200 !$omp do
201 do i = 1, nb
202 k = bo(i)
203 blk_len = b(i)
204 tmp = u(gd(k + 1))
205 do j = 2, blk_len
206 tmp = min(tmp, u(gd(k + j)))
207 end do
208 v(dg(k + 1)) = tmp
209 end do
210 !$omp end do nowait
211
212 if (o .lt. 0) then
213 !$omp do
214 do i = abs(o), m
215 v(dg(i)) = u(gd(i))
216 end do
217 !$omp end do
218 else
219 !$omp do
220 do i = o, m, 2
221 v(dg(i)) = min(u(gd(i)), u(gd(i+1)))
222 end do
223 !$omp end do
224 end if
225
226 end subroutine gs_gather_kernel_min
227
230 subroutine gs_gather_kernel_max(v, m, o, dg, u, n, gd, nb, b, bo)
231 integer, intent(in) :: m
232 integer, intent(in) :: n
233 integer, intent(in) :: nb
234 real(kind=rp), dimension(m), intent(inout) :: v
235 integer, dimension(m), intent(inout) :: dg
236 real(kind=rp), dimension(n), intent(inout) :: u
237 integer, dimension(m), intent(inout) :: gd
238 integer, dimension(nb), intent(inout) :: b
239 integer, dimension(nb), intent(inout) :: bo
240 integer, intent(in) :: o
241 integer :: i, j, k, blk_len
242 real(kind=rp) :: tmp
243
244 !$omp do
245 do i = 1, nb
246 k = bo(i)
247 blk_len = b(i)
248 tmp = u(gd(k + 1))
249 do j = 2, blk_len
250 tmp = max(tmp, u(gd(k + j)))
251 end do
252 v(dg(k + 1)) = tmp
253 end do
254 !$omp end do nowait
255
256 if (o .lt. 0) then
257 !$omp do
258 do i = abs(o), m
259 v(dg(i)) = u(gd(i))
260 end do
261 !$omp end do
262 else
263 !$omp do
264 do i = o, m, 2
265 v(dg(i)) = max(u(gd(i)), u(gd(i+1)))
266 end do
267 !$omp end do
268 end if
269
270 end subroutine gs_gather_kernel_max
271
273 subroutine gs_scatter_cpu(this, v, m, dg, u, n, gd, nb, b, bo, shrd, event)
274 integer, intent(in) :: m
275 integer, intent(in) :: n
276 integer, intent(in) :: nb
277 class(gs_cpu_t), intent(inout) :: this
278 real(kind=rp), dimension(m), intent(inout) :: v
279 integer, dimension(m), intent(inout) :: dg
280 real(kind=rp), dimension(n), intent(inout) :: u
281 integer, dimension(m), intent(inout) :: gd
282 integer, dimension(nb), intent(inout) :: b
283 integer, dimension(nb), intent(inout) :: bo
284 logical, intent(in) :: shrd
285 type(c_ptr) :: event
286
287 call gs_scatter_kernel(v, m, dg, u, n, gd, nb, b, bo)
288
289 end subroutine gs_scatter_cpu
290
292 subroutine gs_scatter_kernel(v, m, dg, u, n, gd, nb, b, bo)
293 integer, intent(in) :: m
294 integer, intent(in) :: n
295 integer, intent(in) :: nb
296 real(kind=rp), dimension(m), intent(inout) :: v
297 integer, dimension(m), intent(inout) :: dg
298 real(kind=rp), dimension(n), intent(inout) :: u
299 integer, dimension(m), intent(inout) :: gd
300 integer, dimension(nb), intent(inout) :: b
301 integer, dimension(nb), intent(inout) :: bo
302 integer :: i, j, k, blk_len
303 real(kind=rp) :: tmp
304
305 !$omp do
306 do i = 1, nb
307 k = bo(i)
308 blk_len = b(i)
309 tmp = v(dg(k + 1))
310 do j = 1, blk_len
311 u(gd(k + j)) = tmp
312 end do
313 end do
314 !$omp end do nowait
315
316 !$omp do
317 do i = (bo(nb) + b(nb) + 1), m
318 u(gd(i)) = v(dg(i))
319 end do
320 !$omp end do
321
322 end subroutine gs_scatter_kernel
323
324end module gs_cpu
Defines a gather-scatter backend.
Definition gs_bcknd.f90:34
Generic Gather-scatter backend for CPUs.
Definition gs_cpu.f90:34
subroutine gs_scatter_kernel(v, m, dg, u, n, gd, nb, b, bo)
Scatter kernel .
Definition gs_cpu.f90:293
subroutine gs_cpu_init(this, nlocal, nshared, nlcl_blks, nshrd_blks)
Dummy backend initialisation.
Definition gs_cpu.f90:55
subroutine gs_gather_kernel_add(v, m, o, dg, u, n, gd, nb, b, bo)
Gather kernel for addition of data .
Definition gs_cpu.f90:99
subroutine gs_gather_kernel_mul(v, m, o, dg, u, n, gd, nb, b, bo)
Gather kernel for multiplication of data .
Definition gs_cpu.f90:143
subroutine gs_cpu_free(this)
Dummy backend deallocation.
Definition gs_cpu.f90:64
subroutine gs_scatter_cpu(this, v, m, dg, u, n, gd, nb, b, bo, shrd, event)
Scatter kernel.
Definition gs_cpu.f90:274
subroutine gs_gather_cpu(this, v, m, o, dg, u, n, gd, nb, b, bo, op, shrd)
Gather kernel.
Definition gs_cpu.f90:69
subroutine gs_gather_kernel_min(v, m, o, dg, u, n, gd, nb, b, bo)
Gather kernel for minimum of data .
Definition gs_cpu.f90:187
subroutine gs_gather_kernel_max(v, m, o, dg, u, n, gd, nb, b, bo)
Gather kernel for maximum of data .
Definition gs_cpu.f90:231
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
integer, parameter, public gs_op_max
Definition gs_ops.f90:36
integer, parameter, public gs_op_min
Definition gs_ops.f90:36
integer, parameter, public gs_op_mul
Definition gs_ops.f90:36
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Gather-scatter backend.
Definition gs_bcknd.f90:44
Gather-scatter backend for CPUs.
Definition gs_cpu.f90:43
#define max(a, b)
Definition tensor.cu:40