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