Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
heaviside_mapping_kernel.h
Go to the documentation of this file.
1
37
#ifndef __NEKO_CUDA_HEAVISIDE_MAPPING_KERNELS__
38
#define __NEKO_CUDA_HEAVISIDE_MAPPING_KERNELS__
39
40
#include <math.h>
41
45
template
<
typename
T>
46
__global__
void
heaviside_mapping_apply_kernel
(
47
const
T
beta,
const
T
eta,
T
*
__restrict__
X_out_d
,
48
T
*
__restrict__
X_in_d
,
const
int
n) {
49
50
const
int
idx
=
blockIdx
.x *
blockDim
.x +
threadIdx
.x;
51
const
int
str
=
blockDim
.x *
gridDim
.x;
52
const
T
den
=
tanh
(beta * eta) +
tanh
(beta * (1.0 - eta));
53
const
T
tanh_beta_eta
=
tanh
(beta * eta);
54
55
for
(
int
i =
idx
; i < n; i +=
str
) {
56
X_out_d
[i] = (
tanh_beta_eta
+
tanh
(beta * (
X_in_d
[i] - eta))) /
den
;
57
}
58
}
59
63
template
<
typename
T>
64
__global__
void
heaviside_mapping_apply_backward_kernel
(
65
const
T
beta,
const
T
eta,
T
*
__restrict__
sens_out_d
,
66
T
*
__restrict__
sens_in_d
,
T
*
__restrict__
X_in_d
,
const
int
n) {
67
68
const
int
idx
=
blockIdx
.x *
blockDim
.x +
threadIdx
.x;
69
const
int
str
=
blockDim
.x *
gridDim
.x;
70
const
T
den
=
tanh
(beta * eta) +
tanh
(beta * (1.0 - eta));
71
72
for
(
int
i =
idx
; i < n; i +=
str
) {
73
const
T
arg
= beta * (
X_in_d
[i] - eta);
74
const
T
tanh_arg
=
tanh
(
arg
);
75
const
T
dproj
= beta * (1.0 -
tanh_arg
*
tanh_arg
) /
den
;
76
sens_out_d
[i] =
dproj
*
sens_in_d
[i];
77
}
78
}
79
80
#endif
// __NEKO_CUDA_HEAVISIDE_MAPPING_KERNELS__
heaviside_mapping_apply_backward_kernel
__global__ void heaviside_mapping_apply_backward_kernel(const T beta, const T eta, T *__restrict__ sens_out_d, T *__restrict__ sens_in_d, T *__restrict__ X_in_d, const int n)
Definition
heaviside_mapping_kernel.h:64
heaviside_mapping_apply_kernel
__global__ void heaviside_mapping_apply_kernel(const T beta, const T eta, T *__restrict__ X_out_d, T *__restrict__ X_in_d, const int n)
Definition
heaviside_mapping_kernel.h:46
sources
mapping_functions
bcknd
device
cuda
heaviside_mapping_kernel.h
Generated by
1.9.8