Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
heaviside_mapping.f90
Go to the documentation of this file.
1
34!
37 use num_types, only: rp
38 use mapping, only: mapping_t
39 use json_module, only: json_file
40 use field, only: field_t
41 use coefs, only: coef_t
42 use neko_config, only: neko_bcknd_device
43 use device_heaviside_mapping, only: &
44 device_heaviside_mapping_apply, &
45 device_heaviside_mapping_apply_backward
46 use heaviside_mapping_cpu, only: &
49 use json_utils, only: json_get_or_default
50 use utils, only: neko_error
51 use continuation_scheduler, only: nekotop_continuation
52
53 implicit none
54 private
55
65 type, public, extends(mapping_t) :: heaviside_mapping_t
67 real(kind=rp) :: beta
69 real(kind=rp) :: eta
70
71 contains
73 procedure, pass(this) :: init => heaviside_mapping_init_from_json
75 procedure, pass(this) :: init_from_attributes => &
76 heaviside_mapping_init_from_attributes
78 procedure, pass(this) :: free => heaviside_mapping_free
80 procedure, pass(this) :: forward_mapping => heaviside_mapping_forward
82 procedure, pass(this) :: backward_mapping => heaviside_mapping_backward
83 end type heaviside_mapping_t
84
85contains
86
88 subroutine heaviside_mapping_init_from_json(this, json, coef)
89 class(heaviside_mapping_t), intent(inout) :: this
90 type(json_file), intent(inout) :: json
91 type(coef_t), intent(inout) :: coef
92 real(kind=rp) :: eta, beta
93
94 ! default value for beta
95 beta = 1.0_rp
96
97 call json_get_or_default(json, 'eta', eta, 0.5_rp)
98
99 call nekotop_continuation%json_get_or_register(json, 'beta', this%beta, &
100 beta)
101
102 ! Initialize base
103 call this%init_base(json, coef, "heaviside_mapping")
104 call this%init_from_attributes(coef, beta, eta)
106
108 subroutine heaviside_mapping_init_from_attributes(this, coef, beta, eta)
109 class(heaviside_mapping_t), intent(inout) :: this
110 type(coef_t), intent(inout) :: coef
111 real(kind=rp), intent(in) :: beta
112 real(kind=rp), intent(in) :: eta
113
114 if (beta .le. 0.0_rp) then
115 call neko_error('"beta" must be > 0 in heaviside mapping')
116 end if
117
118 if (eta .lt. 0.0_rp .or. eta .gt. 1.0_rp) then
119 call neko_error('"eta" must be in [0, 1] in heaviside mapping')
120 end if
121
122 this%beta = beta
123 this%eta = eta
124 end subroutine heaviside_mapping_init_from_attributes
125
127 subroutine heaviside_mapping_free(this)
128 class(heaviside_mapping_t), intent(inout) :: this
129
130 call this%free_base()
131 end subroutine heaviside_mapping_free
132
134 subroutine heaviside_mapping_forward(this, X_out, X_in)
135 class(heaviside_mapping_t), intent(inout) :: this
136 type(field_t), intent(in) :: X_in
137 type(field_t), intent(inout) :: X_out
138 call heaviside_mapping_apply(this%beta, this%eta, x_out, x_in)
139 end subroutine heaviside_mapping_forward
140
142 subroutine heaviside_mapping_backward(this, sens_out, sens_in, X_in)
143 class(heaviside_mapping_t), intent(inout) :: this
144 type(field_t), intent(in) :: X_in
145 type(field_t), intent(in) :: sens_in
146 type(field_t), intent(inout) :: sens_out
147 call heaviside_mapping_apply_backward(this%beta, this%eta, &
148 sens_out, sens_in, x_in)
149 end subroutine heaviside_mapping_backward
150
156 subroutine heaviside_mapping_apply(beta, eta, X_out, X_in)
157 real(kind=rp), intent(in) :: beta, eta
158 type(field_t), intent(in) :: x_in
159 type(field_t), intent(inout) :: X_out
160 integer :: n
161 real(kind=rp) :: den, beta_eta
162
163 beta_eta = beta * eta
164 den = tanh(beta_eta) + tanh(beta * (1.0_rp - eta))
165
166 if (abs(den) .le. tiny(den)) then
167 call neko_error('invalid denominator in heaviside mapping')
168 end if
169
170 n = x_in%dof%size()
171 if (neko_bcknd_device .eq. 1) then
172 call device_heaviside_mapping_apply(beta, eta, &
173 x_out%x_d, x_in%x_d, n)
174 else
175 call heaviside_mapping_apply_cpu(beta, eta, &
176 x_out%x, x_in%x, n)
177 end if
178 end subroutine heaviside_mapping_apply
179
186 subroutine heaviside_mapping_apply_backward(beta, eta, sens_out, sens_in, &
187 X_in)
188 real(kind=rp), intent(in) :: beta, eta
189 type(field_t), intent(in) :: x_in
190 type(field_t), intent(in) :: sens_in
191 type(field_t), intent(inout) :: sens_out
192 integer :: n
193 real(kind=rp) :: den, beta_eta
194
195 beta_eta = beta * eta
196 den = tanh(beta_eta) + tanh(beta * (1.0_rp - eta))
197
198 if (abs(den) .le. tiny(den)) then
199 call neko_error('invalid denominator in heaviside mapping')
200 end if
201
202 n = x_in%dof%size()
203 if (neko_bcknd_device .eq. 1) then
204 call device_heaviside_mapping_apply_backward(beta, eta, &
205 sens_out%x_d, sens_in%x_d, x_in%x_d, n)
206 else
208 sens_out%x, sens_in%x, x_in%x, n)
209 end if
210 end subroutine heaviside_mapping_apply_backward
211
212end module heaviside_mapping
Continuation scheduler for the optimization loop.
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.
Smooth Heaviside mapping.
subroutine heaviside_mapping_init_from_json(this, json, coef)
Constructor from json.
Mappings to be applied to a scalar field.
Definition mapping.f90:36
Base abstract class for mapping.
Definition mapping.f90:46