Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
euler_res_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#ifndef __FLUID_EULER_RES_KERNEL__
35#define __FLUID_EULER_RES_KERNEL__
36
37template< typename T >
39 const T * __restrict__ Binv,
40 const T * __restrict__ lap_sol,
41 const T * __restrict__ h,
42 const T c_avisc,
43 const int n) {
44
45 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
46 const int str = blockDim.x * gridDim.x;
47
48 for (int i = idx; i < n; i += str) {
49 rhs[i] = -rhs[i] - c_avisc * h[i] * Binv[i] * lap_sol[i];
50 }
51}
52
53template< typename T >
55 T * __restrict__ f_y,
56 T * __restrict__ f_z,
57 const T * __restrict__ m_x,
58 const T * __restrict__ m_y,
59 const T * __restrict__ m_z,
60 const T * __restrict__ rho_field,
61 const T * __restrict__ p,
62 const int n) {
63
64 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
65 const int str = blockDim.x * gridDim.x;
66
67 for (int i = idx; i < n; i += str) {
68 f_x[i] = m_x[i] * m_x[i] / rho_field[i] + p[i];
69 f_y[i] = m_x[i] * m_y[i] / rho_field[i];
70 f_z[i] = m_x[i] * m_z[i] / rho_field[i];
71 }
72}
73
74template< typename T >
76 T * __restrict__ f_y,
77 T * __restrict__ f_z,
78 const T * __restrict__ m_x,
79 const T * __restrict__ m_y,
80 const T * __restrict__ m_z,
81 const T * __restrict__ rho_field,
82 const T * __restrict__ p,
83 const int n) {
84
85 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
86 const int str = blockDim.x * gridDim.x;
87
88 for (int i = idx; i < n; i += str) {
89 f_x[i] = m_y[i] * m_x[i] / rho_field[i];
90 f_y[i] = m_y[i] * m_y[i] / rho_field[i] + p[i];
91 f_z[i] = m_y[i] * m_z[i] / rho_field[i];
92 }
93}
94
95template< typename T >
97 T * __restrict__ f_y,
98 T * __restrict__ f_z,
99 const T * __restrict__ m_x,
100 const T * __restrict__ m_y,
101 const T * __restrict__ m_z,
102 const T * __restrict__ rho_field,
103 const T * __restrict__ p,
104 const int n) {
105
106 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
107 const int str = blockDim.x * gridDim.x;
108
109 for (int i = idx; i < n; i += str) {
110 f_x[i] = m_z[i] * m_x[i] / rho_field[i];
111 f_y[i] = m_z[i] * m_y[i] / rho_field[i];
112 f_z[i] = m_z[i] * m_z[i] / rho_field[i] + p[i];
113 }
114}
115
116template< typename T >
118 T * __restrict__ f_y,
119 T * __restrict__ f_z,
120 const T * __restrict__ m_x,
121 const T * __restrict__ m_y,
122 const T * __restrict__ m_z,
123 const T * __restrict__ rho_field,
124 const T * __restrict__ p,
125 const T * __restrict__ E,
126 const int n) {
127
128 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
129 const int str = blockDim.x * gridDim.x;
130
131 for (int i = idx; i < n; i += str) {
132 f_x[i] = (E[i] + p[i]) * m_x[i] / rho_field[i];
133 f_y[i] = (E[i] + p[i]) * m_y[i] / rho_field[i];
134 f_z[i] = (E[i] + p[i]) * m_z[i] / rho_field[i];
135 }
136}
137
138template< typename T >
144 const T * __restrict__ mult,
145 const int n) {
146
147 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
148 const int str = blockDim.x * gridDim.x;
149
150 for (int i = idx; i < n; i += str) {
151 rhs_rho[i] = rhs_rho[i] * mult[i];
152 rhs_m_x[i] = rhs_m_x[i] * mult[i];
153 rhs_m_y[i] = rhs_m_y[i] * mult[i];
154 rhs_m_z[i] = rhs_m_z[i] * mult[i];
155 rhs_E[i] = rhs_E[i] * mult[i];
156 }
157}
158
159template< typename T >
161 T * __restrict__ m_x,
162 T * __restrict__ m_y,
163 T * __restrict__ m_z,
164 T * __restrict__ E,
165 const T * __restrict__ k_rho_i,
166 const T * __restrict__ k_m_x_i,
167 const T * __restrict__ k_m_y_i,
168 const T * __restrict__ k_m_z_i,
169 const T * __restrict__ k_E_i,
170 const T dt,
171 const T c,
172 const int n) {
173
174 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
175 const int str = blockDim.x * gridDim.x;
176
177 for (int i = idx; i < n; i += str) {
178 rho[i] = rho[i] + dt * c * k_rho_i[i];
179 m_x[i] = m_x[i] + dt * c * k_m_x_i[i];
180 m_y[i] = m_y[i] + dt * c * k_m_y_i[i];
181 m_z[i] = m_z[i] + dt * c * k_m_z_i[i];
182 E[i] = E[i] + dt * c * k_E_i[i];
183 }
184}
185
186#endif // __FLUID_EULER_RES_KERNEL__
const int i
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dt
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
__global__ void euler_res_part_visc_kernel(T *__restrict__ rhs, const T *__restrict__ Binv, const T *__restrict__ lap_sol, const T *__restrict__ h, const T c_avisc, const int n)
__global__ void euler_res_part_rk_sum_kernel(T *__restrict__ rho, T *__restrict__ m_x, T *__restrict__ m_y, T *__restrict__ m_z, T *__restrict__ E, const T *__restrict__ k_rho_i, const T *__restrict__ k_m_x_i, const T *__restrict__ k_m_y_i, const T *__restrict__ k_m_z_i, const T *__restrict__ k_E_i, const T dt, const T c, const int n)
__global__ void euler_res_part_my_flux_kernel(T *__restrict__ f_x, T *__restrict__ f_y, T *__restrict__ f_z, const T *__restrict__ m_x, const T *__restrict__ m_y, const T *__restrict__ m_z, const T *__restrict__ rho_field, const T *__restrict__ p, const int n)
__global__ void euler_res_part_E_flux_kernel(T *__restrict__ f_x, T *__restrict__ f_y, T *__restrict__ f_z, const T *__restrict__ m_x, const T *__restrict__ m_y, const T *__restrict__ m_z, const T *__restrict__ rho_field, const T *__restrict__ p, const T *__restrict__ E, const int n)
__global__ void euler_res_part_mz_flux_kernel(T *__restrict__ f_x, T *__restrict__ f_y, T *__restrict__ f_z, const T *__restrict__ m_x, const T *__restrict__ m_y, const T *__restrict__ m_z, const T *__restrict__ rho_field, const T *__restrict__ p, const int n)
__global__ void euler_res_part_coef_mult_kernel(T *__restrict__ rhs_rho, T *__restrict__ rhs_m_x, T *__restrict__ rhs_m_y, T *__restrict__ rhs_m_z, T *__restrict__ rhs_E, const T *__restrict__ mult, const int n)
__global__ void euler_res_part_mx_flux_kernel(T *__restrict__ f_x, T *__restrict__ f_y, T *__restrict__ f_z, const T *__restrict__ m_x, const T *__restrict__ m_y, const T *__restrict__ m_z, const T *__restrict__ rho_field, const T *__restrict__ p, const int n)