Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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
41template< typename T>
43 val = fmax(val, __shfl_down(val, 32));
44 val = fmax(val, __shfl_down(val, 16));
45 val = fmax(val, __shfl_down(val, 8));
46 val = fmax(val, __shfl_down(val, 4));
47 val = fmax(val, __shfl_down(val, 2));
48 val = fmax(val, __shfl_down(val, 1));
49 return val;
50}
51
55template< typename T >
56__global__ void cfl_reduce_kernel(T * bufred, const int n) {
57
58 T cfl = 0;
59 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
60 const int str = blockDim.x * gridDim.x;
61 for (int i = idx; i<n ; i += str)
62 {
63 cfl = fmax(cfl, bufred[i]);
64 }
65
67 unsigned int lane = threadIdx.x % warpSize;
68 unsigned int wid = threadIdx.x / warpSize;
69
70 cfl = cfl_reduce_warp<T>(cfl);
71 if (lane == 0)
72 shared[wid] = cfl;
74
75 cfl = (threadIdx.x < blockDim.x / warpSize) ? shared[lane] : 0;
76 if (wid == 0)
77 cfl = cfl_reduce_warp<T>(cfl);
78
79 if (threadIdx.x == 0)
80 bufred[blockIdx.x] = cfl;
81}
82
87template< typename T, const int LX, const int CHUNKS >
89 const T * __restrict__ u,
90 const T * __restrict__ v,
91 const T * __restrict__ w,
92 const T * __restrict__ drdx,
93 const T * __restrict__ dsdx,
94 const T * __restrict__ dtdx,
95 const T * __restrict__ drdy,
96 const T * __restrict__ dsdy,
97 const T * __restrict__ dtdy,
98 const T * __restrict__ drdz,
99 const T * __restrict__ dsdz,
100 const T * __restrict__ dtdz,
101 const T * __restrict__ dr_inv,
102 const T * __restrict__ ds_inv,
103 const T * __restrict__ dt_inv,
104 const T * __restrict__ jacinv,
105 T * __restrict__ cfl_h) {
106
107 int i,j,k;
108
109 const int e = blockIdx.x;
110 const int iii = threadIdx.x;
111 const int nchunks = (LX * LX * LX - 1) / CHUNKS + 1;
112 const unsigned int lane = threadIdx.x % warpSize;
113 const unsigned int wid = threadIdx.x / warpSize;
114
115 __shared__ T shu[LX * LX * LX];
116 __shared__ T shv[LX * LX * LX];
117 __shared__ T shw[LX * LX * LX];
118
122
124
125 __shared__ T shared[64];
126
127 if (iii < LX) {
128 shdr_inv[iii] = dr_inv[iii];
129 shds_inv[iii] = ds_inv[iii];
130 shdt_inv[iii] = dt_inv[iii];
131 }
132
133 j = iii;
134 while(j < (LX * LX * LX)) {
135 shu[j] = u[j + e * LX * LX * LX];
136 shv[j] = v[j + e * LX * LX * LX];
137 shw[j] = w[j + e * LX * LX * LX];
138
139 shjacinv[j] = jacinv[j + e * LX * LX * LX];
140
141 j = j + CHUNKS;
142 }
143
145
146 T cfl_tmp = 0.0;
147 for (int n = 0; n < nchunks; n++) {
148 const int ijk = iii + n * CHUNKS;
149 const int jk = ijk / LX;
150 i = ijk - jk * LX;
151 k = jk / LX;
152 j = jk - k * LX;
153 if ( i < LX && j < LX && k < LX) {
154 const T cflr = fabs( dt * ( ( shu[ijk] * drdx[ijk + e * LX * LX * LX]
155 + shv[ijk] * drdy[ijk + e * LX * LX * LX]
156 * shw[ijk] * drdz[ijk + e * LX * LX * LX]
157 ) * shjacinv[ijk]) * shdr_inv[i]);
158 const T cfls = fabs( dt * ( ( shu[ijk] * dsdx[ijk + e * LX * LX * LX]
159 + shv[ijk] * dsdy[ijk + e * LX * LX * LX]
160 + shw[ijk] * dsdz[ijk + e * LX * LX * LX]
161 ) * shjacinv[ijk]) * shds_inv[j]);
162 const T cflt = fabs( dt * ( ( shu[ijk] * dtdx[ijk + e * LX * LX * LX]
163 + shv[ijk] * dtdy[ijk + e * LX * LX * LX]
164 + shw[ijk] * dtdz[ijk + e * LX * LX * LX]
165 ) * shjacinv[ijk]) * shdt_inv[k]);
166
168
169 }
170 }
171
173 if (lane == 0)
174 shared[wid] = cfl_tmp;
176
177 cfl_tmp = (threadIdx.x < blockDim.x / warpSize) ? shared[lane] : 0;
178 if (wid == 0)
180
181 if (threadIdx.x == 0)
182 cfl_h[blockIdx.x] = cfl_tmp;
183}
184
185
186#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
__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
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
real * bufred
Definition math.cu:479