Neko 1.99.2
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, bo, 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, dimension(nb), intent(inout) :: bo
104 integer, intent(in) :: o
105 integer, intent(in) :: op
106 logical, intent(in) :: shrd
107
108 if (.not. shrd) then
109 associate(w=>this%local_wrk)
110 select case(op)
111 case (gs_op_add)
112 call gs_gather_kernel_add(v, m, o, dg, u, n, gd, nb, b, w)
113 case (gs_op_mul)
114 call gs_gather_kernel_mul(v, m, o, dg, u, n, gd, nb, b, w)
115 case (gs_op_min)
116 call gs_gather_kernel_min(v, m, o, dg, u, n, gd, nb, b, w)
117 case (gs_op_max)
118 call gs_gather_kernel_max(v, m, o, dg, u, n, gd, nb, b, w)
119 end select
120 end associate
121 else if (shrd) then
122 associate(w=>this%shared_wrk)
123 select case(op)
124 case (gs_op_add)
125 call gs_gather_kernel_add(v, m, o, dg, u, n, gd, nb, b, w)
126 case (gs_op_mul)
127 call gs_gather_kernel_mul(v, m, o, dg, u, n, gd, nb, b, w)
128 case (gs_op_min)
129 call gs_gather_kernel_min(v, m, o, dg, u, n, gd, nb, b, w)
130 case (gs_op_max)
131 call gs_gather_kernel_max(v, m, o, dg, u, n, gd, nb, b, w)
132 end select
133 end associate
134 end if
135
136 end subroutine gs_gather_sx
137
140 subroutine gs_gather_kernel_add(v, m, o, dg, u, n, gd, nb, b, w)
141 integer, intent(in) :: m
142 integer, intent(in) :: n
143 integer, intent(in) :: nb
144 real(kind=rp), dimension(m), intent(inout) :: v
145 real(kind=rp), dimension(m), intent(inout) :: w
146 integer, dimension(m), intent(inout) :: dg
147 real(kind=rp), dimension(n), intent(inout) :: u
148 integer, dimension(m), intent(inout) :: gd
149 integer, dimension(nb), intent(inout) :: b
150 integer, intent(in) :: o
151 integer :: i
152 real(kind=rp) :: tmp
153
154 v = 0d0
155 do i = 1, abs(o) - 1
156 w(i) = u(gd(i))
157 end do
158
159 do i = 1, abs(o) - 1
160 v(dg(i)) = v(dg(i)) + w(i)
161 end do
162
163 if (o .lt. 0) then
164 do i = abs(o), m
165 v(dg(i)) = u(gd(i))
166 end do
167 else
168 do i = o, m, 2
169 tmp = u(gd(i)) + u(gd(i+1))
170 v(dg(i)) = tmp
171 end do
172 end if
173
174 end subroutine gs_gather_kernel_add
175
178 subroutine gs_gather_kernel_mul(v, m, o, dg, u, n, gd, nb, b, w)
179 integer, intent(in) :: m
180 integer, intent(in) :: n
181 integer, intent(in) :: nb
182 real(kind=rp), dimension(m), intent(inout) :: v
183 real(kind=rp), dimension(m), intent(inout) :: w
184 integer, dimension(m), intent(inout) :: dg
185 real(kind=rp), dimension(n), intent(inout) :: u
186 integer, dimension(m), intent(inout) :: gd
187 integer, dimension(nb), intent(inout) :: b
188 integer, intent(in) :: o
189 integer :: i
190 real(kind=rp) :: tmp
191
192 do i = 1, abs(o) - 1
193 w(i) = u(gd(i))
194 end do
195
196 do i = 1, abs(o) - 1
197 v(dg(i)) = v(dg(i)) * w(i)
198 end do
199
200 if (o .lt. 0) then
201 do i = abs(o), m
202 v(dg(i)) = u(gd(i))
203 end do
204 else
205 do i = o, m, 2
206 tmp = u(gd(i)) * u(gd(i+1))
207 v(dg(i)) = tmp
208 end do
209 end if
210
211 end subroutine gs_gather_kernel_mul
212
215 subroutine gs_gather_kernel_min(v, m, o, dg, u, n, gd, nb, b, w)
216 integer, intent(in) :: m
217 integer, intent(in) :: n
218 integer, intent(in) :: nb
219 real(kind=rp), dimension(m), intent(inout) :: v
220 real(kind=rp), dimension(m), intent(inout) :: w
221 integer, dimension(m), intent(inout) :: dg
222 real(kind=rp), dimension(n), intent(inout) :: u
223 integer, dimension(m), intent(inout) :: gd
224 integer, dimension(nb), intent(inout) :: b
225 integer, intent(in) :: o
226 integer :: i
227 real(kind=rp) :: tmp
228
229 do i = 1, abs(o) - 1
230 w(i) = u(gd(i))
231 end do
232
233 do i = 1, abs(o) - 1
234 v(dg(i)) = min(v(dg(i)), w(i))
235 end do
236
237 if (o .lt. 0) then
238 do i = abs(o), m
239 v(dg(i)) = u(gd(i))
240 end do
241 else
242 do i = o, m, 2
243 tmp = min(u(gd(i)), u(gd(i+1)))
244 v(dg(i)) = tmp
245 end do
246 end if
247
248 end subroutine gs_gather_kernel_min
249
252 subroutine gs_gather_kernel_max(v, m, o, dg, u, n, gd, nb, b, w)
253 integer, intent(in) :: m
254 integer, intent(in) :: n
255 integer, intent(in) :: nb
256 real(kind=rp), dimension(m), intent(inout) :: v
257 real(kind=rp), dimension(m), intent(inout) :: w
258 integer, dimension(m), intent(inout) :: dg
259 real(kind=rp), dimension(n), intent(inout) :: u
260 integer, dimension(m), intent(inout) :: gd
261 integer, dimension(nb), intent(inout) :: b
262 integer, intent(in) :: o
263 integer :: i
264 real(kind=rp) :: tmp
265
266 do i = 1, abs(o) - 1
267 w(i) = u(gd(i))
268 end do
269
270 do i = 1, abs(o) - 1
271 v(dg(i)) = max(v(dg(i)), w(i))
272 end do
273
274 if (o .lt. 0) then
275 do i = abs(o), m
276 v(dg(i)) = u(gd(i))
277 end do
278 else
279 do i = o, m, 2
280 tmp = max(u(gd(i)), u(gd(i+1)))
281 v(dg(i)) = tmp
282 end do
283 end if
284
285 end subroutine gs_gather_kernel_max
286
288 subroutine gs_scatter_sx(this, v, m, dg, u, n, gd, nb, b, bo, shrd, event)
289 integer, intent(in) :: m
290 integer, intent(in) :: n
291 integer, intent(in) :: nb
292 class(gs_sx_t), intent(inout) :: this
293 real(kind=rp), dimension(m), intent(inout) :: v
294 integer, dimension(m), intent(inout) :: dg
295 real(kind=rp), dimension(n), intent(inout) :: u
296 integer, dimension(m), intent(inout) :: gd
297 integer, dimension(nb), intent(inout) :: b
298 integer, dimension(nb), intent(inout) :: bo
299 logical, intent(in) :: shrd
300 type(c_ptr) :: event
301
302 if (.not. shrd) then
303 call gs_scatter_kernel(v, m, dg, u, n, gd, nb, b, this%local_wrk)
304 else if (shrd) then
305 call gs_scatter_kernel(v, m, dg, u, n, gd, nb, b, this%shared_wrk)
306 end if
307
308 end subroutine gs_scatter_sx
309
311 subroutine gs_scatter_kernel(v, m, dg, u, n, gd, nb, b, w)
312 integer, intent(in) :: m
313 integer, intent(in) :: n
314 integer, intent(in) :: nb
315 real(kind=rp), dimension(m), intent(inout) :: v
316 real(kind=rp), dimension(m), intent(inout) :: w
317 integer, dimension(m), intent(inout) :: dg
318 real(kind=rp), dimension(n), intent(inout) :: u
319 integer, dimension(m), intent(inout) :: gd
320 integer, dimension(nb), intent(inout) :: b
321 integer :: i
322
323 !NEC$ IVDEP
324 do i = 1, m
325 w(i) = v(dg(i))
326 end do
327
328 !NEC$ IVDEP
329 do i = 1, m
330 u(gd(i)) = w(i)
331 end do
332
333 end subroutine gs_scatter_kernel
334
335end 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_kernel_mul(v, m, o, dg, u, n, gd, nb, b, w)
Gather kernel for multiplication of data .
Definition gs_sx.f90:179
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:253
subroutine gs_scatter_sx(this, v, m, dg, u, n, gd, nb, b, bo, shrd, event)
Scatter kernel.
Definition gs_sx.f90:289
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:216
subroutine gs_gather_sx(this, v, m, o, dg, u, n, gd, nb, b, bo, op, shrd)
Gather kernel.
Definition gs_sx.f90:94
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:141
subroutine gs_scatter_kernel(v, m, dg, u, n, gd, nb, b, w)
Scatter kernel .
Definition gs_sx.f90:312
subroutine gs_sx_init(this, nlocal, nshared, nlcl_blks, nshrd_blks)
SX backend initialisation.
Definition gs_sx.f90:59
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