Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
wale_nut_kernel.h
Go to the documentation of this file.
1#ifndef __COMMON_WALE_NUT_KERNEL_H__
2#define __COMMON_WALE_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__ g12,
45 const T * __restrict__ g13,
46 const T * __restrict__ g21,
47 const T * __restrict__ g22,
48 const T * __restrict__ g23,
49 const T * __restrict__ g31,
50 const T * __restrict__ g32,
51 const T * __restrict__ g33,
52 const T * __restrict__ delta,
53 T * __restrict__ nut,
54 const T * __restrict__ mult,
55 const T c,
56 const T eps,
57 const int n){
58
59 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
60 const int str = blockDim.x * gridDim.x;
61
62 T s11, s22, s33;
63 T s12, s13, s23;
67 T sd11, sd22, sd33;
68 T sd12, sd13, sd23;
70
71 for (int i = idx; i < n; i += str) {
72
73 const T g11_r = g11[i];
74 const T g12_r = g12[i];
75 const T g13_r = g13[i];
76 const T g21_r = g21[i];
77 const T g22_r = g22[i];
78 const T g23_r = g23[i];
79 const T g31_r = g31[i];
80 const T g32_r = g32[i];
81 const T g33_r = g33[i];
82 const T delta_r = delta[i];
83
84 s11 = g11_r;
85 s22 = g22_r;
86 s33 = g33_r;
87 s12 = 0.5 * (g12_r + g21_r);
88 s13 = 0.5 * (g13_r + g31_r);
89 s23 = 0.5 * (g23_r + g32_r);
90
92 + g12_r * g21_r
93 + g13_r * g31_r;
95 + g12_r * g22_r
96 + g13_r * g32_r;
98 + g12_r * g23_r
99 + g13_r * g33_r;
100 gsqr_21 = g21_r * g11_r
101 + g22_r * g21_r
102 + g23_r * g31_r;
103 gsqr_22 = g21_r * g12_r
104 + g22_r * g22_r
105 + g23_r * g32_r;
106 gsqr_23 = g21_r * g13_r
107 + g22_r * g23_r
108 + g23_r * g33_r;
109 gsqr_31 = g31_r * g11_r
110 + g32_r * g21_r
111 + g33_r * g31_r;
112 gsqr_32 = g31_r * g12_r
113 + g32_r * g22_r
114 + g33_r * g32_r;
115 gsqr_33 = g31_r * g13_r
116 + g32_r * g23_r
117 + g33_r * g33_r;
118
119 sd11 = gsqr_11 - ( (gsqr_11 + gsqr_22 + gsqr_33) / 3.0);
120 sd22 = gsqr_22 - ( (gsqr_11 + gsqr_22 + gsqr_33) / 3.0);
121 sd33 = gsqr_33 - ( (gsqr_11 + gsqr_22 + gsqr_33) / 3.0);
122 sd12 = 0.5 * (gsqr_12 + gsqr_21);
123 sd13 = 0.5 * (gsqr_13 + gsqr_31);
124 sd23 = 0.5 * (gsqr_23 + gsqr_32);
125
127 2.0 * (sd12*sd12 + sd13*sd13 + sd23*sd23);
128 Sij_Sij = s11*s11 + s22*s22 + s33*s33 +
129 2.0 * (s12*s12 + s13*s13 + s23*s23);
133 eps);
134
135
136 nut[i] = c * c * delta_r * delta_r
137 * OP_wale * mult[i];
138 }
139}
140#endif // __COMMON_WALE_NUT_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 wale_nut_compute(const T *__restrict__ g11, const T *__restrict__ g12, const T *__restrict__ g13, const T *__restrict__ g21, const T *__restrict__ g22, const T *__restrict__ g23, const T *__restrict__ g31, const T *__restrict__ g32, const T *__restrict__ g33, const T *__restrict__ delta, T *__restrict__ nut, const T *__restrict__ mult, const T c, const T eps, const int n)
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g23
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g22
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g13
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g12
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g33
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ g11
#define max(a, b)
Definition tensor.cu:40