Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
gs_kernels.h
Go to the documentation of this file.
1/*
2 Copyright (c) 2021-2022, The Neko Authors
3 All rights reserved.
4
5 Redistribution and use in source and binary forms, with or without
6 modification, are permitted provided that the following conditions
7 are met:
8
9 * Redistributions of source code must retain the above copyright
10 notice, this list of conditions and the following disclaimer.
11
12 * Redistributions in binary form must reproduce the above
13 copyright notice, this list of conditions and the following
14 disclaimer in the documentation and/or other materials provided
15 with the distribution.
16
17 * Neither the name of the authors nor the names of its
18 contributors may be used to endorse or promote products derived
19 from this software without specific prior written permission.
20
21 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22 "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23 LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
24 FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
25 COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
26 INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
27 BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
29 CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
30 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
31 ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
32 POSSIBILITY OF SUCH DAMAGE.
33*/
34
35#ifndef __GS_GS_KERNELS__
36#define __GS_GS_KERNELS__
37
42template< typename T >
44 const int m,
45 const int o,
46 const int * __restrict__ dg,
47 const T * __restrict__ u,
48 const int n,
49 const int * __restrict__ gd,
50 const int nb,
51 const int * __restrict__ b,
52 const int * __restrict__ bo) {
53
54 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
55 const int str = blockDim.x * gridDim.x;
56
57 for (int i = idx; i < nb; i += str) {
58 const int blk_len = b[i];
59 const int k = bo[i];
60 T tmp = u[gd[k] - 1];
61 for (int j = 1; j < blk_len; j++) {
62 tmp += u[gd[k + j] - 1];
63 }
64 v[dg[k] - 1] = tmp;
65 }
66
67 if (o < 0) {
68 for (int i = ((abs(o) - 1) + idx); i < m ; i += str) {
69 v[dg[i] - 1] = u[gd[i] - 1];
70 }
71 }
72 else {
73 if ((idx%2 == 0)) {
74 for (int i = ((o - 1) + idx); i < m ; i += str) {
75 T tmp = u[gd[i] - 1] + u[gd[i+1] - 1];
76 v[dg[i] - 1] = tmp;
77 }
78 }
79 }
80
81}
82
87template< typename T >
89 const int m,
90 const int o,
91 const int * __restrict__ dg,
92 const T * __restrict__ u,
93 const int n,
94 const int * __restrict__ gd,
95 const int nb,
96 const int * __restrict__ b,
97 const int * __restrict__ bo) {
98
99 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
100 const int str = blockDim.x * gridDim.x;
101
102 for (int i = idx; i < nb; i += str) {
103 const int blk_len = b[i];
104 const int k = bo[i];
105 T tmp = u[gd[k] - 1];
106 for (int j = 1; j < blk_len; j++) {
107 tmp *= u[gd[k + j] - 1];
108 }
109 v[dg[k] - 1] = tmp;
110 }
111
112 if (o < 0) {
113 for (int i = ((abs(o) - 1) + idx); i < m ; i += str) {
114 v[dg[i] - 1] = u[gd[i] - 1];
115 }
116 }
117 else {
118 if ((idx%2 == 0)) {
119 for (int i = ((o - 1) + idx); i < m ; i += str) {
120 T tmp = u[gd[i] - 1] * u[gd[i+1] - 1];
121 v[dg[i] - 1] = tmp;
122 }
123 }
124 }
125
126}
127
132template< typename T >
134 const int m,
135 const int o,
136 const int * __restrict__ dg,
137 const T * __restrict__ u,
138 const int n,
139 const int * __restrict__ gd,
140 const int nb,
141 const int * __restrict__ b,
142 const int * __restrict__ bo) {
143
144 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
145 const int str = blockDim.x * gridDim.x;
146
147 for (int i = idx; i < nb; i += str) {
148 const int blk_len = b[i];
149 const int k = bo[i];
150 T tmp = u[gd[k] - 1];
151 for (int j = 1; j < blk_len; j++) {
152 tmp = min(u[gd[k + j] - 1], tmp);
153 }
154 v[dg[k] - 1] = tmp;
155 }
156
157 if (o < 0) {
158 for (int i = ((abs(o) - 1) + idx); i < m ; i += str) {
159 v[dg[i] - 1] = u[gd[i] - 1];
160 }
161 }
162 else {
163 if ((idx%2 == 0)) {
164 for (int i = ((o - 1) + idx); i < m ; i += str) {
165 T tmp = min(u[gd[i] - 1], u[gd[i+1] - 1]);
166 v[dg[i] - 1] = tmp;
167 }
168 }
169 }
170
171}
172
177template< typename T >
179 const int m,
180 const int o,
181 const int * __restrict__ dg,
182 const T * __restrict__ u,
183 const int n,
184 const int * __restrict__ gd,
185 const int nb,
186 const int * __restrict__ b,
187 const int * __restrict__ bo) {
188
189 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
190 const int str = blockDim.x * gridDim.x;
191
192 for (int i = idx; i < nb; i += str) {
193 const int blk_len = b[i];
194 const int k = bo[i];
195 T tmp = u[gd[k] - 1];
196 for (int j = 1; j < blk_len; j++) {
197 tmp = max(u[gd[k + j] - 1], tmp);
198 }
199 v[dg[k] - 1] = tmp;
200 }
201
202 if (o < 0) {
203 for (int i = ((abs(o) - 1) + idx); i < m ; i += str) {
204 v[dg[i] - 1] = u[gd[i] - 1];
205 }
206 }
207 else {
208 if ((idx%2 == 0)) {
209 for (int i = ((o - 1) + idx); i < m ; i += str) {
210 T tmp = max(u[gd[i] - 1], u[gd[i+1] - 1]);
211 v[dg[i] - 1] = tmp;
212 }
213 }
214 }
215
216}
217
222template< typename T >
224 const int m,
225 const int * __restrict__ dg,
226 T * __restrict__ u,
227 const int n,
228 const int * __restrict__ gd,
229 const int nb,
230 const int *__restrict__ b,
231 const int *__restrict__ bo) {
232
233 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
234 const int str = blockDim.x * gridDim.x;
235
236 for (int i = idx; i < nb; i += str) {
237 const int blk_len = b[i];
238 const int k = bo[i];
239 T tmp = v[dg[k] - 1];
240 for (int j = 0; j < blk_len; j++) {
241 u[gd[k + j] - 1] = tmp;
242 }
243 }
244
245 int facet_offset = 0;
246 if ( nb > 0) {
247 facet_offset = bo[nb - 1] + b[nb - 1];
248 }
249
250
251 for (int i = (facet_offset + idx); i < m; i += str) {
252 u[gd[i] - 1] = v[dg[i] - 1];
253 }
254
255}
256
257template< typename T >
260 const int32_t * __restrict__ dof,
261 const int n) {
262
263 const int j = threadIdx.x + blockDim.x * blockIdx.x;
264
265 if (j >= n)
266 return;
267
268 buf[j] = u[dof[j]-1];
269}
270
271
272template< typename T >
274 const T * __restrict__ buf,
275 const int32_t * __restrict__ dof,
276 const int n) {
277
278 const int j = threadIdx.x + blockDim.x * blockIdx.x;
279
280 if (j >= n)
281 return;
282
283 const int32_t idx = dof[j];
284 const T val = buf[j];
285 if (idx < 0) {
286 atomicAdd(&u[-idx-1], val);
287 } else {
288 u[idx-1] += val;
289 }
290}
291
292template<typename T>
294
295template<>
297 float old;
300
301 return old;
302}
303
304template<>
306 double old;
307 old = !signbit(val) ? __longlong_as_double(atomicMin((unsigned long long*)address,
309 __longlong_as_double(atomicMax((unsigned long long*)address,
311 return old;
312}
313
314template< typename T >
316 const T * __restrict__ buf,
317 const int32_t * __restrict__ dof,
318 const int n) {
319
320 const int j = threadIdx.x + blockDim.x * blockIdx.x;
321
322 if (j >= n)
323 return;
324
325 const int32_t idx = dof[j];
326 const T val = buf[j];
327
328 if (idx < 0) {
329 // Use atomicMin for shared nodal points
330 atomicMinFloat(&u[-idx-1], val);
331 } else {
332 // Directly compute min for nodal points on edges
333 u[idx-1] = min(u[idx-1], val);
334 }
335}
336
337template<typename T>
339
340template<>
342 float old;
345 return old;
346}
347
348template<>
350 double old;
351 old = !signbit(val) ? __longlong_as_double(atomicMax((unsigned long long*)address,
353 __longlong_as_double(atomicMin((unsigned long long*)address,
355 return old;
356}
357
358template< typename T >
360 const T * __restrict__ buf,
361 const int32_t * __restrict__ dof,
362 const int n) {
363
364 const int j = threadIdx.x + blockDim.x * blockIdx.x;
365
366 if (j >= n)
367 return;
368
369 const int32_t idx = dof[j];
370 const T val = buf[j];
371
372 if (idx < 0) {
373 // Use atomicMax for shared nodal points
374 atomicMaxFloat(&u[-idx-1], val);
375 } else {
376 // Directly compute min for nodal points on edges
377 u[idx-1] = max(u[idx-1], val);
378 }
379}
380
381#endif // __GS_GS_KERNELS__
const int i
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ u
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ v
const int j
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
__device__ T atomicMaxFloat(T *address, T val)
__global__ void gs_unpack_max_kernel(T *__restrict__ u, const T *__restrict__ buf, const int32_t *__restrict__ dof, const int n)
Definition gs_kernels.h:363
__global__ void gather_kernel_min(T *__restrict__ v, const int m, const int o, const int *__restrict__ dg, const T *__restrict__ u, const int n, const int *__restrict__ gd, const int nb, const int *__restrict__ b, const int *__restrict__ bo)
Definition gs_kernels.h:133
__global__ void scatter_kernel(T *__restrict__ v, const int m, const int *__restrict__ dg, T *__restrict__ u, const int n, const int *__restrict__ gd, const int nb, const int *__restrict__ b, const int *__restrict__ bo)
Definition gs_kernels.h:223
__global__ void gather_kernel_mul(T *__restrict__ v, const int m, const int o, const int *__restrict__ dg, const T *__restrict__ u, const int n, const int *__restrict__ gd, const int nb, const int *__restrict__ b, const int *__restrict__ bo)
Definition gs_kernels.h:88
__device__ T atomicMinFloat(T *address, T val)
__global__ void gs_pack_kernel(const T *__restrict__ u, T *__restrict__ buf, const int32_t *__restrict__ dof, const int n)
Definition gs_kernels.h:256
__global__ void gs_unpack_min_kernel(T *__restrict__ u, const T *__restrict__ buf, const int32_t *__restrict__ dof, const int n)
Definition gs_kernels.h:317
__global__ void gather_kernel_add(T *__restrict__ v, const int m, const int o, const int *__restrict__ dg, const T *__restrict__ u, const int n, const int *__restrict__ gd, const int nb, const int *__restrict__ b, const int *__restrict__ bo)
Definition gs_kernels.h:43
__device__ double atomicMaxFloat< double >(double *address, double val)
Definition gs_kernels.h:351
__device__ float atomicMaxFloat< float >(float *address, float val)
Definition gs_kernels.h:343
__global__ void gs_unpack_add_kernel(T *__restrict__ u, const T *__restrict__ buf, const int32_t *__restrict__ dof, const int n)
Definition gs_kernels.h:271
__global__ void gather_kernel_max(T *__restrict__ v, const int m, const int o, const int *__restrict__ dg, const T *__restrict__ u, const int n, const int *__restrict__ gd, const int nb, const int *__restrict__ b, const int *__restrict__ bo)
Definition gs_kernels.h:178
__device__ float atomicMinFloat< float >(float *address, float val)
Definition gs_kernels.h:296
__device__ double atomicMinFloat< double >(double *address, double val)
Definition gs_kernels.h:305
real * buf
Definition pipecg_aux.cu:42
#define max(a, b)
Definition tensor.cu:40