Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
heaviside_mapping_cpu.f90
Go to the documentation of this file.
1
34!
37 use num_types, only: rp
38 implicit none
39 private
40
43
44contains
45
52 subroutine heaviside_mapping_apply_cpu(beta, eta, X_out, X_in, n)
53 integer, intent(in) :: n
54 real(kind=rp), intent(in) :: beta, eta
55 real(kind=rp), dimension(n), intent(out) :: x_out
56 real(kind=rp), dimension(n), intent(in) :: x_in
57 real(kind=rp) :: den, tanh_beta_eta
58
59 tanh_beta_eta = tanh(beta * eta)
60 den = tanh_beta_eta + tanh(beta * (1.0_rp - eta))
61
62 x_out = heaviside_mapping_kernel(beta, eta, den, tanh_beta_eta, x_in)
63
64 end subroutine heaviside_mapping_apply_cpu
65
74 sens_out, sens_in, X_in, n)
75 integer, intent(in) :: n
76 real(kind=rp), intent(in) :: beta, eta
77 real(kind=rp), dimension(n), intent(out) :: sens_out
78 real(kind=rp), dimension(n), intent(in) :: sens_in
79 real(kind=rp), dimension(n), intent(in) :: x_in
80 real(kind=rp) :: den
81
82 den = tanh(beta * eta) + tanh(beta * (1.0_rp - eta))
83
84 sens_out = heaviside_mapping_backward_kernel(beta, eta, den, &
85 sens_in, x_in)
86
88
96 elemental function heaviside_mapping_kernel(beta, eta, den, &
97 tanh_beta_eta, X_in) result(X_out)
98 real(kind=rp), intent(in) :: beta, eta, den, tanh_beta_eta
99 real(kind=rp), intent(in) :: x_in
100 real(kind=rp) :: x_out
101
102 x_out = (tanh_beta_eta + tanh(beta * (x_in - eta))) / den
103
104 end function heaviside_mapping_kernel
105
113 elemental function heaviside_mapping_backward_kernel(beta, eta, den, &
114 sens_in, X_in) result(sens_out)
115 real(kind=rp), intent(in) :: beta, eta, den
116 real(kind=rp), intent(in) :: sens_in, x_in
117 real(kind=rp) :: sens_out, arg, tanh_arg
118
119 arg = beta * (x_in - eta)
120 tanh_arg = tanh(arg)
121 sens_out = beta * (1.0_rp - tanh_arg * tanh_arg) / den * sens_in
122
123 end function heaviside_mapping_backward_kernel
124
125end module heaviside_mapping_cpu
CPU backend for smooth Heaviside mapping operations.
subroutine, public heaviside_mapping_apply_cpu(beta, eta, x_out, x_in, n)
Apply smooth Heaviside mapping on CPU.
subroutine, public heaviside_mapping_apply_backward_cpu(beta, eta, sens_out, sens_in, x_in, n)
Apply smooth Heaviside mapping chain rule on CPU.