Neko  0.9.0
A portable framework for high-order spectral element flow simulations
makebdf_kernel.h
Go to the documentation of this file.
1 /*
2  Copyright (c) 2022, 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 __COMMON_MAKEBDF_KERNEL__
36 #define __COMMON_MAKEBDF_KERNEL__
37 
38 template< typename T >
39 __global__ void makebdf_kernel(const T * __restrict__ ulag1,
40  const T * __restrict__ ulag2,
41  const T * __restrict__ vlag1,
42  const T * __restrict__ vlag2,
43  const T * __restrict__ wlag1,
44  const T * __restrict__ wlag2,
45  T * __restrict__ bfx,
46  T * __restrict__ bfy,
47  T * __restrict__ bfz,
48  const T * __restrict__ u,
49  const T * __restrict__ v,
50  const T * __restrict__ w,
51  const T * __restrict__ B,
52  const T rho,
53  const T dt,
54  const T bd2,
55  const T bd3,
56  const T bd4,
57  const int nbd,
58  const int n) {
59 
60  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
61  const int str = blockDim.x * gridDim.x;
62 
63  for (int i = idx; i < n; i += str) {
64  T tb1_val = u[i] * B[i] * bd2;
65  T tb2_val = v[i] * B[i] * bd2;
66  T tb3_val = w[i] * B[i] * bd2;
67 
68  T ta1_val = ulag1[i] * B[i] * bd3;
69  T ta2_val = vlag1[i] * B[i] * bd3;
70  T ta3_val = wlag1[i] * B[i] * bd3;
71 
72  tb1_val += ta1_val;
73  tb2_val += ta2_val;
74  tb3_val += ta3_val;
75 
76  if (nbd == 3) {
77  tb1_val += ulag2[i] * B[i] * bd4;
78  tb2_val += vlag2[i] * B[i] * bd4;
79  tb3_val += wlag2[i] * B[i] * bd4;
80  }
81 
82  bfx[i] = bfx[i] + tb1_val * (rho / dt);
83  bfy[i] = bfy[i] + tb2_val * (rho / dt);
84  bfz[i] = bfz[i] + tb3_val * (rho / dt);
85  }
86 
87 }
88 
89 template< typename T >
90 __global__ void scalar_makebdf_kernel(const T * __restrict__ s_lag,
91  const T * __restrict__ s_laglag,
92  T * __restrict__ fs,
93  const T * __restrict__ s,
94  const T * __restrict__ B,
95  const T rho,
96  const T dt,
97  const T bd2,
98  const T bd3,
99  const T bd4,
100  const int nbd,
101  const int n) {
102 
103  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
104  const int str = blockDim.x * gridDim.x;
105 
106  for (int i = idx; i < n; i += str) {
107  T tb1_val = s[i] * B[i] * bd2;
108 
109  T ta1_val = s_lag[i] * B[i] * bd3;
110 
111  tb1_val += ta1_val;
112 
113  if (nbd == 3) {
114  tb1_val += s_laglag[i] * B[i] * bd4;
115  }
116 
117  fs[i] = fs[i] + tb1_val * (rho / dt);
118  }
119 
120 }
121 
122 #endif // __COMMON_MAKEBDF_KERNEL__
__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
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ v
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dt
Definition: cdtp_kernel.h:109
__global__ void scalar_makebdf_kernel(const T *__restrict__ s_lag, const T *__restrict__ s_laglag, T *__restrict__ fs, const T *__restrict__ s, const T *__restrict__ B, const T rho, const T dt, const T bd2, const T bd3, const T bd4, const int nbd, const int n)
__global__ void makebdf_kernel(const T *__restrict__ ulag1, const T *__restrict__ ulag2, const T *__restrict__ vlag1, const T *__restrict__ vlag2, const T *__restrict__ wlag1, const T *__restrict__ wlag2, T *__restrict__ bfx, T *__restrict__ bfy, T *__restrict__ bfz, const T *__restrict__ u, const T *__restrict__ v, const T *__restrict__ w, const T *__restrict__ B, const T rho, const T dt, const T bd2, const T bd3, const T bd4, const int nbd, const int n)