Neko  0.9.99
A portable framework for high-order spectral element flow simulations
gradient_jump_penalty_kernel.h
Go to the documentation of this file.
1 #ifndef __COMMON_GRADIENT_JUMP_PENALTY_KERNEL_H__
2 #define __COMMON_GRADIENT_JUMP_PENALTY_KERNEL_H__
3 /*
4  Copyright (c) 2024, 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 template< typename T>
41 __global__ void pick_facet_value_hex_kernel(T * __restrict__ b,
42  T *__restrict__ a,
43  const int nx) {
44 
45  const int idx = threadIdx.x;
46  const int nx2 = nx+2;
47  const int el2 = blockIdx.x*nx2*nx2*nx2;
48  const int el = blockIdx.x*nx*nx*nx;
49  for(int ijk = idx; ijk < nx*nx*nx; ijk += blockDim.x){
50  const int jk = ijk / nx;
51  const int i = ijk - jk * nx;
52  const int k = jk / nx;
53  const int j = jk - k * nx;
54  if(i == 0){
55  b[0+(j+1)*nx2+(k+1)*nx2*nx2+el2] = a[ijk+el];
56  }
57  if(i == nx-1){
58  b[nx2-1+(j+1)*nx2+(k+1)*nx2*nx2+el2] = a[ijk+el];
59  }
60  if(j == 0){
61  b[(i+1)+0*nx2+(k+1)*nx2*nx2+el2] = a[ijk+el];
62  }
63  if(j == nx-1){
64  b[(i+1)+(nx2-1)*nx2+(k+1)*nx2*nx2+el2] = a[ijk+el];
65  }
66  if(k == 0){
67  b[(i+1)+(j+1)*nx2+0*nx2*nx2+el2] = a[ijk+el];
68  }
69  if(k == nx-1){
70  b[(i+1)+(j+1)*nx2+(nx2-1)*nx2*nx2+el2] = a[ijk+el];
71  }
72  }
73 }
74 
78 template< typename T>
79 __global__ void gradient_jump_penalty_finalize_kernel(T * __restrict__ penalty_d,
80  T *__restrict__ penalty_facet_d,
81  T *__restrict__ dphidxi_d,
82  const int nx) {
83 
84  const int idx = threadIdx.x;
85  const int nx2 = nx+2;
86  const int el2 = blockIdx.x*nx2*nx2*nx2;
87  const int el = blockIdx.x*nx*nx*nx;
88  for(int ijk = idx; ijk < nx*nx*nx; ijk += blockDim.x){
89  const int jk = ijk / nx;
90  const int i = ijk - jk * nx;
91  const int k = jk / nx;
92  const int j = jk - k * nx;
93  penalty_d[ijk+el] = penalty_facet_d[0+(j+1)*nx2+(k+1)*nx2*nx2+el2] \
94  * dphidxi_d[0+i*nx] + \
95  penalty_facet_d[(nx2-1)+(j+1)*nx2+(k+1)*nx2*nx2+el2] \
96  * dphidxi_d[nx-1+i*nx] + \
97  penalty_facet_d[(i+1)+0*nx2+(k+1)*nx2*nx2+el2] \
98  * dphidxi_d[0+j*nx] + \
99  penalty_facet_d[(i+1)+(nx2-1)*nx2+(k+1)*nx2*nx2+el2] \
100  * dphidxi_d[nx-1+j*nx] + \
101  penalty_facet_d[(i+1)+(j+1)*nx2+0*nx2*nx2+el2] \
102  * dphidxi_d[0+k*nx] + \
103  penalty_facet_d[(i+1)+(j+1)*nx2+(nx2-1)*nx2*nx2+el2] \
104  * dphidxi_d[nx-1+k*nx];
105 
106  }
107 }
108 
109 #endif // __COMMON_GRADIENT_JUMP_PENALTY_KERNEL_H__
110 
const int i
const int j
__global__ void gradient_jump_penalty_finalize_kernel(T *__restrict__ penalty_d, T *__restrict__ penalty_facet_d, T *__restrict__ dphidxi_d, const int nx)
__global__ void pick_facet_value_hex_kernel(T *__restrict__ b, T *__restrict__ a, const int nx)