Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
PDE_filter_mapping.f90
Go to the documentation of this file.
1
34!
37 use num_types, only: rp
38 use json_module, only: json_file
39 use registry, only: neko_registry
40 use field, only: field_t
41 use coefs, only: coef_t
42 use ax_product, only: ax_t, ax_helm_factory
43 use krylov, only: ksp_t, ksp_monitor_t, krylov_solver_factory
44 use precon, only: pc_t, precon_factory, precon_destroy
45 use bc_list, only: bc_list_t
46 use neumann, only: neumann_t
47 use profiler, only: profiler_start_region, profiler_end_region
48 use gather_scatter, only: gs_t, gs_op_add
49 use pnpn_residual, only: pnpn_prs_res_t
50 use mesh, only: mesh_t, neko_msh_max_zlbls, neko_msh_max_zlbl_len
51 use registry, only: neko_registry
52 use mapping, only: mapping_t
53 use scratch_registry, only: neko_scratch_registry
54 use field_math, only: field_copy, field_add3
55 use coefs, only: coef_t
56 use logger, only: neko_log, log_size
57 use neko_config, only: neko_bcknd_device
58 use dofmap, only: dofmap_t
59 use jacobi, only: jacobi_t
60 use device_jacobi, only: device_jacobi_t
61 use sx_jacobi, only: sx_jacobi_t
62 use utils, only: neko_error
63 use device_math, only: device_cfill, device_subcol3, device_cmult
64 use json_utils, only: json_get, json_get_or_default
65 use continuation_scheduler, only: nekotop_continuation
66 implicit none
67 private
68
73 type, public, extends(mapping_t) :: pde_filter_t
75 class(ax_t), allocatable :: ax
77 type(ksp_monitor_t) :: ksp_results(1)
79 class(ksp_t), allocatable :: ksp_filt
81 class(pc_t), allocatable :: pc_filt
83 type(bc_list_t) :: bclst_filt
84
85 ! Inputs from the user
87 real(kind=rp) :: r
89 real(kind=rp) :: abstol_filt
91 integer :: ksp_max_iter
93 character(len=:), allocatable :: ksp_solver
94 ! > preconditioner type
95 character(len=:), allocatable :: precon_type_filt
96 integer :: ksp_n, n, i
97
98
99
100 contains
102 procedure, pass(this) :: init => pde_filter_init_from_json
104 procedure, pass(this) :: init_from_attributes => &
105 pde_filter_init_from_attributes
107 procedure, pass(this) :: free => pde_filter_free
109 procedure, pass(this) :: forward_mapping => pde_filter_forward_mapping
111 ! TODO
112 ! TALK TO NIELS, I think this is correct...
113 ! but it's not exactly "chain ruling back"
114 ! it's filtering the sensitivity
115
116 ! UPDATE:
117 ! After an email with Niels, we should be using the chain rule,
118 ! not a sensitivity filter
119 procedure, pass(this) :: backward_mapping => pde_filter_backward_mapping
120 end type pde_filter_t
121
122contains
123
125 subroutine pde_filter_init_from_json(this, json, coef)
126 class(pde_filter_t), intent(inout) :: this
127 type(json_file), intent(inout) :: json
128 type(coef_t), intent(inout) :: coef
129 real(kind=rp) :: r, tol
130 integer :: max_iter
131 character(len=:), allocatable :: ksp_solver, precon_type
132
133 call nekotop_continuation%json_get_or_register(json, 'r', this%r, r)
134 call json_get_or_default(json, 'tol', tol, 0.0000000001_rp)
135 call json_get_or_default(json, 'max_iter', max_iter, 200)
136 call json_get_or_default(json, 'solver', ksp_solver, "cg")
137 call json_get_or_default(json, 'preconditioner', precon_type, "jacobi")
138
139 call this%init_base(json, coef, "PDE_filter_mapping")
140 call this%init_from_attributes(coef, r, tol, max_iter, &
141 ksp_solver, precon_type)
142
143 end subroutine pde_filter_init_from_json
144
146 subroutine pde_filter_init_from_attributes(this, coef, r, tol, max_iter, &
147 ksp_solver, precon_type)
148 class(pde_filter_t), intent(inout) :: this
149 type(coef_t), intent(inout) :: coef
150 real(kind=rp), intent(in) :: r, tol
151 integer, intent(in) :: max_iter
152 character(len=*), intent(in) :: ksp_solver, precon_type
153 integer :: n
154
155 this%r = r
156 this%abstol_filt = tol
157 this%ksp_max_iter = max_iter
158 this%ksp_solver = ksp_solver
159 this%precon_type_filt = precon_type
160
161 ! set the number of dofs
162 n = this%coef%dof%size()
163
164 ! init the bc list (all Neuman BCs, will remain empty)
165 call this%bclst_filt%init()
166
167 ! Setup backend dependent Ax routines
168 call ax_helm_factory(this%Ax, full_formulation = .false.)
169
170 ! set up krylov solver
171 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
172 this%ksp_max_iter, this%abstol_filt)
173
174 ! set up preconditioner
175 call filter_precon_factory(this%pc_filt, this%ksp_filt, &
176 this%coef, this%coef%dof, this%coef%gs_h, this%bclst_filt, &
177 this%precon_type_filt)
178
179 end subroutine pde_filter_init_from_attributes
180
182 subroutine pde_filter_free(this)
183 class(pde_filter_t), intent(inout) :: this
184
185 if (allocated(this%Ax)) then
186 deallocate(this%Ax)
187 end if
188
189 if (allocated(this%ksp_filt)) then
190 call this%ksp_filt%free()
191 deallocate(this%ksp_filt)
192 end if
193
194 if (allocated(this%pc_filt)) then
195 call precon_destroy(this%pc_filt)
196 deallocate(this%pc_filt)
197 end if
198
199 call this%bclst_filt%free()
200
201 call this%free_base()
202
203 end subroutine pde_filter_free
204
209 subroutine pde_filter_forward_mapping(this, X_out, X_in)
210 class(pde_filter_t), intent(inout) :: this
211 type(field_t), intent(in) :: X_in
212 type(field_t), intent(inout) :: X_out
213 integer :: n, i
214 type(field_t), pointer :: RHS, d_X_out
215 character(len=LOG_SIZE) :: log_buf
216 integer :: temp_indices(2)
217
218 n = this%coef%dof%size()
219 call neko_scratch_registry%request_field(rhs, temp_indices(1), .false.)
220 call neko_scratch_registry%request_field(d_x_out, temp_indices(2), .false.)
221 ! in a similar fasion to pressure/velocity, we will solve for d_X_out.
222
223 ! to improve convergence, we use X_in as an initial guess for X_out.
224 ! so X_out = X_in + d_X_in.
225
226 ! Defining the operator A = -r^2 \nabla^2 + I
227 ! the system changes from:
228 ! A (X_out) = X_in
229 ! to
230 ! A (d_X_out) = X_in - A(X_in)
231
232 ! set up Helmholtz operators and RHS
233 if (neko_bcknd_device .eq. 1) then
234 call device_cfill(this%coef%h1_d, (this%r / (2.0_rp * sqrt(3.0_rp)))**2, n)
235 call device_cfill(this%coef%h2_d, 1.0_rp, n)
236 else
237 ! h1 is already negative in its definition
238 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
239 ! ax_helm includes the mass matrix in h2
240 this%coef%h2 = 1.0_rp
241 end if
242 this%coef%ifh2 = .true.
243
244 ! compute the A(X_in) component of the RHS
245 ! (note, to be safe with the inout intent we first copy X_in to the
246 ! temporary d_X_out)
247 call field_copy(d_x_out, x_in)
248 call this%Ax%compute(rhs%x, d_x_out%x, this%coef, this%coef%msh, &
249 this%coef%Xh)
250
251 if (neko_bcknd_device .eq. 1) then
252 call device_subcol3(rhs%x_d, x_in%x_d, this%coef%B_d, n)
253 call device_cmult(rhs%x_d, -1.0_rp, n)
254 else
255 do i = 1, n
256 ! mass matrix should be included here
257 rhs%x(i,1,1,1) = x_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
258 - rhs%x(i,1,1,1)
259 end do
260 end if
261
262 ! gather scatter
263 call this%coef%gs_h%op(rhs, gs_op_add)
264
265 ! set BCs
266 call this%bclst_filt%apply_scalar(rhs%x, n)
267
268 ! Solve Helmholtz equation
269 call profiler_start_region('filter solve')
270 this%ksp_results(1) = &
271 this%ksp_filt%solve(this%Ax, d_x_out, rhs%x, n, this%coef, &
272 this%bclst_filt, this%coef%gs_h)
273
274 call profiler_end_region
275
276 ! add result
277 call field_add3(x_out, x_in, d_x_out)
278 ! update preconditioner (needed?)
279 call this%pc_filt%update()
280
281 ! write it all out
282 call neko_log%message('Filter')
283
284 write(log_buf, '(A,A,A)') 'Iterations: ',&
285 'Start residual: ', 'Final residual:'
286 call neko_log%message(log_buf)
287 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
288 this%ksp_results%res_start, this%ksp_results%res_final
289 call neko_log%message(log_buf)
290
291 call neko_scratch_registry%relinquish_field(temp_indices)
292
293
294
295 end subroutine pde_filter_forward_mapping
296
302 ! TODO
303 ! this really confuses me!
304 ! it's not really a chain rule back, it's just a filtering of the sensitivity
305 ! ?
306 !
307 ! Update:
308 ! After an email exchange with Niels:
309 ! We DON'T want to be filtering the sensitivity, this IS just the chain rule.
310 ! So to the best of my knowledge, we're just applying the same filter
311 ! on the sensitivity field.
312 !
313 ! Niels did mention the order of the RHS assembly should be reversed however.
314 ! I'm not exactly sure how this applies to us, but it should be brought up
315 ! in the next group meeting.
316 subroutine pde_filter_backward_mapping(this, sens_out, sens_in, X_in)
317 class(pde_filter_t), intent(inout) :: this
318 type(field_t), intent(in) :: X_in
319 type(field_t), intent(in) :: sens_in
320 type(field_t), intent(inout) :: sens_out
321 integer :: n, i
322 type(field_t), pointer :: RHS, delta ! I'm so sorry for this notation..
323 integer :: temp_indices(2)
324 character(len=LOG_SIZE) :: log_buf
325
326 n = this%coef%dof%size()
327
328 call neko_scratch_registry%request_field(rhs, temp_indices(1), .false.)
329 call neko_scratch_registry%request_field(delta, temp_indices(2), .false.)
330
331 ! set up Helmholtz operators and RHS
332 if (neko_bcknd_device .eq. 1) then
333 call device_cfill(this%coef%h1_d, (this%r / (2.0_rp * sqrt(3.0_rp)))**2, n)
334 call device_cfill(this%coef%h2_d, 1.0_rp, n)
335 else
336 ! h1 is already negative in its definition
337 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
338 ! ax_helm includes the mass matrix in h2
339 this%coef%h2 = 1.0_rp
340 end if
341 this%coef%ifh2 = .true.
342
343 ! compute the A(sens_in) component of the RHS
344 ! (note, to be safe with the inout intent we first copy sens_in to the
345 ! temporary delta)
346 call field_copy(delta, sens_in)
347 call this%Ax%compute(rhs%x, delta%x, this%coef, this%coef%msh, &
348 this%coef%Xh)
349
350 if (neko_bcknd_device .eq. 1) then
351 call device_subcol3(rhs%x_d, sens_in%x_d, this%coef%B_d, n)
352 call device_cmult(rhs%x_d, -1.0_rp, n)
353 else
354 do i = 1, n
355 ! mass matrix should be included here
356 rhs%x(i,1,1,1) = sens_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
357 - rhs%x(i,1,1,1)
358 end do
359 end if
360
361 ! gather scatter
362 call this%coef%gs_h%op(rhs, gs_op_add)
363
364 ! set BCs
365 call this%bclst_filt%apply_scalar(rhs%x, n)
366
367 ! Solve Helmholtz equation
368 call profiler_start_region('filter solve')
369 this%ksp_results(1) = &
370 this%ksp_filt%solve(this%Ax, delta, rhs%x, n, this%coef, &
371 this%bclst_filt, this%coef%gs_h)
372
373 ! add result
374 call field_add3(sens_out, sens_in, delta)
375
376 call profiler_end_region
377
378 ! update preconditioner (needed?)
379 call this%pc_filt%update()
380
381 ! write it all out
382 call neko_log%message('Filter')
383
384 write(log_buf, '(A,A,A)') 'Iterations: ',&
385 'Start residual: ', 'Final residual:'
386 call neko_log%message(log_buf)
387 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
388 this%ksp_results%res_start, this%ksp_results%res_final
389 call neko_log%message(log_buf)
390
391 call neko_scratch_registry%relinquish_field(temp_indices)
392
393 end subroutine pde_filter_backward_mapping
394
395 subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype)
396
397 implicit none
398 class(pc_t), allocatable, target, intent(inout) :: pc
399 class(ksp_t), target, intent(inout) :: ksp
400 type(coef_t), target, intent(in) :: coef
401 type(dofmap_t), target, intent(in) :: dof
402 type(gs_t), target, intent(inout) :: gs
403 type(bc_list_t), target, intent(inout) :: bclst
404 character(len=*) :: pctype
405
406 call precon_factory(pc, pctype)
407
408 select type (pcp => pc)
409 type is (jacobi_t)
410 call pcp%init(coef, dof, gs)
411 type is (sx_jacobi_t)
412 call pcp%init(coef, dof, gs)
413 type is (device_jacobi_t)
414 call pcp%init(coef, dof, gs)
415 end select
416
417 call ksp%set_pc(pc)
418
419 end subroutine filter_precon_factory
420
421end module pde_filter_mapping
Continuation scheduler for the optimization loop.
Mappings to be applied to a scalar field.
Definition mapping.f90:36
A PDE based filter.
subroutine pde_filter_init_from_json(this, json, coef)
Constructor from json.
Base abstract class for mapping.
Definition mapping.f90:47
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .