Neko  0.9.99
A portable framework for high-order spectral element flow simulations
projection.cu
Go to the documentation of this file.
1 /*
2  Copyright (c) 2022, The Neko Authors
3  All rights reserved.
4 
5  Redistribution and use in source and binary forms, with or without
6  modification, are permitted provided that the following conditions
7  are met:
8 
9  * Redistributions of source code must retain the above copyright
10  notice, this list of conditions and the following disclaimer.
11 
12  * Redistributions in binary form must reproduce the above
13  copyright notice, this list of conditions and the following
14  disclaimer in the documentation and/or other materials provided
15  with the distribution.
16 
17  * Neither the name of the authors nor the names of its
18  contributors may be used to endorse or promote products derived
19  from this software without specific prior written permission.
20 
21  THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
22  "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
23  LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
24  FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
25  COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
26  INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
27  BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28  LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
29  CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
30  LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
31  ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
32  POSSIBILITY OF SUCH DAMAGE.
33 */
34 
35 #include <device/device_config.h>
36 #include <device/cuda/check.h>
37 
38 #include "projection_kernel.h"
40 
41 
42 /*
43  * Reduction buffer
44  */
45 int proj_red_s = 0;
47 
48 extern "C" {
49 
52 
53  void cuda_project_on(void *alpha, void * b, void *xx, void *bb, void *mult,
54  void *xbar, int *j, int *n){
55 
56  const cudaStream_t stream = (cudaStream_t) glb_cmd_queue;
57  int pow2 = 1;
58  while(pow2 < (*j)){
59  pow2 = 2*pow2;
60  }
61  const int nt = 1024/pow2;
62  const dim3 glsc3_nthrds(pow2, nt, 1);
63  const dim3 glsc3_nblcks(((*n)+nt - 1)/nt, 1, 1);
64  const int glsc3_nb = ((*n) + nt - 1)/nt;
65  if((*j)*glsc3_nb>proj_red_s){
66  proj_red_s = (*j)*glsc3_nb;
67  if (proj_bufred_d != NULL) {
68  CUDA_CHECK(cudaFree(proj_bufred_d));
69  }
70  CUDA_CHECK(cudaMalloc(&proj_bufred_d, (*j)*glsc3_nb*sizeof(real)));
71  }
72 
73  /* First glsc3_many call */
74  glsc3_many_kernel<real>
75  <<<glsc3_nblcks, glsc3_nthrds, 0, stream>>>((const real *) b,
76  (const real **) xx,
77  (const real *) mult,
78  proj_bufred_d, *j, *n);
79  CUDA_CHECK(cudaGetLastError());
80  glsc3_reduce_kernel<real>
81  <<<(*j), 1024, 0, stream>>>(proj_bufred_d, glsc3_nb, *j);
82  CUDA_CHECK(cudaGetLastError());
83  CUDA_CHECK(cudaMemcpyAsync(alpha, proj_bufred_d, (*j) * sizeof(real),
84  cudaMemcpyDeviceToDevice, stream));
85  CUDA_CHECK(cudaMemsetAsync(xbar, 0, (*n) * sizeof(real), stream));
86 
87  cudaStreamSynchronize(stream);
89 
90  const dim3 vec_nthrds(1024, 1, 1);
91  const dim3 vec_nblcks(((*n)+1024 - 1)/ 1024, 1, 1);
92 
93  /* First vector operation block */
94  project_on_vec_kernel<real>
95  <<<vec_nblcks, vec_nthrds, 0, stream>>>((real *) xbar,
96  (const real **) xx,
97  (real *) b,
98  (const real **) bb,
99  (const real *) alpha,
100  *j, *n);
101  /* Second glsc3_many call */
102  glsc3_many_kernel<real>
103  <<<glsc3_nblcks, glsc3_nthrds, 0, stream>>>((const real *) b,
104  (const real **) xx,
105  (const real *) mult,
106  proj_bufred_d, *j, *n);
107  CUDA_CHECK(cudaGetLastError());
108  glsc3_reduce_kernel<real>
109  <<<(*j), 1024, 0, stream>>>(proj_bufred_d, glsc3_nb, *j);
110  CUDA_CHECK(cudaGetLastError());
111  CUDA_CHECK(cudaMemcpyAsync(alpha, proj_bufred_d, (*j) * sizeof(real),
112  cudaMemcpyDeviceToDevice, stream));
113 
114  cudaStreamSynchronize(stream);
115  device_mpi_allreduce_inplace(alpha, (*j), sizeof(real), DEVICE_MPI_SUM);
116 
117  /* Second vector operation block */
118  project_on_vec_kernel<real>
119  <<<vec_nblcks, vec_nthrds, 0, stream>>>((real *) xbar,
120  (const real **) xx,
121  (real *) b,
122  (const real **) bb,
123  (const real *) alpha,
124  *j, *n);
125  }
126 
127  void cuda_project_ortho(void *alpha, void * b, void *xx, void *bb,
128  void *w, void *xm, int *j, int *n, real *nrm){
129 
130  const cudaStream_t stream = (cudaStream_t) glb_cmd_queue;
131 
132  int pow2 = 1;
133  while(pow2 < (*j)){
134  pow2 = 2*pow2;
135  }
136  const int nt = 1024/pow2;
137  const dim3 glsc3_nthrds(pow2, nt, 1);
138  const dim3 glsc3_nblcks(((*n)+nt - 1)/nt, 1, 1);
139  const int glsc3_nb = ((*n) + nt - 1)/nt;
140  if((*j)*glsc3_nb>proj_red_s){
141  proj_red_s = (*j)*glsc3_nb;
142  if (proj_bufred_d != NULL) {
143  CUDA_CHECK(cudaFree(proj_bufred_d));
144  }
145  CUDA_CHECK(cudaMalloc(&proj_bufred_d, (*j)*glsc3_nb*sizeof(real)));
146  }
147 
148  /* First glsc3_many call */
149  glsc3_many_kernel<real>
150  <<<glsc3_nblcks, glsc3_nthrds, 0, stream>>>((const real *) b,
151  (const real **) xx,
152  (const real *) w,
153  proj_bufred_d, *j, *n);
154  CUDA_CHECK(cudaGetLastError());
155  glsc3_reduce_kernel<real>
156  <<<(*j), 1024, 0, stream>>>(proj_bufred_d, glsc3_nb, *j);
157  CUDA_CHECK(cudaGetLastError());
158  CUDA_CHECK(cudaMemcpyAsync(alpha, proj_bufred_d, (*j) * sizeof(real),
159  cudaMemcpyDeviceToDevice, stream));
160 
161  cudaStreamSynchronize(stream);
162  device_mpi_allreduce_inplace(alpha, (*j), sizeof(real), DEVICE_MPI_SUM);
163 
164  CUDA_CHECK(cudaMemcpyAsync(nrm, (real *) alpha + (*j - 1),
165  sizeof(real), cudaMemcpyDeviceToHost, stream));
166  (*nrm) = sqrt(*nrm);
167 
168 
169  const dim3 vec_nthrds(1024, 1, 1);
170  const dim3 vec_nblcks(((*n)+1024 - 1)/ 1024, 1, 1);
171 
172  /* First vector operation block */
173  project_ortho_vec_kernel<real>
174  <<<vec_nblcks, vec_nthrds, 0, stream>>>((real *) xm, (const real **) xx,
175  (real *) b, (const real **) bb,
176  (const real *) alpha, *j, *n);
177 
178  /* Second glsc3_many call */
179  glsc3_many_kernel<real>
180  <<<glsc3_nblcks, glsc3_nthrds, 0, stream>>>((const real *) b,
181  (const real **) xx,
182  (const real *) w,
183  proj_bufred_d, *j, *n);
184  CUDA_CHECK(cudaGetLastError());
185  glsc3_reduce_kernel<real>
186  <<<(*j), 1024, 0, stream>>> (proj_bufred_d, glsc3_nb, *j);
187  CUDA_CHECK(cudaGetLastError());
188  CUDA_CHECK(cudaMemcpyAsync(alpha, proj_bufred_d, (*j) * sizeof(real),
189  cudaMemcpyDeviceToDevice, stream));
190 
191  cudaStreamSynchronize(stream);
192  device_mpi_allreduce_inplace(alpha, (*j), sizeof(real), DEVICE_MPI_SUM);
193 
194  /* Second vector operation block */
195  project_ortho_vec_kernel<real>
196  <<<vec_nblcks, vec_nthrds, 0, stream>>>((real *) xm, (const real **) xx,
197  (real *) b, (const real **) bb,
198  (const real *) alpha, *j, *n);
199 
200  }
201 
202 }
203 
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ w
const int j
#define CUDA_CHECK(err)
Definition: check.h:6
double real
Definition: device_config.h:12
void * glb_cmd_queue
#define DEVICE_MPI_SUM
Definition: device_mpi_op.h:9
void device_mpi_allreduce_inplace(void *buf_d, int count, int nbytes, int op)
void cuda_project_ortho(void *alpha, void *b, void *xx, void *bb, void *w, void *xm, int *j, int *n, real *nrm)
Definition: projection.cu:127
int proj_red_s
Definition: projection.cu:45
real * proj_bufred_d
Definition: projection.cu:46
void cuda_project_on(void *alpha, void *b, void *xx, void *bb, void *mult, void *xbar, int *j, int *n)
Definition: projection.cu:53