Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
entropy_viscosity.c
Go to the documentation of this file.
1/*
2 Copyright (c) 2025, 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>
43#include <device/opencl/jit.h>
45#include <device/opencl/check.h>
46
47#include "entropy_viscosity_kernel.cl.h"
48
49void opencl_entropy_visc_compute_residual(void *entropy_residual,
50 void *S, void *S_lag1,
51 void *S_lag2, void *S_lag3,
54 real dt, int n) {
55 cl_int err;
56
60
62 "entropy_visc_compute_residual_kernel", &err);
64
65 CL_CHECK(clSetKernelArg(kernel, 0, sizeof(cl_mem), (void *) &entropy_residual));
66 CL_CHECK(clSetKernelArg(kernel, 1, sizeof(cl_mem), (void *) &S));
67 CL_CHECK(clSetKernelArg(kernel, 2, sizeof(cl_mem), (void *) &S_lag1));
68 CL_CHECK(clSetKernelArg(kernel, 3, sizeof(cl_mem), (void *) &S_lag2));
69 CL_CHECK(clSetKernelArg(kernel, 4, sizeof(cl_mem), (void *) &S_lag3));
70 CL_CHECK(clSetKernelArg(kernel, 5, sizeof(real), &bdf1));
71 CL_CHECK(clSetKernelArg(kernel, 6, sizeof(real), &bdf2));
72 CL_CHECK(clSetKernelArg(kernel, 7, sizeof(real), &bdf3));
73 CL_CHECK(clSetKernelArg(kernel, 8, sizeof(real), &bdf4));
74 CL_CHECK(clSetKernelArg(kernel, 9, sizeof(real), &dt));
75 CL_CHECK(clSetKernelArg(kernel, 10, sizeof(int), &n));
76
77 const int nb = (n + 256 - 1) / 256;
78 const size_t global_item_size = 256 * nb;
79 const size_t local_item_size = 256;
80
83 0, NULL, NULL));
85}
86
88 void *entropy_residual,
89 void *h,
90 real c_avisc_entropy,
91 real n_S,
92 int n) {
93 cl_int err;
94
98
100 "entropy_visc_compute_viscosity_kernel", &err);
101 CL_CHECK(err);
102
103 CL_CHECK(clSetKernelArg(kernel, 0, sizeof(cl_mem), (void *) &reg_coeff));
104 CL_CHECK(clSetKernelArg(kernel, 1, sizeof(cl_mem), (void *) &entropy_residual));
105 CL_CHECK(clSetKernelArg(kernel, 2, sizeof(cl_mem), (void *) &h));
106 CL_CHECK(clSetKernelArg(kernel, 3, sizeof(real), &c_avisc_entropy));
107 CL_CHECK(clSetKernelArg(kernel, 4, sizeof(real), &n_S));
108 CL_CHECK(clSetKernelArg(kernel, 5, sizeof(int), &n));
109
110 const int nb = (n + 256 - 1) / 256;
111 const size_t global_item_size = 256 * nb;
112 const size_t local_item_size = 256;
113
116 0, NULL, NULL));
118}
119
121 int lx,
122 int nelv) {
123 cl_int err;
124 const int lx3 = lx * lx * lx;
125
129
131 "entropy_visc_apply_element_max_kernel", &err);
132 CL_CHECK(err);
133
134 CL_CHECK(clSetKernelArg(kernel, 0, sizeof(cl_mem), (void *) &reg_coeff));
135 CL_CHECK(clSetKernelArg(kernel, 1, sizeof(int), &lx3));
136 CL_CHECK(clSetKernelArg(kernel, 2, sizeof(int), &nelv));
137
138 const int local_size = (lx3 < 256) ? lx3 : 256;
139 const size_t global_item_size = local_size * nelv;
140 const size_t local_item_size = local_size;
141
144 0, NULL, NULL));
146}
147
149 void *h,
150 void *max_wave_speed,
151 real c_avisc_low,
152 int n) {
153 cl_int err;
154
158
160 "entropy_visc_clamp_to_low_order_kernel", &err);
161 CL_CHECK(err);
162
163 CL_CHECK(clSetKernelArg(kernel, 0, sizeof(cl_mem), (void *) &reg_coeff));
164 CL_CHECK(clSetKernelArg(kernel, 1, sizeof(cl_mem), (void *) &h));
165 CL_CHECK(clSetKernelArg(kernel, 2, sizeof(cl_mem), (void *) &max_wave_speed));
166 CL_CHECK(clSetKernelArg(kernel, 3, sizeof(real), &c_avisc_low));
167 CL_CHECK(clSetKernelArg(kernel, 4, sizeof(int), &n));
168
169 const int nb = (n + 256 - 1) / 256;
170 const size_t global_item_size = 256 * nb;
171 const size_t local_item_size = 256;
172
175 0, NULL, NULL));
177}
178
180 void *temp_field,
181 void *mult_field,
182 int n) {
183 cl_int err;
184
188
190 "entropy_visc_smooth_divide_kernel", &err);
191 CL_CHECK(err);
192
193 CL_CHECK(clSetKernelArg(kernel, 0, sizeof(cl_mem), (void *) &reg_coeff));
194 CL_CHECK(clSetKernelArg(kernel, 1, sizeof(cl_mem), (void *) &temp_field));
195 CL_CHECK(clSetKernelArg(kernel, 2, sizeof(cl_mem), (void *) &mult_field));
196 CL_CHECK(clSetKernelArg(kernel, 3, sizeof(int), &n));
197
198 const int nb = (n + 256 - 1) / 256;
199 const size_t global_item_size = 256 * nb;
200 const size_t local_item_size = 256;
201
204 0, NULL, NULL));
206}
207
__global__ void const T *__restrict__ const T *__restrict__ const T *__restrict__ const T *__restrict__ dt
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
double real
void opencl_entropy_visc_compute_viscosity(void *reg_coeff, void *entropy_residual, void *h, real c_avisc_entropy, real n_S, int n)
void opencl_entropy_visc_compute_residual(void *entropy_residual, void *S, void *S_lag1, void *S_lag2, void *S_lag3, real bdf1, real bdf2, real bdf3, real bdf4, real dt, int n)
void opencl_entropy_visc_smooth_divide(void *reg_coeff, void *temp_field, void *mult_field, int n)
void opencl_entropy_visc_clamp_to_low_order(void *reg_coeff, void *h, void *max_wave_speed, real c_avisc_low, int n)
void opencl_entropy_visc_apply_element_max(void *reg_coeff, int lx, int nelv)
void opencl_kernel_jit(const char *kernel, cl_program *program)
Definition jit.c:50
#define CL_CHECK(err)
Definition check.h:12
void * entropy_viscosity_program