Neko  0.9.0
A portable framework for high-order spectral element flow simulations
pc_jacobi.cu
Go to the documentation of this file.
1 /*
2  Copyright (c) 2021, 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 #include <stdio.h>
36 #include <device/device_config.h>
37 
38 template< typename T, const int LX >
39 __global__ void jacobi_kernel(T * __restrict__ du,
40  const T * __restrict__ dxt,
41  const T * __restrict__ dyt,
42  const T * __restrict__ dzt,
43  const T * __restrict__ G11,
44  const T * __restrict__ G22,
45  const T * __restrict__ G33,
46  const T * __restrict__ G12,
47  const T * __restrict__ G13,
48  const T * __restrict__ G23,
49  const int nel) {
50  const int idx = threadIdx.x + blockIdx.x * blockDim.x;
51  const int e = idx / (LX*LX*LX);
52  const int ijk = idx - e*LX*LX*LX;
53  const int jk = ijk / LX;
54  const int i = ijk - jk * LX;
55  const int k = jk / LX;
56  const int j = jk - k * LX;
57 
58  if (e >= nel)
59  return;
60 
61  T d = 0.0;
62 
63  for (int l = 0; l < LX; l++) {
64  T g = G11[l + LX*j + LX*LX*k + LX*LX*LX*e];
65  T t = dxt[i + LX*l];
66  d += g*t*t;
67  }
68 
69  for (int l = 0; l < LX; l++) {
70  T g = G22[i + LX*l + LX*LX*k + LX*LX*LX*e];
71  T t = dyt[j + LX*l];
72  d += g*t*t;
73  }
74 
75  for (int l = 0; l < LX; l++) {
76  T g = G33[i + LX*j + LX*LX*l + LX*LX*LX*e];
77  T t = dzt[k + LX*l];
78  d += g*t*t;
79  }
80 
81  // Corrections for deformed elements
82  if (i == 0 || i == LX-1) {
83  d += G12[i + LX*j + LX*LX*k + LX*LX*LX*e] * dxt[i + LX*i] * dyt[j + LX*j];
84  d += G13[i + LX*j + LX*LX*k + LX*LX*LX*e] * dxt[i + LX*i] * dzt[k + LX*k];
85  }
86 
87  if (j == 0 || j == LX-1) {
88  d += G12[i + LX*j + LX*LX*k + LX*LX*LX*e] * dyt[j + LX*j] * dxt[i + LX*i];
89  d += G23[i + LX*j + LX*LX*k + LX*LX*LX*e] * dyt[j + LX*j] * dzt[k + LX*k];
90  }
91 
92  if (k == 0 || k == LX-1) {
93  d += G13[i + LX*j + LX*LX*k + LX*LX*LX*e] * dzt[k + LX*k] * dxt[i + LX*i];
94  d += G23[i + LX*j + LX*LX*k + LX*LX*LX*e] * dzt[k + LX*k] * dyt[j + LX*j];
95  }
96 
97  du[idx] = d;
98 }
99 
100 extern "C" {
101  void cuda_jacobi_update(void *d,
102  void *dxt, void *dyt, void *dzt,
103  void *G11, void *G22, void *G33,
104  void *G12, void *G13, void *G23,
105  int *nel, int *lxp) {
106 
107  const int lx = *lxp;
108  const int threads = 1024;
109  const int blocks = ((*nel * lx*lx*lx) + threads - 1) / threads;
110  const cudaStream_t stream = (cudaStream_t) glb_cmd_queue;
111 
112 #define CASE(N)\
113  case N:\
114  jacobi_kernel<real, N><<<blocks, threads, 0, stream>>>( \
115  (real*)d,\
116  (real*)dxt, (real*)dyt, (real*)dzt,\
117  (real*)G11, (real*)G22, (real*)G33,\
118  (real*)G12, (real*)G13, (real*)G23,\
119  *nel);\
120  break
121 
122  switch (lx) {
123  CASE(1);
124  CASE(2);
125  CASE(3);
126  CASE(4);
127  CASE(5);
128  CASE(6);
129  CASE(7);
130  CASE(8);
131  CASE(9);
132  CASE(10);
133  CASE(11);
134  CASE(12);
135  CASE(13);
136  CASE(14);
137  CASE(15);
138  default:
139  fprintf(stderr, __FILE__ ": size not supported: %d\n", lx);
140  }
141 
142  cudaError_t err = cudaGetLastError();
143  if (err != cudaSuccess) {
144  fprintf(stderr, __FILE__ ": %s\n", cudaGetErrorString(err));
145  }
146 
147  }
148 }
const int i
const int e
const int j
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dzt
Definition: cdtp_kernel.h:112
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dyt
Definition: cdtp_kernel.h:111
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dxt
Definition: cdtp_kernel.h:110
void * glb_cmd_queue
#define CASE(N)
__global__ void jacobi_kernel(T *__restrict__ du, const T *__restrict__ dxt, const T *__restrict__ dyt, const T *__restrict__ dzt, const T *__restrict__ G11, const T *__restrict__ G22, const T *__restrict__ G33, const T *__restrict__ G12, const T *__restrict__ G13, const T *__restrict__ G23, const int nel)
Definition: pc_jacobi.cu:39
void cuda_jacobi_update(void *d, void *dxt, void *dyt, void *dzt, void *G11, void *G22, void *G33, void *G12, void *G13, void *G23, int *nel, int *lxp)
Definition: pc_jacobi.cu:101