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 implicit none
49 private
50
78
79 type, public, extends(mapping_t) :: ramp_mapping_t
81 real(kind=rp) :: f_min
83 real(kind=rp) :: f_max
85 real(kind=rp) :: q
88 logical :: convex_up
89
90 contains
92 procedure, pass(this) :: init => ramp_mapping_init_from_json
94 procedure, pass(this) :: init_from_attributes => &
95 ramp_mapping_init_from_attributes
97 procedure, pass(this) :: free => ramp_mapping_free
99 procedure, pass(this) :: forward_mapping => ramp_forward_mapping
101 procedure, pass(this) :: backward_mapping => ramp_backward_mapping
102 end type ramp_mapping_t
103
104contains
105
107 subroutine ramp_mapping_init_from_json(this, json, coef)
108 class(ramp_mapping_t), intent(inout) :: this
109 type(json_file), intent(inout) :: json
110 type(coef_t), intent(inout) :: coef
111 real(kind=rp) :: f_min, f_max, q
112 logical :: convex_up
113
114 call json_get_or_default(json, 'f_min', f_min, 0.0_rp)
115 call json_get(json, 'f_max', f_max)
116 call json_get_or_default(json, 'q', q, 1.0_rp)
117 call json_get_or_default(json, 'convex_up', convex_up, .false.)
118
119 call this%init_base(json, coef)
120 call this%init_from_attributes(coef, f_min, f_max, q, &
121 convex_up)
122
123 end subroutine ramp_mapping_init_from_json
124
126 subroutine ramp_mapping_init_from_attributes(this, coef, f_min, f_max, q, &
127 convex_up)
128 class(ramp_mapping_t), intent(inout) :: this
129 type(coef_t), intent(inout) :: coef
130 real(kind=rp), intent(in) :: f_min, f_max, q
131 logical, intent(in) :: convex_up
132
133 this%f_min = f_min
134 this%f_max = f_max
135 this%q = q
136 this%convex_up = convex_up
137
138 end subroutine ramp_mapping_init_from_attributes
139
141 subroutine ramp_mapping_free(this)
142 class(ramp_mapping_t), intent(inout) :: this
143
144 call this%free_base()
145
146 end subroutine ramp_mapping_free
147
152 subroutine ramp_forward_mapping(this, X_out, X_in)
153 class(ramp_mapping_t), intent(inout) :: this
154 type(field_t), intent(in) :: X_in
155 type(field_t), intent(inout) :: X_out
156
157 if (this%convex_up .eqv. .true.) then
158 call convex_up_ramp_mapping_apply(this%f_min, this%f_max, &
159 this%q, x_out, x_in)
160 else
161 call convex_down_ramp_mapping_apply(this%f_min, this%f_max, &
162 this%q, x_out, x_in)
163 end if
164
165 end subroutine ramp_forward_mapping
166
167
173 subroutine ramp_backward_mapping(this, sens_out, sens_in, X_in)
174 class(ramp_mapping_t), intent(inout) :: this
175 type(field_t), intent(in) :: X_in
176 type(field_t), intent(in) :: sens_in
177 type(field_t), intent(inout) :: sens_out
178
179 if (this%convex_up .eqv. .true.) then
180 call convex_up_ramp_mapping_apply_backward(this%f_min, this%f_max, &
181 this%q, sens_out, sens_in, x_in)
182 else
183 call convex_down_ramp_mapping_apply_backward(this%f_min, this%f_max, &
184 this%q, sens_out, sens_in, x_in)
185 end if
186
187 end subroutine ramp_backward_mapping
188
195 subroutine convex_down_ramp_mapping_apply(f_min, f_max, q, X_out, X_in)
196 real(kind=rp), intent(in) :: q, f_min, f_max
197 type(field_t), intent(in) :: x_in
198 type(field_t), intent(inout) :: X_out
199 integer :: n, i
200
201 ! x_out = f_min + (f_max - f_min) * x_in / (1 + q * (1 - x_in) )
202
203 n = x_in%dof%size()
204 if (neko_bcknd_device .eq. 1) then
205 call device_convex_down_ramp_mapping_apply(f_min, f_max, q, &
206 x_out%x_d, x_in%x_d, n)
207 else
208 do i = 1, n
209 x_out%x(i,1,1,1) = f_min + (f_max - f_min) * &
210 x_in%x(i,1,1,1) / (1.0_rp + q * (1.0_rp - x_in%x(i,1,1,1) ) )
211 end do
212 end if
213
214 end subroutine convex_down_ramp_mapping_apply
215
216
224 subroutine convex_down_ramp_mapping_apply_backward(f_min, f_max, q, &
225 sens_out, sens_in, X_in)
226 real(kind=rp), intent(in) :: f_min, f_max, q
227 type(field_t), intent(in) :: x_in
228 type(field_t), intent(in) :: sens_in
229 type(field_t), intent(inout) :: sens_out
230 integer :: n, i
231
232 ! df/dx_in = df/dx_out * dx_out/dx_in
233
234 ! dx_out/dx_in = (f_min - f_max) * (q + 1) / (1 - q*(x - 1))**2
235
236 n = x_in%dof%size()
237
238 if (neko_bcknd_device .eq. 1) then
239 call device_convex_down_ramp_mapping_apply_backward(f_min, f_max, q, &
240 sens_out%x_d, sens_in%x_d, x_in%x_d, n)
241 else
242 do i = 1, n
243 sens_out%x(i,1,1,1) = (f_max - f_min) * (q + 1.0_rp) / &
244 ((1.0_rp - q * (x_in%x(i,1,1,1) - 1.0_rp))**2) * &
245 sens_in%x(i,1,1,1)
246 end do
247 end if
248
249 end subroutine convex_down_ramp_mapping_apply_backward
250
257 subroutine convex_up_ramp_mapping_apply(f_min, f_max, q, X_out, X_in)
258 real(kind=rp), intent(in) :: f_min, f_max, q
259 type(field_t), intent(in) :: x_in
260 type(field_t), intent(inout) :: X_out
261 integer :: n, i
262
263 ! x_out = f_min + (f_max - f_min) * x_in * (q + 1) / (x_in + q)
264
265 n = x_in%dof%size()
266
267 if (neko_bcknd_device .eq. 1) then
268 call device_convex_up_ramp_mapping_apply(f_min, f_max, q, &
269 x_out%x_d, x_in%x_d, n)
270 else
271 do i = 1, n
272 x_out%x(i,1,1,1) = f_min + (f_max - f_min) * &
273 x_in%x(i,1,1,1) * (1.0_rp + q ) / (x_in%x(i,1,1,1) + q)
274 end do
275 end if
276
277
278 end subroutine convex_up_ramp_mapping_apply
279
280
288 subroutine convex_up_ramp_mapping_apply_backward(f_min, f_max, q, &
289 sens_out, sens_in, X_in)
290 real(kind=rp), intent(in) :: f_min, f_max, q
291 type(field_t), intent(in) :: x_in
292 type(field_t), intent(in) :: sens_in
293 type(field_t), intent(inout) :: sens_out
294 integer :: n, i
295
296 ! df/dx_in = df/dx_out * dx_out/dx_in
297
298 ! dx_out/dx_in = (f_min - f_max) * (q + 1) * q / (q + x)**2
299
300 n = x_in%dof%size()
301
302 if (neko_bcknd_device .eq. 1) then
303 call device_convex_up_ramp_mapping_apply_backward(f_min, f_max, q, &
304 sens_out%x_d, sens_in%x_d, x_in%x_d, n)
305 else
306 do i = 1, n
307 sens_out%x(i,1,1,1) = (f_max - f_min) * (q + 1.0_rp) * q / &
308 ( (x_in%x(i,1,1,1) + q)**2) * &
309 sens_in%x(i,1,1,1)
310 end do
311 end if
312
313 end subroutine convex_up_ramp_mapping_apply_backward
314end 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....