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