Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
neumann_kernel.h
Go to the documentation of this file.
1/*
2 Copyright (c) 2025, 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_NEUMANN_KERNEL__
36#define __BC_NEUMANN_KERNEL__
37
38#include "bc_utils.h"
39
45#define coef_normal_area_idx(i, j, k, l, lx, nf) \
46 (((i) + (lx) * (((j) - 1) + (lx) * (((k) - 1) + (nf) * (((l) - 1))))) - 1)
47
48
52template< typename T >
55 const int * __restrict__ facet,
57 const T * __restrict__ flux,
58 const T * __restrict__ area,
59 const int lx,
60 const int m) {
61 int index[4];
62 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
63 const int str = blockDim.x * gridDim.x;
64
65 for (int i = (idx + 1); i < m; i += str) {
66 const int k = (msk[i] - 1);
67 const int f = (facet[i]);
68 nonlinear_index(msk[i], lx, index);
69
70 switch(f) {
71 case 1:
72 case 2:
73 {
74 const int na_idx = coef_normal_area_idx(index[1], index[2],
75 f, index[3], lx, 6);
76 x[k] += flux[i-1] * area[na_idx];
77 break;
78 }
79 case 3:
80 case 4:
81 {
82 const int na_idx = coef_normal_area_idx(index[0], index[2],
83 f, index[3], lx, 6);
84 x[k] += flux[i-1] * area[na_idx];
85 break;
86 }
87 case 5:
88 case 6:
89 {
90 const int na_idx = coef_normal_area_idx(index[0], index[1],
91 f, index[3], lx, 6);
92 x[k] += flux[i-1] * area[na_idx];
93 break;
94 }
95 }
96 }
97}
98
102template< typename T >
105 const int * __restrict__ facet,
106 T * __restrict__ x,
107 T * __restrict__ y,
108 T * __restrict__ z,
109 const T * __restrict__ flux_x,
110 const T * __restrict__ flux_y,
111 const T * __restrict__ flux_z,
112 const T * __restrict__ area,
113 const int lx,
114 const int m) {
115 int index[4];
116 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
117 const int str = blockDim.x * gridDim.x;
118
119 for (int i = (idx + 1); i < m; i += str) {
120 const int k = msk[i] - 1;
121 const int f = facet[i];
122 nonlinear_index(msk[i], lx, index);
123
124 switch(f) {
125 case 1:
126 case 2:
127 {
128 const int na_idx = coef_normal_area_idx(index[1], index[2],
129 f, index[3], lx, 6);
130 x[k] += flux_x[i-1] * area[na_idx];
131 y[k] += flux_y[i-1] * area[na_idx];
132 z[k] += flux_z[i-1] * area[na_idx];
133 break;
134 }
135 case 3:
136 case 4:
137 {
138 const int na_idx = coef_normal_area_idx(index[0], index[2],
139 f, index[3], lx, 6);
140 x[k] += flux_x[i-1] * area[na_idx];
141 y[k] += flux_y[i-1] * area[na_idx];
142 z[k] += flux_z[i-1] * area[na_idx];
143 break;
144 }
145 case 5:
146 case 6:
147 {
148 const int na_idx = coef_normal_area_idx(index[0], index[1],
149 f, index[3], lx, 6);
150 x[k] += flux_x[i-1] * area[na_idx];
151 y[k] += flux_y[i-1] * area[na_idx];
152 z[k] += flux_z[i-1] * area[na_idx];
153 break;
154 }
155 }
156 }
157}
158
159#endif // __BC_NEUMANN_KERNEL__
const int i
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
__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 neumann_apply_vector_kernel(const int *__restrict__ msk, const int *__restrict__ facet, T *__restrict__ x, T *__restrict__ y, T *__restrict__ z, const T *__restrict__ flux_x, const T *__restrict__ flux_y, const T *__restrict__ flux_z, const T *__restrict__ area, const int lx, const int m)
__global__ void neumann_apply_scalar_kernel(const int *__restrict__ msk, const int *__restrict__ facet, T *__restrict__ x, const T *__restrict__ flux, const T *__restrict__ area, const int lx, const int m)
#define coef_normal_area_idx(i, j, k, l, lx, nf)