Neko  0.8.99
A portable framework for high-order spectral element flow simulations
fusedcg_kernel.h
Go to the documentation of this file.
1 #ifndef __KRYLOV_FUSEDCG_KERNEL_H__
2 #define __KRYLOV_FUSEDCG_KERNEL_H__
3 /*
4  Copyright (c) 2021-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 
38 
42 template< typename T >
43 __global__ void fusedcg_update_p_kernel(T * __restrict__ p,
44  const T * __restrict__ z,
45  const T * __restrict__ po,
46  const T beta,
47  const int n) {
48 
49  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
50  const int str = blockDim.x * gridDim.x;
51 
52  for (int i = idx; i < n; i+= str) {
53  p[i] = beta*po[i] + z[i];
54  }
55 
56 }
57 
61 template< typename T >
62 __global__ void fusedcg_update_x_kernel(T * __restrict__ x,
63  const T ** p,
64  const T * __restrict__ alpha,
65  const int p_cur,
66  const int n) {
67 
68  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
69  const int str = blockDim.x * gridDim.x;
70 
71  for (int i = idx; i < n; i+= str) {
72  T tmp = 0.0;
73  for (int j = 0; j < p_cur; j ++) {
74  tmp += p[j][i] * alpha[j];
75  }
76  x[i] += tmp;
77  }
78 
79 }
80 
84 template< typename T>
85 __global__ void fusedcg_part2_kernel(T * __restrict__ a,
86  const T * __restrict__ b,
87  const T * __restrict__ c,
88  const T alpha,
89  T * buf_h,
90  const int n) {
91 
92  const int idx = blockIdx.x * blockDim.x + threadIdx.x;
93  const int str = blockDim.x * gridDim.x;
94 
95  const unsigned int lane = threadIdx.x % warpSize;
96  const unsigned int wid = threadIdx.x / warpSize;
97 
98  __shared__ T buf[64];
99  T tmp = 0.0;
100 
101  for (int i = idx; i < n; i+= str) {
102  T rt = a[i] - alpha * c[i];
103  tmp = tmp + rt * b[i] * rt;
104  a[i] = rt;
105  }
106 
107  tmp = reduce_warp<T>(tmp);
108  if (lane == 0) {
109  buf[wid] = tmp;
110  }
111  __syncthreads();
112 
113  tmp = (threadIdx.x < blockDim.x / warpSize) ? buf[lane] : 0;
114  if (wid == 0) {
115  tmp = reduce_warp<T>(tmp);
116  }
117 
118  if (threadIdx.x == 0) {
119  buf_h[blockIdx.x] = tmp;
120  }
121 }
122 
123 
124 #endif // __KRYLOV_FUSEDCG_KERNEL_H__
const int i
const int j
__syncthreads()
__global__ void const T *__restrict__ x
Definition: cdtp_kernel.h:106
__global__ void fusedcg_update_p_kernel(T *__restrict__ p, const T *__restrict__ z, const T *__restrict__ po, const T beta, const int n)
__global__ void fusedcg_part2_kernel(T *__restrict__ a, const T *__restrict__ b, const T *__restrict__ c, const T alpha, T *buf_h, const int n)
__global__ void fusedcg_update_x_kernel(T *__restrict__ x, const T **p, const T *__restrict__ alpha, const int p_cur, const int n)
real * buf
Definition: pipecg_aux.cu:42