Neko  0.9.99
A portable framework for high-order spectral element flow simulations
pipecg_kernel.h
Go to the documentation of this file.
1 #ifndef __KRYLOV_PIPECG_KERNEL_H__
2 #define __KRYLOV_PIPECG_KERNEL_H__
3 /*
4  Copyright (c) 2021-2022, 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 
38 
42 template< typename T >
43 __global__ void cg_update_xp_kernel(T * __restrict__ x,
44  T * __restrict__ p,
45  T ** __restrict__ u,
46  const T * alpha,
47  const T * beta,
48  const int p_cur,
49  const int p_space,
50  const int n) {
51 
52  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
53  const int str = blockDim.x * gridDim.x;
54 
55 
56  for (int i = idx; i < n; i+= str) {
57  T tmp = 0.0;
58  int p_prev = p_space;
59  for (int j = 0; j < p_cur; j ++) {
60  p[i] = beta[j]*p[i] + u[p_prev][i];
61  tmp += alpha[j]*p[i];
62  p_prev = j;
63  }
64  x[i] += tmp;
65  u[p_space][i] = u[p_space-1][i];
66  }
67 }
71 template< typename T >
72 __global__ void pipecg_vecops_kernel(T * __restrict__ p,
73  T * __restrict__ q,
74  T * __restrict__ r,
75  T * __restrict__ s,
76  T * __restrict__ u1,
77  T * __restrict__ u2,
78  T * __restrict__ w,
79  T * __restrict__ z,
80  T * __restrict__ ni,
81  T * __restrict__ mi,
82  const T alpha,
83  const T beta,
84  const T * mult,
85  T * buf_h1,
86  T * buf_h2,
87  T * buf_h3,
88  const int n) {
89 
90  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
91  const int str = blockDim.x * gridDim.x;
92 
93  const unsigned int lane = threadIdx.x % warpSize;
94  const unsigned int wid = threadIdx.x / warpSize;
95 
96  __shared__ T buf1[32];
97  __shared__ T buf2[32];
98  __shared__ T buf3[32];
99  T tmp1 = 0.0;
100  T tmp2 = 0.0;
101  T tmp3 = 0.0;
102 
103  for (int i = idx; i < n; i+= str) {
104  z[i] = beta * z[i] + ni[i];
105  q[i] = beta * q[i] + mi[i];
106  s[i] = beta * s[i] + w[i];
107  r[i] = r[i] - alpha * s[i];
108  u2[i] = u1[i] - alpha * q[i];
109  w[i] = w[i] - alpha * z[i];
110  tmp1 = tmp1 + r[i] * mult[i] * u2[i];
111  tmp2 = tmp2 + w[i] * mult[i] * u2[i];
112  tmp3 = tmp3 + r[i] * mult[i] * r[i];
113 
114  }
115 
116  tmp1 = reduce_warp<T>(tmp1);
117  tmp2 = reduce_warp<T>(tmp2);
118  tmp3 = reduce_warp<T>(tmp3);
119  if (lane == 0) {
120  buf1[wid] = tmp1;
121  buf2[wid] = tmp2;
122  buf3[wid] = tmp3;
123  }
124  __syncthreads();
125 
126  tmp1 = (threadIdx.x < blockDim.x / warpSize) ? buf1[lane] : 0;
127  tmp2 = (threadIdx.x < blockDim.x / warpSize) ? buf2[lane] : 0;
128  tmp3 = (threadIdx.x < blockDim.x / warpSize) ? buf3[lane] : 0;
129  if (wid == 0) {
130  tmp1 = reduce_warp<T>(tmp1);
131  tmp2 = reduce_warp<T>(tmp2);
132  tmp3 = reduce_warp<T>(tmp3);
133  }
134 
135  if (threadIdx.x == 0) {
136  buf_h1[blockIdx.x] = tmp1;
137  buf_h2[blockIdx.x] = tmp2;
138  buf_h3[blockIdx.x] = tmp3;
139  }
140 
141 }
142 
143 #endif // __KRYLOV_PIPECG_KERNEL_H__
__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
const int j
__syncthreads()
__global__ void const T *__restrict__ x
Definition: cdtp_kernel.h:106
__global__ void pipecg_vecops_kernel(T *__restrict__ p, T *__restrict__ q, T *__restrict__ r, T *__restrict__ s, T *__restrict__ u1, T *__restrict__ u2, T *__restrict__ w, T *__restrict__ z, T *__restrict__ ni, T *__restrict__ mi, const T alpha, const T beta, const T *mult, T *buf_h1, T *buf_h2, T *buf_h3, const int n)
Definition: pipecg_kernel.h:72
__global__ void cg_update_xp_kernel(T *__restrict__ x, T *__restrict__ p, T **__restrict__ u, const T *alpha, const T *beta, const int p_cur, const int p_space, const int n)
Definition: pipecg_kernel.h:43