Neko  0.8.99
A portable framework for high-order spectral element flow simulations
fusedcg_cpld_aux.cu
Go to the documentation of this file.
1 /*
2  Copyright (c) 2024, 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 "fusedcg_cpld_kernel.h"
36 #include <device/device_config.h>
37 #include <device/cuda/check.h>
38 
45 
46 extern "C" {
47 
50 
51  void cuda_fusedcg_cpld_part1(void *a1, void *a2, void *a3,
52  void *b1, void *b2, void *b3,
53  void *tmp, int *n) {
54 
55  const dim3 nthrds(1024, 1, 1);
56  const dim3 nblcks(((*n)+1024 - 1)/ 1024, 1, 1);
57  const cudaStream_t stream = (cudaStream_t) glb_cmd_queue;
58 
59  fusedcg_cpld_part1_kernel<real>
60  <<<nblcks, nthrds, 0, stream>>>((real *) a1, (real *) a2, (real *) a3,
61  (real *) b1, (real *) b2, (real *) b3,
62  (real *) tmp, *n);
63  CUDA_CHECK(cudaGetLastError());
64  }
65 
66  void cuda_fusedcg_cpld_update_p(void *p1, void *p2, void *p3,
67  void *z1, void *z2, void *z3,
68  void *po1, void *po2, void *po3,
69  real *beta, int *n) {
70 
71  const dim3 nthrds(1024, 1, 1);
72  const dim3 nblcks(((*n)+1024 - 1)/ 1024, 1, 1);
73  const cudaStream_t stream = (cudaStream_t) glb_cmd_queue;
74 
75  fusedcg_cpld_update_p_kernel<real>
76  <<<nblcks, nthrds, 0, stream>>>((real *) p1, (real *) p2, (real *) p3,
77  (real *) z1, (real *) z2, (real *) z3,
78  (real *) po1, (real *) po2, (real *) po3,
79  *beta, *n);
80  CUDA_CHECK(cudaGetLastError());
81 
82  }
83 
84  void cuda_fusedcg_cpld_update_x(void *x1, void *x2, void *x3,
85  void *p1, void *p2, void *p3,
86  void *alpha, int *p_cur, int *n) {
87 
88  const dim3 nthrds(1024, 1, 1);
89  const dim3 nblcks(((*n)+1024 - 1)/ 1024, 1, 1);
90  const cudaStream_t stream = (cudaStream_t) glb_cmd_queue;
91 
92  fusedcg_cpld_update_x_kernel<real>
93  <<<nblcks, nthrds, 0, stream>>>((real *) x1, (real *) x2, (real *) x3,
94  (const real **) p1, (const real **) p2,
95  (const real **) p3, (const real *) alpha,
96  *p_cur, *n);
97  CUDA_CHECK(cudaGetLastError());
98  }
99 
100  real cuda_fusedcg_cpld_part2(void *a1, void *a2, void *a3, void *b,
101  void *c1, void *c2, void *c3, void *alpha_d ,
102  real *alpha, int *p_cur, int * n) {
103 
104  const dim3 nthrds(1024, 1, 1);
105  const dim3 nblcks(((*n)+1024 - 1)/ 1024, 1, 1);
106  const int nb = ((*n) + 1024 - 1)/ 1024;
107  const cudaStream_t stream = (cudaStream_t) glb_cmd_queue;
108 
109  if (fusedcg_cpld_buf != NULL && fusedcg_cpld_buf_len < nb) {
110  CUDA_CHECK(cudaFreeHost(fusedcg_cpld_buf));
111  CUDA_CHECK(cudaFree(fusedcg_cpld_buf_d));
112  fusedcg_cpld_buf = NULL;
113  }
114 
115  if (fusedcg_cpld_buf == NULL){
116  CUDA_CHECK(cudaMallocHost(&fusedcg_cpld_buf, 2*sizeof(real)));
117  CUDA_CHECK(cudaMalloc(&fusedcg_cpld_buf_d, nb*sizeof(real)));
119  }
120 
121  /* Store alpha(p_cur) in pinned memory */
122  fusedcg_cpld_buf[1] = (*alpha);
123 
124  /* Update alpha_d(p_cur) = alpha(p_cur) */
125  real *alpha_d_p_cur = ((real *) alpha_d) + ((*p_cur - 1));
126  CUDA_CHECK(cudaMemcpyAsync(alpha_d_p_cur, &fusedcg_cpld_buf[1],
127  sizeof(real), cudaMemcpyHostToDevice,
128  stream));
129 
130 
131  fusedcg_cpld_part2_kernel<real>
132  <<<nblcks, nthrds, 0, stream>>>((real *) a1, (real *) a2, (real *) a3,
133  (real *) b, (real *) c1, (real *) c2,
134  (real *) c3, *alpha, fusedcg_cpld_buf_d,
135  *n);
136  CUDA_CHECK(cudaGetLastError());
137 
138  reduce_kernel<real><<<1, 1024, 0, stream>>>(fusedcg_cpld_buf_d, nb);
139  CUDA_CHECK(cudaGetLastError());
140 
141 #ifdef HAVE_DEVICE_MPI
142  cudaStreamSynchronize(stream);
144  sizeof(real), DEVICE_MPI_SUM);
145 #else
146  CUDA_CHECK(cudaMemcpyAsync(fusedcg_cpld_buf, fusedcg_cpld_buf_d, sizeof(real),
147  cudaMemcpyDeviceToHost, stream));
148  cudaStreamSynchronize(stream);
149 #endif
150 
151  return fusedcg_cpld_buf[0];
152  }
153 }
154 
#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(void *buf_d, void *buf, int count, int nbytes, int op)
void cuda_fusedcg_cpld_update_x(void *x1, void *x2, void *x3, void *p1, void *p2, void *p3, void *alpha, int *p_cur, int *n)
int fusedcg_cpld_buf_len
void cuda_fusedcg_cpld_part1(void *a1, void *a2, void *a3, void *b1, void *b2, void *b3, void *tmp, int *n)
real cuda_fusedcg_cpld_part2(void *a1, void *a2, void *a3, void *b, void *c1, void *c2, void *c3, void *alpha_d, real *alpha, int *p_cur, int *n)
real * fusedcg_cpld_buf
real * fusedcg_cpld_buf_d
void cuda_fusedcg_cpld_update_p(void *p1, void *p2, void *p3, void *z1, void *z2, void *z3, void *po1, void *po2, void *po3, real *beta, int *n)