Neko  0.8.99
A portable framework for high-order spectral element flow simulations
opr_cfl.hip
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 <stdio.h>
36 #include <hip/hip_runtime.h>
37 #include <device/device_config.h>
38 #include <device/hip/check.h>
39 #include "cfl_kernel.h"
40 
41 extern "C" {
42 
45 
49  real *cfl_d = NULL;
50 
54  real hip_cfl(real *dt, void *u, void *v, void *w,
55  void *drdx, void *dsdx, void *dtdx,
56  void *drdy, void *dsdy, void *dtdy,
57  void *drdz, void *dsdz, void *dtdz,
58  void *dr_inv, void *ds_inv, void *dt_inv,
59  void *jacinv, int *nel, int *lx) {
60 
61  const dim3 nthrds(1024, 1, 1);
62  const dim3 nblcks((*nel), 1, 1);
63 
64  if (cfl_d == NULL) {
65  HIP_CHECK(hipMalloc(&cfl_d, (*nel) * sizeof(real)));
66  }
67 
68 #define CASE(LX) \
69  case LX: \
70  hipLaunchKernelGGL(HIP_KERNEL_NAME( cfl_kernel<real, LX, 1024> ), \
71  nblcks, nthrds, 0, (hipStream_t) glb_cmd_queue, \
72  *dt, (real *) u, (real *) v, (real *) w, \
73  (real *) drdx, (real *) dsdx, (real *) dtdx, \
74  (real *) drdy, (real *) dsdy, (real *) dtdy, \
75  (real *) drdz, (real *) dsdz, (real *) dtdz, \
76  (real *) dr_inv, (real *) ds_inv, (real *) dt_inv, \
77  (real *) jacinv, (real *) cfl_d); \
78  HIP_CHECK(hipGetLastError()); \
79  break
80 
81  switch(*lx) {
82  CASE(2);
83  CASE(3);
84  CASE(4);
85  CASE(5);
86  CASE(6);
87  CASE(7);
88  CASE(8);
89  CASE(9);
90  CASE(10);
91  CASE(11);
92  CASE(12);
93  default:
94  {
95  fprintf(stderr, __FILE__ ": size not supported: %d\n", *lx);
96  exit(1);
97  }
98  }
99 
100 
101  hipLaunchKernelGGL(HIP_KERNEL_NAME(cfl_reduce_kernel<real>),
102  1, 1024, 0, (hipStream_t) glb_cmd_queue,
103  cfl_d, (*nel));
104  HIP_CHECK(hipGetLastError());
105 
106  real cfl = 0.0;
107 #ifdef HAVE_DEVICE_MPI
108  HIP_CHECK(hipStreamSynchronize((hipStream_t) glb_cmd_queue));
110 #else
111  HIP_CHECK(hipMemcpyAsync(&cfl, cfl_d, sizeof(real),
112  hipMemcpyDeviceToHost,
113  (hipStream_t) glb_cmd_queue));
114  HIP_CHECK(hipStreamSynchronize((hipStream_t) glb_cmd_queue));
115 #endif
116 
117 
118  return cfl;
119  }
120 }
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dt
Definition: cdtp_kernel.h:109
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ jacinv
Definition: conv1_kernel.h:148
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ drdx
Definition: conv1_kernel.h:139
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dtdx
Definition: conv1_kernel.h:141
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dtdz
Definition: conv1_kernel.h:147
__global__ void const T *__restrict__ u
Definition: conv1_kernel.h:132
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dsdx
Definition: conv1_kernel.h:140
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dtdy
Definition: conv1_kernel.h:144
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dsdz
Definition: conv1_kernel.h:146
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dsdy
Definition: conv1_kernel.h:143
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ drdy
Definition: conv1_kernel.h:142
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ drdz
Definition: conv1_kernel.h:145
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ w
__global__ void const T *__restrict__ const T *__restrict__ v
double real
Definition: device_config.h:12
void * glb_cmd_queue
#define DEVICE_MPI_MAX
Definition: device_mpi_op.h:10
void device_mpi_allreduce(void *buf_d, void *buf, int count, int nbytes, int op)
#define HIP_CHECK(err)
Definition: check.h:8
real(kind=rp) function, public cfl(dt, u, v, w, Xh, coef, nelv, gdim)
Definition: operators.f90:392
real * cfl_d
Definition: opr_cfl.hip:49
#define CASE(LX)
real hip_cfl(real *dt, void *u, void *v, void *w, void *drdx, void *dsdx, void *dtdx, void *drdy, void *dsdy, void *dtdy, void *drdz, void *dsdz, void *dtdz, void *dr_inv, void *ds_inv, void *dt_inv, void *jacinv, int *nel, int *lx)
Definition: opr_cfl.hip:54