Neko  0.9.0
A portable framework for high-order spectral element flow simulations
cfl_kernel.h
Go to the documentation of this file.
1 #ifndef __MATH_CFL_KERNEL_H__
2 #define __MATH_CFL_KERNEL_H__
3 /*
4  Copyright (c) 2022, The Neko Authors
5  All rights reserved.
6 
7  Redistribution and use in source and binary forms, with or without
8  modification, are permitted provided that the following conditions
9  are met:
10 
11  * Redistributions of source code must retain the above copyright
12  notice, this list of conditions and the following disclaimer.
13 
14  * Redistributions in binary form must reproduce the above
15  copyright notice, this list of conditions and the following
16  disclaimer in the documentation and/or other materials provided
17  with the distribution.
18 
19  * Neither the name of the authors nor the names of its
20  contributors may be used to endorse or promote products derived
21  from this software without specific prior written permission.
22 
23  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
24  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
25  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
26  FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
27  COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
28  INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
29  BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
30  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
31  CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
32  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
33  ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
34  POSSIBILITY OF SUCH DAMAGE.
35 */
36 
37 
41 template< typename T>
42 __inline__ __device__ T cfl_reduce_warp(T val) {
43  val = fmax(val, __shfl_down_sync(0xffffffff, val, 16));
44  val = fmax(val, __shfl_down_sync(0xffffffff, val, 8));
45  val = fmax(val, __shfl_down_sync(0xffffffff, val, 4));
46  val = fmax(val, __shfl_down_sync(0xffffffff, val, 2));
47  val = fmax(val, __shfl_down_sync(0xffffffff, val, 1));
48  return val;
49 }
50 
54 template< typename T >
55 __global__ void cfl_reduce_kernel(T * bufred, const int n) {
56 
57  T cfl = 0;
58  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
59  const int str = blockDim.x * gridDim.x;
60  for (int i = idx; i<n ; i += str)
61  {
62  cfl = fmax(cfl, bufred[i]);
63  }
64 
65  __shared__ T shared[32];
66  unsigned int lane = threadIdx.x % warpSize;
67  unsigned int wid = threadIdx.x / warpSize;
68 
69  cfl = cfl_reduce_warp<T>(cfl);
70  if (lane == 0)
71  shared[wid] = cfl;
72  __syncthreads();
73 
74  cfl = (threadIdx.x < blockDim.x / warpSize) ? shared[lane] : 0;
75  if (wid == 0)
76  cfl = cfl_reduce_warp<T>(cfl);
77 
78  if (threadIdx.x == 0)
79  bufred[blockIdx.x] = cfl;
80 }
81 
86 template< typename T, const int LX, const int CHUNKS >
87 __global__ void cfl_kernel(const T dt,
88  const T * __restrict__ u,
89  const T * __restrict__ v,
90  const T * __restrict__ w,
91  const T * __restrict__ drdx,
92  const T * __restrict__ dsdx,
93  const T * __restrict__ dtdx,
94  const T * __restrict__ drdy,
95  const T * __restrict__ dsdy,
96  const T * __restrict__ dtdy,
97  const T * __restrict__ drdz,
98  const T * __restrict__ dsdz,
99  const T * __restrict__ dtdz,
100  const T * __restrict__ dr_inv,
101  const T * __restrict__ ds_inv,
102  const T * __restrict__ dt_inv,
103  const T * __restrict__ jacinv,
104  T * __restrict__ cfl_h) {
105 
106  int i,j,k;
107 
108  const int e = blockIdx.x;
109  const int iii = threadIdx.x;
110  const int nchunks = (LX * LX * LX - 1) / CHUNKS + 1;
111  const unsigned int lane = threadIdx.x % warpSize;
112  const unsigned int wid = threadIdx.x / warpSize;
113 
114  __shared__ T shu[LX * LX * LX];
115  __shared__ T shv[LX * LX * LX];
116  __shared__ T shw[LX * LX * LX];
117 
118  __shared__ T shdr_inv[LX];
119  __shared__ T shds_inv[LX];
120  __shared__ T shdt_inv[LX];
121 
122  __shared__ T shjacinv[LX * LX * LX];
123 
124  __shared__ T shared[32];
125 
126  if (iii < LX) {
127  shdr_inv[iii] = dr_inv[iii];
128  shds_inv[iii] = ds_inv[iii];
129  shdt_inv[iii] = dt_inv[iii];
130  }
131 
132  j = iii;
133  while(j < (LX * LX * LX)) {
134  shu[j] = u[j + e * LX * LX * LX];
135  shv[j] = v[j + e * LX * LX * LX];
136  shw[j] = w[j + e * LX * LX * LX];
137 
138  shjacinv[j] = jacinv[j + e * LX * LX * LX];
139 
140  j = j + CHUNKS;
141  }
142 
143  __syncthreads();
144 
145  T cfl_tmp = 0.0;
146  for (int n = 0; n < nchunks; n++) {
147  const int ijk = iii + n * CHUNKS;
148  const int jk = ijk / LX;
149  i = ijk - jk * LX;
150  k = jk / LX;
151  j = jk - k * LX;
152  if ( i < LX && j < LX && k < LX) {
153  const T cflr = fabs( dt * ( ( shu[ijk] * drdx[ijk + e * LX * LX * LX]
154  + shv[ijk] * drdy[ijk + e * LX * LX * LX]
155  * shw[ijk] * drdz[ijk + e * LX * LX * LX]
156  ) * shjacinv[ijk]) * shdr_inv[i]);
157  const T cfls = fabs( dt * ( ( shu[ijk] * dsdx[ijk + e * LX * LX * LX]
158  + shv[ijk] * dsdy[ijk + e * LX * LX * LX]
159  + shw[ijk] * dsdz[ijk + e * LX * LX * LX]
160  ) * shjacinv[ijk]) * shds_inv[j]);
161  const T cflt = fabs( dt * ( ( shu[ijk] * dtdx[ijk + e * LX * LX * LX]
162  + shv[ijk] * dtdy[ijk + e * LX * LX * LX]
163  + shw[ijk] * dtdz[ijk + e * LX * LX * LX]
164  ) * shjacinv[ijk]) * shdt_inv[k]);
165 
166  cfl_tmp = fmax(cflr + cfls + cflt, cfl_tmp);
167 
168  }
169  }
170 
171  cfl_tmp = cfl_reduce_warp<T>(cfl_tmp);
172  if (lane == 0)
173  shared[wid] = cfl_tmp;
174  __syncthreads();
175 
176  cfl_tmp = (threadIdx.x < blockDim.x / warpSize) ? shared[lane] : 0;
177  if (wid == 0)
178  cfl_tmp = cfl_reduce_warp<T>(cfl_tmp);
179 
180  if (threadIdx.x == 0)
181  cfl_h[blockIdx.x] = cfl_tmp;
182 }
183 
184 
185 #endif // __MATH_CFL_KERNEL_H__
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ drdy
__shared__ T shu[LX *LX]
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ drdz
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dsdz
__shared__ T shw[LX *LX]
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dsdy
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ w
const int i
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dtdy
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ u
const int e
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ drdx
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dtdz
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dsdx
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ v
const int j
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dtdx
__syncthreads()
__shared__ T shv[LX *LX]
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ jacinv
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dt
Definition: cdtp_kernel.h:109
__global__ void cfl_kernel(const T dt, const T *__restrict__ u, const T *__restrict__ v, const T *__restrict__ w, const T *__restrict__ drdx, const T *__restrict__ dsdx, const T *__restrict__ dtdx, const T *__restrict__ drdy, const T *__restrict__ dsdy, const T *__restrict__ dtdy, const T *__restrict__ drdz, const T *__restrict__ dsdz, const T *__restrict__ dtdz, const T *__restrict__ dr_inv, const T *__restrict__ ds_inv, const T *__restrict__ dt_inv, const T *__restrict__ jacinv, T *__restrict__ cfl_h)
Definition: cfl_kernel.h:87
__global__ void cfl_reduce_kernel(T *bufred, const int n)
Definition: cfl_kernel.h:55
__inline__ __device__ T cfl_reduce_warp(T val)
Definition: cfl_kernel.h:42
real * bufred
Definition: math.cu:479
real(kind=rp) function, public cfl(dt, u, v, w, Xh, coef, nelv, gdim)
Definition: operators.f90:394