Neko  0.9.0
A portable framework for high-order spectral element flow simulations
fdm_kernel.h
Go to the documentation of this file.
1 #ifndef __MATH_FDM_KERNEL_H__
2 #define __MATH_FDM_KERNEL_H__
3 /*
4  Copyright (c) 2021-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 template< typename T, const int NL >
38 __global__ void fdm_do_fast_kernel(T * __restrict__ e,
39  T * __restrict__ r,
40  T * __restrict__ s,
41  T * __restrict__ d) {
42  __shared__ T shwork[NL*NL*NL];
43  __shared__ T shwork2[NL*NL*NL];
44  __shared__ T A[NL*NL];
45  __shared__ T Bt[NL*NL];
46  __shared__ T Ct[NL*NL];
47 
48  const int idx = threadIdx.x;
49  const int str = blockDim.x;
50  const int el = blockIdx.x;
51  if( idx < NL*NL){
52  A[idx] = s[idx+NL*NL+el*NL*NL*3*2];
53  Bt[idx] = s[idx+2*NL*NL+el*NL*NL*3*2];
54  Ct[idx] = s[idx+2*2*NL*NL+el*NL*NL*3*2];
55  }
56  __syncthreads();
57 
58  for (int ii = idx; ii< NL*NL*NL; ii += str){
59  T tmp = 0.0;
60  int j = ii/NL;
61  int i = ii - j*NL;
62  for( int l = 0; l < NL; l++){
63  tmp += A[i+l*NL]*r[l+NL*j+el*NL*NL*NL];
64  }
65  shwork[ii] = tmp;
66  }
67  __syncthreads();
68  for (int ijk = idx; ijk< NL*NL*NL; ijk += str){
69  const int jk = ijk / NL;
70  const int i = ijk - jk * NL;
71  const int k = jk / NL;
72  const int j = jk - k * NL;
73  T tmp = 0.0;
74  const int ik2 = i + k*NL*NL;
75  for( int l = 0; l < NL; l++){
76  tmp += Bt[l+j*NL]*shwork[l*NL+ik2];
77  }
78  shwork2[ijk] = tmp;
79  }
80  __syncthreads();
81  for (int ijk = idx; ijk< NL*NL*NL; ijk += str){
82  const int jk = ijk / NL;
83  const int i = ijk - jk * NL;
84  const int k = jk / NL;
85  const int j = jk - k * NL;
86  T tmp = 0.0;
87  const int ij2 = i + j*NL;
88  for( int l = 0; l < NL; l++){
89  tmp += Ct[l+k*NL]*shwork2[ij2 + l*NL*NL];
90  }
91  r[ijk+el*NL*NL*NL] = tmp*d[ijk+el*NL*NL*NL];
92  }
93  __syncthreads();
94  if( idx < NL*NL){
95  A[idx] = s[idx+el*NL*NL*3*2];
96  Bt[idx] = s[idx+NL*NL+2*NL*NL+el*NL*NL*3*2];
97  Ct[idx] = s[idx+NL*NL+2*2*NL*NL+el*NL*NL*3*2];
98  }
99  __syncthreads();
100 
101  for (int ii = idx; ii< NL*NL*NL; ii += str){
102  T tmp = 0.0;
103  int j = ii/NL;
104  int i = ii - j*NL;
105  for( int l = 0; l < NL; l++){
106  tmp += A[i+l*NL]*r[l+NL*j+el*NL*NL*NL];
107  }
108  shwork[ii] = tmp;
109  }
110  __syncthreads();
111  for (int ijk = idx; ijk< NL*NL*NL; ijk += str){
112  const int jk = ijk / NL;
113  const int i = ijk - jk * NL;
114  const int k = jk / NL;
115  const int j = jk - k * NL;
116  T tmp = 0.0;
117  const int ik2 = i + k*NL*NL;
118  for( int l = 0; l < NL; l++){
119  tmp += Bt[l+j*NL]*shwork[l*NL+ik2];
120  }
121  shwork2[ijk] = tmp;
122  }
123  __syncthreads();
124  for (int ijk = idx; ijk< NL*NL*NL; ijk += str){
125  const int jk = ijk / NL;
126  const int i = ijk - jk * NL;
127  const int k = jk / NL;
128  const int j = jk - k * NL;
129  T tmp = 0.0;
130  const int ij2 = i + j*NL;
131  for( int l = 0; l < NL; l++){
132  tmp += Ct[l+k*NL]*shwork2[ij2 + l*NL*NL];
133  }
134  e[ijk+el*NL*NL*NL] = tmp;
135  }
136 
137 }
138 
139 
140 
141 #endif // __MATH_FDM_KERNEL_H__
const int i
const int e
const int j
__syncthreads()
__global__ void fdm_do_fast_kernel(T *__restrict__ e, T *__restrict__ r, T *__restrict__ s, T *__restrict__ d)
Definition: fdm_kernel.h:38