Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
deardorff_nut_kernel.h
Go to the documentation of this file.
1#ifndef __COMMON_DEARDORFF_NUT_KERNEL_H__
2#define __COMMON_DEARDORFF_NUT_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>
44 const T *__restrict__ dTdx,
45 const T *__restrict__ dTdy,
46 const T *__restrict__ dTdz,
47 const T * __restrict__ a11,
48 const T * __restrict__ a12,
49 const T * __restrict__ a13,
50 const T * __restrict__ a21,
51 const T * __restrict__ a22,
52 const T * __restrict__ a23,
53 const T * __restrict__ a31,
54 const T * __restrict__ a32,
55 const T * __restrict__ a33,
56 const T * __restrict__ delta,
57 T * __restrict__ nut,
58 T * __restrict__ temperature_alphat,
61 const T c_k,
62 const T T0,
63 const T g1,
64 const T g2,
65 const T g3,
66 const T eps,
67 const int n){
68
69 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
70 const int str = blockDim.x * gridDim.x;
71
72 T s11, s22, s33;
73 T s12, s13, s23;
74 T N2, l;
76
77 for (int i = idx; i < n; i += str) {
78 const T dTdx_r = dTdx[i];
79 const T dTdy_r = dTdy[i];
80 const T dTdz_r = dTdz[i];
81 const T a11_r = a11[i];
82 const T a12_r = a12[i];
83 const T a13_r = a13[i];
84 const T a21_r = a21[i];
85 const T a22_r = a22[i];
86 const T a23_r = a23[i];
87 const T a31_r = a31[i];
88 const T a32_r = a32[i];
89 const T a33_r = a33[i];
90 const T delta_r = delta[i];
91
92 if (TKE[i] < eps) {
93 TKE[i] = eps;
94 }
95 const T TKE_r = TKE[i];
96
97 N2 = (dTdx_r * g1 + dTdy_r * g2 + dTdz_r * g3) / T0;
98 if (N2 > 0.0) {
99 l = 0.76 * sqrt(TKE_r / N2);
100 l = min(l, delta_r);
101 }
102 else {
103 l = delta_r;
104 }
105
106 nut[i] = c_k * l * sqrt(TKE_r);
107 temperature_alphat[i] = (1.0 + 2.0 * l/delta_r) * nut[i];
108 TKE_alphat[i] = 2.0 * nut[i];
109
110 s11 = a11_r + a11_r;
111 s22 = a22_r + a22_r;
112 s33 = a33_r + a33_r;
113 s12 = a12_r + a21_r;
114 s13 = a13_r + a31_r;
115 s23 = a23_r + a32_r;
116
117 shear = nut[i] * (s11*a11_r + s12*a12_r + s13*a13_r
118 + s12*a21_r + s22*a22_r + s23*a23_r
119 + s13*a31_r + s23*a32_r + s33*a33_r);
120 buoyancy = - (dTdx_r * g1 +
121 dTdy_r * g2 +
122 dTdz_r * g3) * temperature_alphat[i] / T0;
123 dissipation = - (0.19 + 0.74 * l/delta_r) * sqrt(TKE_r*TKE_r*TKE_r) / l;
124
126 }
127}
128#endif // __COMMON_DEARDORFF_NUT_KERNEL_H__
const int i
__global__ void deardorff_nut_compute(T *__restrict__ TKE, const T *__restrict__ dTdx, const T *__restrict__ dTdy, const T *__restrict__ dTdz, const T *__restrict__ a11, const T *__restrict__ a12, const T *__restrict__ a13, const T *__restrict__ a21, const T *__restrict__ a22, const T *__restrict__ a23, const T *__restrict__ a31, const T *__restrict__ a32, const T *__restrict__ a33, const T *__restrict__ delta, T *__restrict__ nut, T *__restrict__ temperature_alphat, T *__restrict__ TKE_alphat, T *__restrict__ TKE_source, const T c_k, const T T0, const T g1, const T g2, const T g3, const T eps, const int n)
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)