Neko  0.8.99
A portable framework for high-order spectral element flow simulations
gmres_kernel.h
Go to the documentation of this file.
1 #ifndef __KRYLOV_GMRES_KERNEL_H__
2 #define __KRYLOV_GMRES_KERNEL_H__
3 /*
4  Copyright (c) 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 
42 
43 template< typename T >
44 __global__ void gmres_part2_kernel(T * __restrict__ w,
45  T * const * __restrict__ v,
46  const T * __restrict__ mult,
47  const T * __restrict__ h,
48  T * __restrict__ buf_h1,
49  const int j,
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  const unsigned int lane = threadIdx.x % warpSize;
56  const unsigned int wid = threadIdx.x / warpSize;
57 
58  __shared__ T shared[64];
59  T tmp1 = 0.0;
60 
61  for (int i = idx; i < n; i+= str) {
62  T tmp = 0.0;
63  for (int k = 0; k < j; k ++) {
64  tmp += -h[k]*v[k][i];
65  }
66  w[i] += tmp;
67  tmp1 += w[i]*w[i]*mult[i];
68  }
69 
70  tmp1 = reduce_warp<T>(tmp1);
71  if (lane == 0)
72  shared[wid] = tmp1;
73  __syncthreads();
74 
75  tmp1 = (threadIdx.x < blockDim.x / warpSize) ? shared[lane] : 0;
76  if (wid == 0)
77  tmp1 = reduce_warp<T>(tmp1);
78 
79  if (threadIdx.x == 0)
80  buf_h1[blockIdx.x] = tmp1;
81 
82 }
83 
84 
85 
86 #endif // __KRYLOV_GMRES_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__ const T *__restrict__ v
const int j
__syncthreads()
__global__ void gmres_part2_kernel(T *__restrict__ w, T *const *__restrict__ v, const T *__restrict__ mult, const T *__restrict__ h, T *__restrict__ buf_h1, const int j, const int n)
Definition: gmres_kernel.h:43