Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
opr_cfl.c
Go to the documentation of this file.
1/*
2 Copyright (c) 2022-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#ifdef __APPLE__
36#include <OpenCL/cl.h>
37#else
38#include <CL/cl.h>
39#endif
40
41#include <stdio.h>
42#include <stdlib.h>
43#include <math.h>
45#include <device/opencl/jit.h>
47#include <device/opencl/check.h>
48
49#include "cfl_kernel.cl.h"
50
54real opencl_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 cl_int err;
62 int i;
63 if (cfl_program == NULL)
65
66 const size_t global_item_size = 256 * (*nel);
67 const size_t local_item_size = 256;
68
69 real * cfl = (real *) malloc((*nel) * sizeof(real));
71 (*nel) * sizeof(real), NULL, &err);
73
74#define STR(X) #X
75#define CASE(LX) \
76 case LX: \
77 { \
78 cl_kernel kernel = clCreateKernel(cfl_program, \
79 STR(cfl_kernel_lx##LX), &err); \
80 CL_CHECK(err); \
81 \
82 CL_CHECK(clSetKernelArg(kernel, 0, sizeof(real), dt)); \
83 CL_CHECK(clSetKernelArg(kernel, 1, sizeof(cl_mem), (void *) &u)); \
84 CL_CHECK(clSetKernelArg(kernel, 2, sizeof(cl_mem), (void *) &v)); \
85 CL_CHECK(clSetKernelArg(kernel, 3, sizeof(cl_mem), (void *) &w)); \
86 CL_CHECK(clSetKernelArg(kernel, 4, sizeof(cl_mem), (void *) &drdx)); \
87 CL_CHECK(clSetKernelArg(kernel, 5, sizeof(cl_mem), (void *) &dsdx)); \
88 CL_CHECK(clSetKernelArg(kernel, 6, sizeof(cl_mem), (void *) &dtdx)); \
89 CL_CHECK(clSetKernelArg(kernel, 7, sizeof(cl_mem), (void *) &drdy)); \
90 CL_CHECK(clSetKernelArg(kernel, 8, sizeof(cl_mem), (void *) &dsdy)); \
91 CL_CHECK(clSetKernelArg(kernel, 9, sizeof(cl_mem), (void *) &dtdy)); \
92 CL_CHECK(clSetKernelArg(kernel, 10, sizeof(cl_mem), (void *) &drdz)); \
93 CL_CHECK(clSetKernelArg(kernel, 11, sizeof(cl_mem), (void *) &dsdz)); \
94 CL_CHECK(clSetKernelArg(kernel, 12, sizeof(cl_mem), (void *) &dtdz)); \
95 CL_CHECK(clSetKernelArg(kernel, 13, sizeof(cl_mem), (void *) &dr_inv)); \
96 CL_CHECK(clSetKernelArg(kernel, 14, sizeof(cl_mem), (void *) &ds_inv)); \
97 CL_CHECK(clSetKernelArg(kernel, 15, sizeof(cl_mem), (void *) &dt_inv)); \
98 CL_CHECK(clSetKernelArg(kernel, 16, sizeof(cl_mem), (void *) &jacinv)); \
99 CL_CHECK(clSetKernelArg(kernel, 17, sizeof(cl_mem), (void *) &cfl_d)); \
100 \
101 CL_CHECK(clEnqueueNDRangeKernel((cl_command_queue) glb_cmd_queue, \
102 kernel, 1, NULL, &global_item_size, \
103 &local_item_size, 0, NULL, &kern_wait)); \
104 } \
105 break
106
107 switch(*lx) {
108 CASE(2);
109 CASE(3);
110 CASE(4);
111 CASE(5);
112 CASE(6);
113 CASE(7);
114 CASE(8);
115 CASE(9);
116 CASE(10);
117 CASE(11);
118 CASE(12);
119 }
120
122 cfl_d, CL_TRUE, 0, (*nel) * sizeof(real),
123 cfl, 1, &kern_wait, NULL));
124
125 real cfl_max = 0.0;
126 for (i = 0; i < (*nel); i++) {
127 cfl_max = fmax(cfl_max, cfl[i]);
128 }
129
130 free(cfl);
132
133 return cfl_max;
134
135}
__global__ void T *__restrict__ 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
__global__ void T *__restrict__ 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
__global__ void T *__restrict__ 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
__global__ void T *__restrict__ 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__ dsdy
__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__ 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
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ u
__global__ void T *__restrict__ 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__ drdx
__global__ void T *__restrict__ 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__ dtdz
__global__ void T *__restrict__ 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__ dsdx
__global__ void T *__restrict__ T *__restrict__ const T *__restrict__ const T *__restrict__ v
__global__ void T *__restrict__ 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__ dtdx
__global__ void T *__restrict__ 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__ const T *__restrict__ jacinv
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dt
__global__ void cfl_kernel(const T dt, const T *__restrict__ u, const T *__restrict__ v, const T *__restrict__ w, const T *__restrict__ drdx, const T *__restrict__ dsdx, const T *__restrict__ dtdx, const T *__restrict__ drdy, const T *__restrict__ dsdy, const T *__restrict__ dtdy, const T *__restrict__ drdz, const T *__restrict__ dsdz, const T *__restrict__ dtdz, const T *__restrict__ dr_inv, const T *__restrict__ ds_inv, const T *__restrict__ dt_inv, const T *__restrict__ jacinv, T *__restrict__ cfl_h)
Definition cfl_kernel.h:87
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
double real
void * glb_ctx
void opencl_kernel_jit(const char *kernel, cl_program *program)
Definition jit.c:50
#define CL_CHECK(err)
Definition check.h:12
#define CASE(LX)
real opencl_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.c:54
real * cfl_d
Definition opr_cfl.cu:49
void * cfl_program