Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
RAMP_mapping_kernel.h
Go to the documentation of this file.
1/*
2 Copyright (c) 2021-2023, 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#ifndef __NEKO_CUDA_RAMP_MAPPING_KERNELS__
36#define __NEKO_CUDA_RAMP_MAPPING_KERNELS__
37
41template <typename T>
43 const T f_min, const T f_max, const T q, T* __restrict__ X_out_d,
44 T* __restrict__ X_in_d, const int n) {
45
46 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
47 const int str = blockDim.x * gridDim.x;
48
49 for (int i = idx; i < n; i += str) {
50 X_out_d[i] = f_min
51 + (f_max - f_min) * X_in_d[i] / (1.0 + q * (1.0 + X_in_d[i]));
52 }
53}
54
58template <typename T>
60 const T f_min, const T f_max, const T q, T* __restrict__ dF_dX_in_d,
62
63 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
64 const int str = blockDim.x * gridDim.x;
65
66 for (int i = idx; i < n; i += str) {
67 dF_dX_in_d[i] = (f_max - f_min) * (q + 1.0) / (
68 (1.0 - q * (X_in_d[i] - 1.0)) * (1.0 - q * (X_in_d[i] - 1.0))
69 ) * dF_dX_out_d[i];
70 }
71}
75template <typename T>
77 const T f_min, const T f_max, const T q, T* __restrict__ X_out_d,
78 T* __restrict__ X_in_d, const int n) {
79
80 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
81 const int str = blockDim.x * gridDim.x;
82
83 for (int i = idx; i < n; i += str) {
84 X_out_d[i] = f_min
85 + (f_max - f_min) * X_in_d[i] * (1.0 + q) / (X_in_d[i] + q);
86 }
87}
88
92template <typename T>
94 const T f_min, const T f_max, const T q, T* __restrict__ dF_dX_in_d,
96
97 const int idx = blockIdx.x * blockDim.x + threadIdx.x;
98 const int str = blockDim.x * gridDim.x;
99
100 for (int i = idx; i < n; i += str) {
101 dF_dX_in_d[i] = (f_max - f_min)
102 * (q + 1.0) / ( (X_in_d[i] + q) * (X_in_d[i] + q)) * dF_dX_out_d[i];
103 }
104}
105#endif // __NEKO_CUDA_RAMP_MAPPING_KERNELS__
__global__ void convex_down_RAMP_mapping_apply_backward_kernel(const T f_min, const T f_max, const T q, T *__restrict__ dF_dX_in_d, T *__restrict__ dF_dX_out_d, T *__restrict__ X_in_d, const int n)
__global__ void convex_down_RAMP_mapping_apply_kernel(const T f_min, const T f_max, const T q, T *__restrict__ X_out_d, T *__restrict__ X_in_d, const int n)
__global__ void convex_up_RAMP_mapping_apply_kernel(const T f_min, const T f_max, const T q, T *__restrict__ X_out_d, T *__restrict__ X_in_d, const int n)
__global__ void convex_up_RAMP_mapping_apply_backward_kernel(const T f_min, const T f_max, const T q, T *__restrict__ dF_dX_in_d, T *__restrict__ dF_dX_out_d, T *__restrict__ X_in_d, const int n)