Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
RAMP_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_ramp_mapping, only: device_convex_down_ramp_mapping_apply, &
44 device_convex_down_ramp_mapping_apply_backward, &
45 device_convex_up_ramp_mapping_apply, &
46 device_convex_up_ramp_mapping_apply_backward
47 use json_utils, only: json_get, json_get_or_default
48 use logger, only: neko_log
49 implicit none
50 private
51
79
80 type, public, extends(mapping_t) :: ramp_mapping_t
82 real(kind=rp) :: f_min
84 real(kind=rp) :: f_max
86 real(kind=rp) :: q
89 logical :: convex_up
90
91 contains
93 procedure, pass(this) :: init => ramp_mapping_init_from_json
95 procedure, pass(this) :: init_from_attributes => &
96 ramp_mapping_init_from_attributes
98 procedure, pass(this) :: free => ramp_mapping_free
100 procedure, pass(this) :: forward_mapping => ramp_forward_mapping
102 procedure, pass(this) :: backward_mapping => ramp_backward_mapping
103 end type ramp_mapping_t
104
105contains
106
108 subroutine ramp_mapping_init_from_json(this, json, coef)
109 class(ramp_mapping_t), intent(inout) :: this
110 type(json_file), intent(inout) :: json
111 type(coef_t), intent(inout) :: coef
112 real(kind=rp) :: f_min, f_max, q
113 logical :: convex_up
114
115 call json_get_or_default(json, 'f_min', f_min, 0.0_rp)
116 call json_get(json, 'f_max', f_max)
117 call json_get_or_default(json, 'q', q, 1.0_rp)
118 call json_get_or_default(json, 'convex_up', convex_up, .false.)
119
120 call this%init_base(json, coef)
121 call this%init_from_attributes(coef, f_min, f_max, q, &
122 convex_up)
123
124 end subroutine ramp_mapping_init_from_json
125
127 subroutine ramp_mapping_init_from_attributes(this, coef, f_min, f_max, q, &
128 convex_up)
129 class(ramp_mapping_t), intent(inout) :: this
130 type(coef_t), intent(inout) :: coef
131 real(kind=rp), intent(in) :: f_min, f_max, q
132 logical, intent(in) :: convex_up
133 character(len=256) :: msg
134
135 this%f_min = f_min
136 this%f_max = f_max
137 this%q = q
138 this%convex_up = convex_up
139
140 call neko_log%section('RAMP Mapping')
141 write(msg, '(A,F8.4)') ' f_min: ', this%f_min
142 call neko_log%message(msg)
143 write(msg, '(A,F8.4)') ' f_max: ', this%f_max
144 call neko_log%message(msg)
145 write(msg, '(A,F8.4)') ' q: ', this%q
146 call neko_log%message(msg)
147 if (this%convex_up .eqv. .true.) then
148 call neko_log%message(' convexity: up (Borrvall & Peterson)')
149 else
150 call neko_log%message(' convexity: down (standard RAMP)')
151 end if
152 call neko_log%end_section()
153
154 end subroutine ramp_mapping_init_from_attributes
155
157 subroutine ramp_mapping_free(this)
158 class(ramp_mapping_t), intent(inout) :: this
159
160 call this%free_base()
161
162 end subroutine ramp_mapping_free
163
168 subroutine ramp_forward_mapping(this, X_out, X_in)
169 class(ramp_mapping_t), intent(inout) :: this
170 type(field_t), intent(in) :: X_in
171 type(field_t), intent(inout) :: X_out
172
173 if (this%convex_up .eqv. .true.) then
174 call convex_up_ramp_mapping_apply(this%f_min, this%f_max, &
175 this%q, x_out, x_in)
176 else
177 call convex_down_ramp_mapping_apply(this%f_min, this%f_max, &
178 this%q, x_out, x_in)
179 end if
180
181 end subroutine ramp_forward_mapping
182
183
189 subroutine ramp_backward_mapping(this, sens_out, sens_in, X_in)
190 class(ramp_mapping_t), intent(inout) :: this
191 type(field_t), intent(in) :: X_in
192 type(field_t), intent(in) :: sens_in
193 type(field_t), intent(inout) :: sens_out
194
195 if (this%convex_up .eqv. .true.) then
196 call convex_up_ramp_mapping_apply_backward(this%f_min, this%f_max, &
197 this%q, sens_out, sens_in, x_in)
198 else
199 call convex_down_ramp_mapping_apply_backward(this%f_min, this%f_max, &
200 this%q, sens_out, sens_in, x_in)
201 end if
202
203 end subroutine ramp_backward_mapping
204
211 subroutine convex_down_ramp_mapping_apply(f_min, f_max, q, X_out, X_in)
212 real(kind=rp), intent(in) :: q, f_min, f_max
213 type(field_t), intent(in) :: x_in
214 type(field_t), intent(inout) :: X_out
215 integer :: n, i
216
217 ! x_out = f_min + (f_max - f_min) * x_in / (1 + q * (1 - x_in) )
218
219 n = x_in%dof%size()
220 if (neko_bcknd_device .eq. 1) then
221 call device_convex_down_ramp_mapping_apply(f_min, f_max, q, &
222 x_out%x_d, x_in%x_d, n)
223 else
224 do i = 1, n
225 x_out%x(i,1,1,1) = f_min + (f_max - f_min) * &
226 x_in%x(i,1,1,1) / (1.0_rp + q * (1.0_rp - x_in%x(i,1,1,1) ) )
227 end do
228 end if
229
230 end subroutine convex_down_ramp_mapping_apply
231
232
240 subroutine convex_down_ramp_mapping_apply_backward(f_min, f_max, q, &
241 sens_out, sens_in, X_in)
242 real(kind=rp), intent(in) :: f_min, f_max, q
243 type(field_t), intent(in) :: x_in
244 type(field_t), intent(in) :: sens_in
245 type(field_t), intent(inout) :: sens_out
246 integer :: n, i
247
248 ! df/dx_in = df/dx_out * dx_out/dx_in
249
250 ! dx_out/dx_in = (f_min - f_max) * (q + 1) / (1 - q*(x - 1))**2
251
252 n = x_in%dof%size()
253
254 if (neko_bcknd_device .eq. 1) then
255 call device_convex_down_ramp_mapping_apply_backward(f_min, f_max, q, &
256 sens_out%x_d, sens_in%x_d, x_in%x_d, n)
257 else
258 do i = 1, n
259 sens_out%x(i,1,1,1) = (f_max - f_min) * (q + 1.0_rp) / &
260 ((1.0_rp - q * (x_in%x(i,1,1,1) - 1.0_rp))**2) * &
261 sens_in%x(i,1,1,1)
262 end do
263 end if
264
265 end subroutine convex_down_ramp_mapping_apply_backward
266
273 subroutine convex_up_ramp_mapping_apply(f_min, f_max, q, X_out, X_in)
274 real(kind=rp), intent(in) :: f_min, f_max, q
275 type(field_t), intent(in) :: x_in
276 type(field_t), intent(inout) :: X_out
277 integer :: n, i
278
279 ! x_out = f_min + (f_max - f_min) * x_in * (q + 1) / (x_in + q)
280
281 n = x_in%dof%size()
282
283 if (neko_bcknd_device .eq. 1) then
284 call device_convex_up_ramp_mapping_apply(f_min, f_max, q, &
285 x_out%x_d, x_in%x_d, n)
286 else
287 do i = 1, n
288 x_out%x(i,1,1,1) = f_min + (f_max - f_min) * &
289 x_in%x(i,1,1,1) * (1.0_rp + q ) / (x_in%x(i,1,1,1) + q)
290 end do
291 end if
292
293
294 end subroutine convex_up_ramp_mapping_apply
295
296
304 subroutine convex_up_ramp_mapping_apply_backward(f_min, f_max, q, &
305 sens_out, sens_in, X_in)
306 real(kind=rp), intent(in) :: f_min, f_max, q
307 type(field_t), intent(in) :: x_in
308 type(field_t), intent(in) :: sens_in
309 type(field_t), intent(inout) :: sens_out
310 integer :: n, i
311
312 ! df/dx_in = df/dx_out * dx_out/dx_in
313
314 ! dx_out/dx_in = (f_min - f_max) * (q + 1) * q / (q + x)**2
315
316 n = x_in%dof%size()
317
318 if (neko_bcknd_device .eq. 1) then
319 call device_convex_up_ramp_mapping_apply_backward(f_min, f_max, q, &
320 sens_out%x_d, sens_in%x_d, x_in%x_d, n)
321 else
322 do i = 1, n
323 sens_out%x(i,1,1,1) = (f_max - f_min) * (q + 1.0_rp) * q / &
324 ( (x_in%x(i,1,1,1) + q)**2) * &
325 sens_in%x(i,1,1,1)
326 end do
327 end if
328
329 end subroutine convex_up_ramp_mapping_apply_backward
330end module ramp_mapping
Mappings to be applied to a scalar field.
Definition mapping.f90:36
A RAMP mapping of coefficients.
subroutine ramp_mapping_init_from_json(this, json, coef)
Constructor from json.
Base abstract class for mapping.
Definition mapping.f90:46
A RAMP mapping of coefficients This is the standard RAMP described in https://doi....