Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
compressible_ops_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 __FLUID_COMPRESSIBLE_OPS_KERNEL_H__
36#define __FLUID_COMPRESSIBLE_OPS_KERNEL_H__
37
39
43template< typename T>
45 const T * __restrict__ u,
46 const T * __restrict__ v,
47 const T * __restrict__ w,
48 const T gamma,
49 const T * __restrict__ p,
50 const T * __restrict__ rho,
51 const int n) {
52
53 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
54 const int str = blockDim.x * gridDim.x;
55
56 for (int i = idx; i < n; i += str) {
57 T vel_mag = sqrt(u[i]*u[i] + v[i]*v[i] + w[i]*w[i]);
58 T sound_speed = sqrt(gamma * p[i] / rho[i]);
59 max_wave_speed[i] = vel_mag + sound_speed;
60 }
61
62}
63
67template< typename T>
69 const T * __restrict__ p,
70 const T * __restrict__ rho,
71 const T gamma,
72 const int n) {
73
74 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
75 const int str = blockDim.x * gridDim.x;
76
77 for (int i = idx; i < n; i += str) {
78 // S = 1/(gamma-1) * rho * (log(p) - gamma * log(rho))
79 T log_p = log(p[i]);
80 T log_rho = log(rho[i]);
81 S[i] = (1.0 / (gamma - 1.0)) * rho[i] * (log_p - gamma * log_rho);
82 }
83
84}
85
89template< typename T>
93 const T* __restrict__ m_x,
94 const T* __restrict__ m_y,
95 const T* __restrict__ m_z,
96 const T* __restrict__ rho,
97 const int n) {
98
99
100 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
101 const int str = blockDim.x * gridDim.x;
102
103 for (int i = idx; i < n; i += str) {
104 u[i] = m_x[i] / rho[i];
105 v[i] = m_y[i] / rho[i];
106 w[i] = m_z[i] / rho[i];
107 }
108
109}
110
111
115template< typename T>
117 T* __restrict__ m_y,
118 T* __restrict__ m_z,
119 T* __restrict__ p,
121 const T* __restrict__ u,
122 const T* __restrict__ v,
123 const T* __restrict__ w,
124 const T* __restrict__ E,
125 const T* __restrict__ rho,
126 const T gamma,
127 const int n) {
128
129 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
130 const int str = blockDim.x * gridDim.x;
131
132 for (int i = idx; i < n; i += str) {
133 m_x[i] = u[i] * rho[i];
134 m_y[i] = v[i] * rho[i];
135 m_z[i] = w[i] * rho[i];
136
137 /* Update p = (gamma - 1) * (E - 0.5 * rho * (u^2 + v^2 + w^2)) */
138 const real tmp = 0.5 * rho[i] * (u[i]*u[i] + v[i]*v[i] + w[i]*w[i]);
139 p[i] = (gamma - 1.0) * (E[i] - tmp);
140 ruvw[i] = tmp;
141 }
142}
143
144#define MAX(a,b) (((a)>(b))?(a):(b))
145
149template< typename T>
151 T* __restrict__ p,
152 const T* __restrict__ ruvw,
153 const T gamma,
154 const int n) {
155
156 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
157 const int str = blockDim.x * gridDim.x;
158
159 for (int i = idx; i < n; i += str) {
160 /* Ensure pressure is positive */
161 p[i] = MAX(p[i], 1e-12);
162 /* E = p / (gamma - 1) + 0.5 * rho * (u^2 + v^2 + w^2) */
163 E[i] = p[i] * (1.0 / (gamma - 1.0)) + ruvw[i];
164 }
165}
166
167#endif // __FLUID_COMPRESSIBLE_OPS_KERNEL_H__
__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
const int e
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ v
__global__ void compute_max_wave_speed_kernel(T *__restrict__ max_wave_speed, const T *__restrict__ u, const T *__restrict__ v, const T *__restrict__ w, const T gamma, const T *__restrict__ p, const T *__restrict__ rho, const int n)
__global__ void compute_entropy_kernel(T *__restrict__ S, const T *__restrict__ p, const T *__restrict__ rho, const T gamma, const int n)
__global__ void update_mxyz_p_ruvw_kernel(T *__restrict__ m_x, T *__restrict__ m_y, T *__restrict__ m_z, T *__restrict__ p, T *__restrict__ ruvw, const T *__restrict__ u, const T *__restrict__ v, const T *__restrict__ w, const T *__restrict__ E, const T *__restrict__ rho, const T gamma, const int n)
__global__ void update_e_kernel(T *__restrict__ E, T *__restrict__ p, const T *__restrict__ ruvw, const T gamma, const int n)
__global__ void update_uvw_kernel(T *__restrict__ u, T *__restrict__ v, T *__restrict__ w, const T *__restrict__ m_x, const T *__restrict__ m_y, const T *__restrict__ m_z, const T *__restrict__ rho, const int n)
#define MAX(a, b)
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
double real