Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
set_convect_rst_kernel.h
Go to the documentation of this file.
1#ifndef __MATH_SET_CONVECT_RST_KERNEL_H__
2#define __MATH_SET_CONVECT_RST_KERNEL_H__
3/*
4 Copyright (c) 2021-2023, 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
41template< typename T, const int LX, const int CHUNKS >
45 const T * __restrict__ cx,
46 const T * __restrict__ cy,
47 const T * __restrict__ cz,
48 const T * __restrict__ drdx,
49 const T * __restrict__ dsdx,
50 const T * __restrict__ dtdx,
51 const T * __restrict__ drdy,
52 const T * __restrict__ dsdy,
53 const T * __restrict__ dtdy,
54 const T * __restrict__ drdz,
55 const T * __restrict__ dsdz,
56 const T * __restrict__ dtdz,
57 const T * __restrict__ w3) {
58
59 int i,j,k;
60
61 const int e = blockIdx.x;
62 const int iii = threadIdx.x;
63 const int nchunks = (LX * LX * LX - 1) / CHUNKS + 1;
64
65 for (int n = 0; n < nchunks; n++) {
66 const int ijk = iii + n * CHUNKS;
67 const int jk = ijk / LX;
68 i = ijk - jk * LX;
69 k = jk / LX;
70 j = jk - k * LX;
71 if ( i < LX && j < LX && k < LX ) {
72
73 cr[ijk + e * LX * LX * LX] = w3[ijk]
74 * (drdx[ijk + e * LX * LX * LX] * cx[ijk + e * LX * LX * LX]
75 + drdy[ijk + e * LX * LX * LX] * cy[ijk + e * LX * LX * LX]
76 + drdz[ijk + e * LX * LX * LX] * cz[ijk + e * LX * LX * LX]);
77
78 cs[ijk + e * LX * LX * LX] = w3[ijk]
79 * (dsdx[ijk + e * LX * LX * LX] * cx[ijk + e * LX * LX * LX]
80 + dsdy[ijk + e * LX * LX * LX] * cy[ijk + e * LX * LX * LX]
81 + dsdz[ijk + e * LX * LX * LX] * cz[ijk + e * LX * LX * LX]);
82
83 ct[ijk + e * LX * LX * LX] = w3[ijk]
84 * (dtdx[ijk + e * LX * LX * LX] * cx[ijk + e * LX * LX * LX]
85 + dtdy[ijk + e * LX * LX * LX] * cy[ijk + e * LX * LX * LX]
86 + dtdz[ijk + e * LX * LX * LX] * cz[ijk + e * LX * LX * LX]);
87
88 }
89 }
90
91}
92
93template< typename T, const int LX >
111
112 const int e = blockIdx.x;
113 const int j = threadIdx.y;
114 const int i = threadIdx.x;
115 const int ij = i + j * LX;
116 const int ele = e*LX*LX*LX;
117
118 #pragma unroll
119 for (int k = 0; k < LX; ++k) {
120 const int ijk = ij + k*LX*LX;
121 const T W3 = w3[ijk];
122
123 cr[ijk + ele] = W3 * (drdx[ijk + ele] * cx[ijk + ele]
124 + drdy[ijk + ele] * cy[ijk + ele]
125 + drdz[ijk + ele] * cz[ijk + ele]);
126
127 cs[ijk + ele] = W3 * (dsdx[ijk + ele] * cx[ijk + ele]
128 + dsdy[ijk + ele] * cy[ijk + ele]
129 + dsdz[ijk + ele] * cz[ijk + ele]);
130
131 ct[ijk + ele] = W3 * (dtdx[ijk + ele] * cx[ijk + ele]
132 + dtdy[ijk + ele] * cy[ijk + ele]
133 + dtdz[ijk + ele] * cz[ijk + ele]);
134 }
135}
136
137
138#endif // __MATH_SET_CONVECT_RST_KERNEL_H__
const int e
__global__ void const T *__restrict__ const T *__restrict__ cr
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
__global__ void T *__restrict__ 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__ dtdz
__global__ void T *__restrict__ T *__restrict__ ct
__global__ void __launch_bounds__(LX *LX, 3) set_convect_rst_kernel_kstep(T *__restrict__ cr
__global__ void set_convect_rst_kernel_1d(T *__restrict__ cr, T *__restrict__ cs, T *__restrict__ ct, const T *__restrict__ cx, const T *__restrict__ cy, const T *__restrict__ cz, 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__ w3)
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ cz
__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__ dtdy
__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__ dsdy
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ cx
const int i
__global__ void T *__restrict__ cs
const int ij
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ cy
__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__ w3
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dsdx
const int ele
__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__ dtdx
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__ drdy
__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__ dsdz