Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
spalding_kernel.h
Go to the documentation of this file.
1#ifndef __COMMON_SPALDING_KERNEL_H__
2#define __COMMON_SPALDING_KERNEL_H__
3/*
4 Copyright (c) 2025, 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
40#include <cmath>
41#include <algorithm>
42template<typename T>
43__device__ T solve(const T u, const T y, const T guess, const T nu,
44 const T kappa, const T B);
45
49template<typename T>
51 const T * __restrict__ v_d,
52 const T * __restrict__ w_d,
53 const int * __restrict__ ind_r_d,
54 const int * __restrict__ ind_s_d,
55 const int * __restrict__ ind_t_d,
56 const int * __restrict__ ind_e_d,
57 const T * __restrict__ n_x_d,
58 const T * __restrict__ n_y_d,
59 const T * __restrict__ n_z_d,
60 const T * __restrict__ nu_d,
61 const T * __restrict__ rho_w_d,
62 const T * __restrict__ h_d,
66 const int n_nodes,
67 const int lx,
68 const T kappa,
69 const T B,
70 const int tstep) {
71
72 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
73 const int str = blockDim.x * gridDim.x;
74 for (int i = idx; i < n_nodes; i += str) {
75 // Sample the velocity
76 const int index = (ind_e_d[i] - 1) * lx * lx * lx +
77 (ind_t_d[i] - 1) * lx * lx +
78 (ind_s_d[i] - 1) * lx +
79 (ind_r_d[i] - 1);
80
81 T ui = u_d[index];
82 T vi = v_d[index];
83 T wi = w_d[index];
84 T rho = rho_w_d[i];
85
86 // Load normal vectors and wall shear stress values once
87 T nx = n_x_d[i];
88 T ny = n_y_d[i];
89 T nz = n_z_d[i];
90 T h = h_d[i];
91
92 // Project on tangential direction
93 T normu = ui * nx + vi * ny + wi * nz;
94
95 ui -= normu * nx;
96 vi -= normu * ny;
97 wi -= normu * nz;
98
99 T magu = sqrt(ui * ui + vi * vi + wi * wi);
100
101 // Get initial guess for Newton solver
102 T guess;
103 if (tstep == 1) {
104 guess = sqrt(magu * nu_d[i] / h);
105 } else {
106 guess = tau_x_d[i] * tau_x_d[i] +
107 tau_y_d[i] * tau_y_d[i] +
108 tau_z_d[i] * tau_z_d[i];
109 guess = sqrt(sqrt(guess));
110 }
111
112 // Solve for utau using Newton's method
113 T utau = solve(magu, h, guess, nu_d[i], kappa, B);
114
115 // Distribute according to the velocity vector
116 tau_x_d[i] = -rho * utau * utau * ui / magu;
117 tau_y_d[i] = -rho * utau * utau * vi / magu;
118 tau_z_d[i] = -rho * utau * utau * wi / magu;
119 }
120}
121
125template<typename T>
126__device__ T solve(const T u, const T y, const T guess, const T nu,
127 const T kappa, const T B) {
128 T utau = guess;
129 T yp, up, f, df, old, error;
130 const int maxiter = 100;
131
132 for (int k = 0; k < maxiter; ++k) {
133 up = u / utau;
134 yp = y * utau / nu;
135 old = utau;
136
137 // Evaluate function and its derivative
138 f = (up + exp(-kappa * B) *
139 (exp(kappa * up) - 1.0 - kappa * up -
140 0.5 * (kappa * up) * (kappa * up) -
141 (1.0 / 6.0) * (kappa * up) * (kappa * up) * (kappa * up)) -
142 yp);
143
144 df = (-y / nu - u / (utau * utau) -
145 kappa * up / utau * exp(-kappa * B) *
146 (exp(kappa * up) - 1.0 - kappa * up -
147 0.5 * (kappa * up) * (kappa * up)));
148
149 // Update solution
150 utau -= f / df;
151
152 error = fabs((old - utau) / old);
153
154 if (error < 1e-3) {
155 break;
156 }
157 }
158
159 return utau;
160}
161
162#endif // __COMMON_SPALDING_KERNEL_H__
const int i
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ u
const int e
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
__global__ void spalding_compute(const T *__restrict__ u_d, const T *__restrict__ v_d, const T *__restrict__ w_d, const int *__restrict__ ind_r_d, const int *__restrict__ ind_s_d, const int *__restrict__ ind_t_d, const int *__restrict__ ind_e_d, const T *__restrict__ n_x_d, const T *__restrict__ n_y_d, const T *__restrict__ n_z_d, const T *__restrict__ nu_d, const T *__restrict__ rho_w_d, const T *__restrict__ h_d, T *__restrict__ tau_x_d, T *__restrict__ tau_y_d, T *__restrict__ tau_z_d, const int n_nodes, const int lx, const T kappa, const T B, const int tstep)
solve(case_descr)
Definition intf.py:323