Neko  0.8.99
A portable framework for high-order spectral element flow simulations
prs_res_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 __FLUID_PRS_RES_KERNEL__
36 #define __FLUID_PRS_RES_KERNEL__
37 
38 template< typename T >
39 __global__ void prs_res_part1_kernel(T * __restrict__ ta1,
40  T * __restrict__ ta2,
41  T * __restrict__ ta3,
42  const T * __restrict__ wa1,
43  const T * __restrict__ wa2,
44  const T * __restrict__ wa3,
45  const T * __restrict__ f_u,
46  const T * __restrict__ f_v,
47  const T * __restrict__ f_w,
48  const T * __restrict__ B,
49  T * __restrict__ h1,
50  const T mu,
51  const T rho,
52  const int n) {
53 
54  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
55  const int str = blockDim.x * gridDim.x;
56  const T inv_rho = 1.0 / rho;
57 
58  for (int i = idx; i < n; i += str) {
59  h1[i] = inv_rho;
60  ta1[i] = (f_u[i] / rho) - ((wa1[i] * (mu / rho)) * B[i]);
61  ta2[i] = (f_v[i] / rho) - ((wa2[i] * (mu / rho)) * B[i]);
62  ta3[i] = (f_w[i] / rho) - ((wa3[i] * (mu / rho)) * B[i]);
63  }
64 
65 }
66 
67 
68 template< typename T >
69 __global__ void prs_res_part2_kernel(T * __restrict__ p_res,
70  const T * __restrict__ wa1,
71  const T * __restrict__ wa2,
72  const T * __restrict__ wa3,
73  const int n) {
74 
75  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
76  const int str = blockDim.x * gridDim.x;
77 
78  for (int i = idx; i < n; i += str) {
79  p_res[i] = (-p_res[i]) + (wa1[i] + wa2[i] + wa3[i]);
80  }
81 }
82 
83 template< typename T >
84 __global__ void prs_res_part3_kernel(T * __restrict__ p_res,
85  const T * __restrict__ ta1,
86  const T * __restrict__ ta2,
87  const T * __restrict__ ta3,
88  const T dtbd,
89  const int n) {
90 
91  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
92  const int str = blockDim.x * gridDim.x;
93 
94  for (int i = idx; i < n; i += str) {
95  p_res[i] = p_res[i] - (dtbd * (ta1[i] + ta2[i] + ta3[i]));
96  }
97 }
98 
99 #endif // __FLUID_PRS_RES_KERNEL__
const int i
Definition: cdtp_kernel.h:128
__global__ void prs_res_part3_kernel(T *__restrict__ p_res, const T *__restrict__ ta1, const T *__restrict__ ta2, const T *__restrict__ ta3, const T dtbd, const int n)
__global__ void prs_res_part2_kernel(T *__restrict__ p_res, const T *__restrict__ wa1, const T *__restrict__ wa2, const T *__restrict__ wa3, const int n)
__global__ void prs_res_part1_kernel(T *__restrict__ ta1, T *__restrict__ ta2, T *__restrict__ ta3, const T *__restrict__ wa1, const T *__restrict__ wa2, const T *__restrict__ wa3, const T *__restrict__ f_u, const T *__restrict__ f_v, const T *__restrict__ f_w, const T *__restrict__ B, T *__restrict__ h1, const T mu, const T rho, const int n)
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ h1