Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
facet_normal_kernel.h
Go to the documentation of this file.
1/*
2 Copyright (c) 2021-2022, The Neko Authors
3 All rights reserved.
4
5 Redistribution and use in source and binary forms, with or without
6 modification, are permitted provided that the following conditions
7 are met:
8
9 * Redistributions of source code must retain the above copyright
10 notice, this list of conditions and the following disclaimer.
11
12 * Redistributions in binary form must reproduce the above
13 copyright notice, this list of conditions and the following
14 disclaimer in the documentation and/or other materials provided
15 with the distribution.
16
17 * Neither the name of the authors nor the names of its
18 contributors may be used to endorse or promote products derived
19 from this software without specific prior written permission.
20
21 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22 "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23 LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
24 FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
25 COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
26 INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
27 BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
29 CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
30 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
31 ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
32 POSSIBILITY OF SUCH DAMAGE.
33*/
34
35#ifndef __BC_FACET_NORMAL_KERNEL__
36#define __BC_FACET_NORMAL_KERNEL__
37
43#define coef_normal_area_idx(i, j, k, l, lx, nf) \
44 (((i) + (lx) * (((j) - 1) + (lx) * (((k) - 1) + (nf) * (((l) - 1))))) - 1)
45
51void nonlinear_index(const int idx, const int lx, int *index) {
52 const int idx2 = idx -1;
53 index[3] = idx2/(lx * lx * lx) ;
54 index[2] = (idx2 - (lx*lx*lx)*index[3])/(lx * lx);
55 index[1] = (idx2 - (lx*lx*lx)*index[3] - (lx*lx) * index[2]) / lx;
56 index[0] = (idx2 - (lx*lx*lx)*index[3] - (lx*lx) * index[2]) - lx*index[1];
57 index[0]++;
58 index[1]++;
59 index[2]++;
60 index[3]++;
61}
62
63
67template< typename T >
70 const int * __restrict__ facet,
72 T * __restrict__ y,
73 T * __restrict__ z,
74 const T * __restrict__ u,
75 const T * __restrict__ v,
76 const T * __restrict__ w,
77 const T * __restrict__ nx,
78 const T * __restrict__ ny,
79 const T * __restrict__ nz,
80 const T * __restrict__ area,
81 const int lx,
82 const int m) {
83 int index[4];
84 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
85 const int str = blockDim.x * gridDim.x;
86
87 for (int i = (idx + 1); i < m; i += str) {
88 const int k = (msk[i] - 1);
89 const int f = (facet[i]);
90 nonlinear_index(msk[i], lx, index);
91
92
93 switch(f) {
94 case 1:
95 case 2:
96 {
97 const int na_idx = coef_normal_area_idx(index[1], index[2],
98 f, index[3], lx, 6);
99 x[k] = u[k] * nx[na_idx] * area[na_idx];
100 y[k] = v[k] * ny[na_idx] * area[na_idx];
101 z[k] = w[k] * nz[na_idx] * area[na_idx];
102 break;
103 }
104 case 3:
105 case 4:
106 {
107 const int na_idx = coef_normal_area_idx(index[0], index[2],
108 f, index[3], lx, 6);
109 x[k] = u[k] * nx[na_idx] * area[na_idx];
110 y[k] = v[k] * ny[na_idx] * area[na_idx];
111 z[k] = w[k] * nz[na_idx] * area[na_idx];
112 break;
113 }
114 case 5:
115 case 6:
116 {
117 const int na_idx = coef_normal_area_idx(index[0], index[1],
118 f, index[3], lx, 6);
119 x[k] = u[k] * nx[na_idx] * area[na_idx];
120 y[k] = v[k] * ny[na_idx] * area[na_idx];
121 z[k] = w[k] * nz[na_idx] * area[na_idx];
122 break;
123 }
124 }
125 }
126}
127
128#endif // __BC_FACET_NORMAL_KERNEL__
__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__ u
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ v
__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)
__device__ void nonlinear_index(const int idx, const int lx, int *index)
__global__ void facet_normal_apply_surfvec_kernel(const int *__restrict__ msk, const int *__restrict__ facet, T *__restrict__ x, T *__restrict__ y, T *__restrict__ z, const T *__restrict__ u, const T *__restrict__ v, const T *__restrict__ w, const T *__restrict__ nx, const T *__restrict__ ny, const T *__restrict__ nz, const T *__restrict__ area, const int lx, const int m)
#define coef_normal_area_idx(i, j, k, l, lx, nf)