Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
schwarz_kernel.h
Go to the documentation of this file.
1#ifndef __MATH_SCHWARZ_KERNEL_H__
2#define __MATH_SCHWARZ_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
45template< typename T, const int NX>
47 const int l1,
48 const T f1,
49 T * a2,
50 const int l2,
51 const T f2){
52
53 const int idx = threadIdx.x;
54 const int el = blockIdx.x*NX*NX*NX;
55 const int x = idx%(NX-2) + 1;
56 const int y = idx/(NX-2) + 1;
57 int idx1,idx2;
58
59 idx1 = l1 + x*NX + y*NX*NX + el;
60 idx2 = l2 + x*NX + y*NX*NX + el;
61 a1[idx1] = f1*a1[idx1] + f2*a2[idx2];
62
63 idx1 = (NX-1-l1) + x*NX + y*NX*NX + el;
64 idx2 = (NX-1-l2) + x*NX + y*NX*NX + el;
65 a1[idx1] = f1*a1[idx1] + f2*a2[idx2];
66
68
69 idx1 = x + l1*NX + y*NX*NX + el;
70 idx2 = x + l2*NX + y*NX*NX + el;
71 a1[idx1] = f1*a1[idx1] + f2*a2[idx2];
72
73 idx1 = x + (NX-1-l1)*NX + y*NX*NX + el;
74 idx2 = x + (NX-1-l2)*NX + y*NX*NX + el;
75 a1[idx1] = f1*a1[idx1] + f2*a2[idx2];
76
78
79 idx1 = x + y*NX + l1*NX*NX + el;
80 idx2 = x + y*NX + l2*NX*NX + el;
81 a1[idx1] = f1*a1[idx1] + f2*a2[idx2];
82
83 idx1 = x + y*NX + (NX-1-l1)*NX*NX + el;
84 idx2 = x + y*NX + (NX-1-l2)*NX*NX + el;
85 a1[idx1] = f1*a1[idx1] + f2*a2[idx2];
86}
87
88
92template< typename T>
94 T *__restrict__ b,
95 const int nx) {
96
97 const int idx = threadIdx.x;
98 const int nx2 = nx+2;
99 const int el2 = blockIdx.x*nx2*nx2*nx2;
100 const int el = blockIdx.x*nx*nx*nx;
101 for(int i = idx; i<nx2*nx2*nx2; i+=blockDim.x){
102 a[i+el2] = 0.0;
103 }
105 for(int ijk = idx; ijk<nx*nx*nx; ijk+=blockDim.x){
106 const int jk = ijk / nx;
107 const int i = ijk - jk * nx;
108 const int k = jk / nx;
109 const int j = jk - k * nx;
110 a[(i+1)+(j+1)*nx2+(k+1)*nx2*nx2+el2] = b[ijk+el];
111 }
112}
113
114template< typename T>
117 const int nx) {
118
119 const int idx = threadIdx.x;
120 const int nx2 = nx+2;
121 const int el2 = blockIdx.x*nx2*nx2*nx2;
122 const int el = blockIdx.x*nx*nx*nx;
123 for(int ijk = idx; ijk<nx*nx*nx; ijk+=blockDim.x){
124 const int jk = ijk / nx;
125 const int i = ijk - jk * nx;
126 const int k = jk / nx;
127 const int j = jk - k * nx;
128 b[ijk+el] = a[(i+1)+(j+1)*nx2+(k+1)*nx2*nx2+el2];
129 }
130}
131
132#endif // __MATH_SCHWARZ_KERNEL_H__
const int i
const int j
__syncthreads()
__global__ void const T *__restrict__ x
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
__global__ void schwarz_extrude_kernel(T *a1, const int l1, const T f1, T *a2, const int l2, const T f2)
__global__ void schwarz_toreg3d_kernel(T *__restrict__ b, T *__restrict__ a, const int nx)
__global__ void schwarz_toext3d_kernel(T *__restrict__ a, T *__restrict__ b, const int nx)