Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fusedcg_cpld_kernel.h
Go to the documentation of this file.
1#ifndef __KRYLOV_FUSEDCG_CPLD_KERNEL_H__
2#define __KRYLOV_FUSEDCG_CPLD_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
42template< typename T >
49 T * __restrict__ tmp,
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 for (int i = idx; i < n; i+= str) {
56 tmp[i] = a1[i]*b1[i] + a2[i]*b2[i] + a3[i]*b3[i];
57 }
58
59}
60
64template< typename T >
66 T * __restrict__ p2,
67 T * __restrict__ p3,
68 const T * __restrict__ z1,
69 const T * __restrict__ z2,
70 const T * __restrict__ z3,
71 const T * __restrict__ po1,
72 const T * __restrict__ po2,
73 const T * __restrict__ po3,
74 const T beta,
75 const int n) {
76
77 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
78 const int str = blockDim.x * gridDim.x;
79
80 for (int i = idx; i < n; i+= str) {
81 p1[i] = beta*po1[i] + z1[i];
82 p2[i] = beta*po2[i] + z2[i];
83 p3[i] = beta*po3[i] + z3[i];
84 }
85
86}
87
91template< typename T >
95 const T ** p1,
96 const T ** p2,
97 const T ** p3,
98 const T * __restrict__ alpha,
99 const int p_cur,
100 const int n) {
101
102 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
103 const int str = blockDim.x * gridDim.x;
104
105 for (int i = idx; i < n; i+= str) {
106 T tmp1 = 0.0;
107 T tmp2 = 0.0;
108 T tmp3 = 0.0;
109 for (int j = 0; j < p_cur; j ++) {
110 tmp1 += p1[j][i] * alpha[j];
111 tmp2 += p2[j][i] * alpha[j];
112 tmp3 += p3[j][i] * alpha[j];
113 }
114 x1[i] += tmp1;
115 x2[i] += tmp2;
116 x3[i] += tmp3;
117 }
118
119}
120
124template< typename T>
126 T * __restrict__ a2,
127 T * __restrict__ a3,
128 const T * __restrict__ b,
129 const T * __restrict__ c1,
130 const T * __restrict__ c2,
131 const T * __restrict__ c3,
132 const T alpha,
133 T * buf_h,
134 const int n) {
135
136 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
137 const int str = blockDim.x * gridDim.x;
138
139 const unsigned int lane = threadIdx.x % warpSize;
140 const unsigned int wid = threadIdx.x / warpSize;
141
142 __shared__ T buf[32];
143 T tmp = 0.0;
144
145 for (int i = idx; i < n; i+= str) {
146 T rt1 = a1[i] - alpha * c1[i];
147 T rt2 = a2[i] - alpha * c2[i];
148 T rt3 = a3[i] - alpha * c3[i];
149 tmp = tmp + ((rt1*rt1 + rt2*rt2 + rt3*rt3) * b[i]);
150 a1[i] = rt1;
151 a2[i] = rt2;
152 a3[i] = rt3;
153 }
154
155 tmp = reduce_warp<T>(tmp);
156 if (lane == 0) {
157 buf[wid] = tmp;
158 }
160
161 tmp = (threadIdx.x < blockDim.x / warpSize) ? buf[lane] : 0;
162 if (wid == 0) {
163 tmp = reduce_warp<T>(tmp);
164 }
165
166 if (threadIdx.x == 0) {
167 buf_h[blockIdx.x] = tmp;
168 }
169}
170
171#endif // __KRYLOV_FUSEDCG_CPLD_KERNEL_H__
const int i
const int j
__syncthreads()
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
__global__ void fusedcg_cpld_update_p_kernel(T *__restrict__ p1, T *__restrict__ p2, T *__restrict__ p3, const T *__restrict__ z1, const T *__restrict__ z2, const T *__restrict__ z3, const T *__restrict__ po1, const T *__restrict__ po2, const T *__restrict__ po3, const T beta, const int n)
__global__ void fusedcg_cpld_part2_kernel(T *__restrict__ a1, T *__restrict__ a2, T *__restrict__ a3, const T *__restrict__ b, const T *__restrict__ c1, const T *__restrict__ c2, const T *__restrict__ c3, const T alpha, T *buf_h, const int n)
__global__ void fusedcg_cpld_update_x_kernel(T *__restrict__ x1, T *__restrict__ x2, T *__restrict__ x3, const T **p1, const T **p2, const T **p3, const T *__restrict__ alpha, const int p_cur, const int n)
__global__ void fusedcg_cpld_part1_kernel(T *__restrict__ a1, T *__restrict__ a2, T *__restrict__ a3, T *__restrict__ b1, T *__restrict__ b2, T *__restrict__ b3, T *__restrict__ tmp, const int n)
real * buf
Definition pipecg_aux.cu:42