Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
rough_log_law_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 ROUGH_LOG_LAW_KERNEL_H
36#define ROUGH_LOG_LAW_KERNEL_H
37
38#include <cmath>
39#include <algorithm>
40
44template<typename T>
46 const T* __restrict__ v_d,
47 const T* __restrict__ w_d,
48 const int* __restrict__ ind_r_d,
49 const int* __restrict__ ind_s_d,
50 const int* __restrict__ ind_t_d,
51 const int* __restrict__ ind_e_d,
52 const T* __restrict__ n_x_d,
53 const T* __restrict__ n_y_d,
54 const T* __restrict__ n_z_d,
55 const T* __restrict__ h_d,
59 const int n_nodes,
60 const int lx,
61 const T kappa,
62 const T * __restrict__ rho_w_d,
63 const T B,
64 const T z0) {
65
66 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
67 const int str = blockDim.x * gridDim.x;
68 for (int i = idx; i < n_nodes; i += str) {
69 // Sample the velocity
70 const int index = (ind_e_d[i] - 1) * lx * lx * lx +
71 (ind_t_d[i] - 1) * lx * lx +
72 (ind_s_d[i] - 1) * lx +
73 (ind_r_d[i] - 1);
74
75 T ui = u_d[index];
76 T vi = v_d[index];
77 T wi = w_d[index];
78 T rho = rho_w_d[i];
79
80 // Load normal vectors and wall shear stress values once
81 T nx = n_x_d[i];
82 T ny = n_y_d[i];
83 T nz = n_z_d[i];
84 T h = h_d[i];
85
86 // Project on tangential direction
87 T normu = ui * nx + vi * ny + wi * nz;
88
89 ui -= normu * nx;
90 vi -= normu * ny;
91 wi -= normu * nz;
92
93 T magu = sqrt(ui * ui + vi * vi + wi * wi);
94
95 // Compute the wall shear stress using the rough log-law
96 T utau = 0.0;
97 if (h > z0) {
98 utau = (magu - B) * kappa / log(h / z0);
99 }
100
101 // Distribute according to the velocity vector
102 tau_x_d[i] = -rho *utau * utau * ui / magu;
103 tau_y_d[i] = -rho *utau * utau * vi / magu;
104 tau_z_d[i] = -rho *utau * utau * wi / magu;
105 }
106}
107
108#endif // ROUGH_LOG_LAW_KERNEL_H
const int i
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
__global__ void rough_log_law_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__ 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 *__restrict__ rho_w_d, const T B, const T z0)