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