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 r = 1.0_rp
134 call nekotop_continuation%json_get_or_register(json, 'r', this%r, r)
135 call json_get_or_default(json, 'tol', tol, 0.0000000001_rp)
136 call json_get_or_default(json, 'max_iter', max_iter, 200)
137 call json_get_or_default(json, 'solver', ksp_solver, "cg")
138 call json_get_or_default(json, 'preconditioner', precon_type, "jacobi")
139
140 call this%init_base(json, coef, "PDE_filter_mapping")
141 call this%init_from_attributes(coef, r, tol, max_iter, &
142 ksp_solver, precon_type)
143
144 end subroutine pde_filter_init_from_json
145
147 subroutine pde_filter_init_from_attributes(this, coef, r, tol, max_iter, &
148 ksp_solver, precon_type)
149 class(pde_filter_t), intent(inout) :: this
150 type(coef_t), intent(inout) :: coef
151 real(kind=rp), intent(in) :: r, tol
152 integer, intent(in) :: max_iter
153 character(len=*), intent(in) :: ksp_solver, precon_type
154 integer :: n
155
156 this%r = r
157 this%abstol_filt = tol
158 this%ksp_max_iter = max_iter
159 this%ksp_solver = ksp_solver
160 this%precon_type_filt = precon_type
161
162 ! set the number of dofs
163 n = this%coef%dof%size()
164
165 ! init the bc list (all Neuman BCs, will remain empty)
166 call this%bclst_filt%init()
167
168 ! Setup backend dependent Ax routines
169 call ax_helm_factory(this%Ax, full_formulation = .false.)
170
171 ! set up krylov solver
172 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
173 this%ksp_max_iter, this%abstol_filt)
174
175 ! set up preconditioner
176 call filter_precon_factory(this%pc_filt, this%ksp_filt, &
177 this%coef, this%coef%dof, this%coef%gs_h, this%bclst_filt, &
178 this%precon_type_filt)
179
180 end subroutine pde_filter_init_from_attributes
181
183 subroutine pde_filter_free(this)
184 class(pde_filter_t), intent(inout) :: this
185
186 if (allocated(this%Ax)) then
187 deallocate(this%Ax)
188 end if
189
190 if (allocated(this%ksp_filt)) then
191 call this%ksp_filt%free()
192 deallocate(this%ksp_filt)
193 end if
194
195 if (allocated(this%pc_filt)) then
196 call precon_destroy(this%pc_filt)
197 deallocate(this%pc_filt)
198 end if
199
200 call this%bclst_filt%free()
201
202 call this%free_base()
203
204 end subroutine pde_filter_free
205
210 subroutine pde_filter_forward_mapping(this, X_out, X_in)
211 class(pde_filter_t), intent(inout) :: this
212 type(field_t), intent(in) :: X_in
213 type(field_t), intent(inout) :: X_out
214 integer :: n, i
215 type(field_t), pointer :: RHS, d_X_out
216 character(len=LOG_SIZE) :: log_buf
217 integer :: temp_indices(2)
218
219 n = this%coef%dof%size()
220 call neko_scratch_registry%request_field(rhs, temp_indices(1), .false.)
221 call neko_scratch_registry%request_field(d_x_out, temp_indices(2), .false.)
222 ! in a similar fasion to pressure/velocity, we will solve for d_X_out.
223
224 ! to improve convergence, we use X_in as an initial guess for X_out.
225 ! so X_out = X_in + d_X_in.
226
227 ! Defining the operator A = -r^2 \nabla^2 + I
228 ! the system changes from:
229 ! A (X_out) = X_in
230 ! to
231 ! A (d_X_out) = X_in - A(X_in)
232
233 ! set up Helmholtz operators and RHS
234 if (neko_bcknd_device .eq. 1) then
235 call device_cfill(this%coef%h1_d, (this%r / (2.0_rp * sqrt(3.0_rp)))**2, n)
236 call device_cfill(this%coef%h2_d, 1.0_rp, n)
237 else
238 ! h1 is already negative in its definition
239 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
240 ! ax_helm includes the mass matrix in h2
241 this%coef%h2 = 1.0_rp
242 end if
243 this%coef%ifh2 = .true.
244
245 ! compute the A(X_in) component of the RHS
246 ! (note, to be safe with the inout intent we first copy X_in to the
247 ! temporary d_X_out)
248 call field_copy(d_x_out, x_in)
249 call this%Ax%compute(rhs%x, d_x_out%x, this%coef, this%coef%msh, &
250 this%coef%Xh)
251
252 if (neko_bcknd_device .eq. 1) then
253 call device_subcol3(rhs%x_d, x_in%x_d, this%coef%B_d, n)
254 call device_cmult(rhs%x_d, -1.0_rp, n)
255 else
256 do i = 1, n
257 ! mass matrix should be included here
258 rhs%x(i,1,1,1) = x_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
259 - rhs%x(i,1,1,1)
260 end do
261 end if
262
263 ! gather scatter
264 call this%coef%gs_h%op(rhs, gs_op_add)
265
266 ! set BCs
267 call this%bclst_filt%apply_scalar(rhs%x, n)
268
269 ! Solve Helmholtz equation
270 call profiler_start_region('filter solve')
271 this%ksp_results(1) = &
272 this%ksp_filt%solve(this%Ax, d_x_out, rhs%x, n, this%coef, &
273 this%bclst_filt, this%coef%gs_h)
274
275 call profiler_end_region
276
277 ! add result
278 call field_add3(x_out, x_in, d_x_out)
279 ! update preconditioner (needed?)
280 call this%pc_filt%update()
281
282 ! write it all out
283 call neko_log%message('Filter')
284
285 write(log_buf, '(A,A,A)') 'Iterations: ',&
286 'Start residual: ', 'Final residual:'
287 call neko_log%message(log_buf)
288 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
289 this%ksp_results%res_start, this%ksp_results%res_final
290 call neko_log%message(log_buf)
291
292 call neko_scratch_registry%relinquish_field(temp_indices)
293
294
295
296 end subroutine pde_filter_forward_mapping
297
303 ! TODO
304 ! this really confuses me!
305 ! it's not really a chain rule back, it's just a filtering of the sensitivity
306 ! ?
307 !
308 ! Update:
309 ! After an email exchange with Niels:
310 ! We DON'T want to be filtering the sensitivity, this IS just the chain rule.
311 ! So to the best of my knowledge, we're just applying the same filter
312 ! on the sensitivity field.
313 !
314 ! Niels did mention the order of the RHS assembly should be reversed however.
315 ! I'm not exactly sure how this applies to us, but it should be brought up
316 ! in the next group meeting.
317 subroutine pde_filter_backward_mapping(this, sens_out, sens_in, X_in)
318 class(pde_filter_t), intent(inout) :: this
319 type(field_t), intent(in) :: X_in
320 type(field_t), intent(in) :: sens_in
321 type(field_t), intent(inout) :: sens_out
322 integer :: n, i
323 type(field_t), pointer :: RHS, delta ! I'm so sorry for this notation..
324 integer :: temp_indices(2)
325 character(len=LOG_SIZE) :: log_buf
326
327 n = this%coef%dof%size()
328
329 call neko_scratch_registry%request_field(rhs, temp_indices(1), .false.)
330 call neko_scratch_registry%request_field(delta, temp_indices(2), .false.)
331
332 ! set up Helmholtz operators and RHS
333 if (neko_bcknd_device .eq. 1) then
334 call device_cfill(this%coef%h1_d, (this%r / (2.0_rp * sqrt(3.0_rp)))**2, n)
335 call device_cfill(this%coef%h2_d, 1.0_rp, n)
336 else
337 ! h1 is already negative in its definition
338 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
339 ! ax_helm includes the mass matrix in h2
340 this%coef%h2 = 1.0_rp
341 end if
342 this%coef%ifh2 = .true.
343
344 ! compute the A(sens_in) component of the RHS
345 ! (note, to be safe with the inout intent we first copy sens_in to the
346 ! temporary delta)
347 call field_copy(delta, sens_in)
348 call this%Ax%compute(rhs%x, delta%x, this%coef, this%coef%msh, &
349 this%coef%Xh)
350
351 if (neko_bcknd_device .eq. 1) then
352 call device_subcol3(rhs%x_d, sens_in%x_d, this%coef%B_d, n)
353 call device_cmult(rhs%x_d, -1.0_rp, n)
354 else
355 do i = 1, n
356 ! mass matrix should be included here
357 rhs%x(i,1,1,1) = sens_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
358 - rhs%x(i,1,1,1)
359 end do
360 end if
361
362 ! gather scatter
363 call this%coef%gs_h%op(rhs, gs_op_add)
364
365 ! set BCs
366 call this%bclst_filt%apply_scalar(rhs%x, n)
367
368 ! Solve Helmholtz equation
369 call profiler_start_region('filter solve')
370 this%ksp_results(1) = &
371 this%ksp_filt%solve(this%Ax, delta, rhs%x, n, this%coef, &
372 this%bclst_filt, this%coef%gs_h)
373
374 ! add result
375 call field_add3(sens_out, sens_in, delta)
376
377 call profiler_end_region
378
379 ! update preconditioner (needed?)
380 call this%pc_filt%update()
381
382 ! write it all out
383 call neko_log%message('Filter')
384
385 write(log_buf, '(A,A,A)') 'Iterations: ',&
386 'Start residual: ', 'Final residual:'
387 call neko_log%message(log_buf)
388 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
389 this%ksp_results%res_start, this%ksp_results%res_final
390 call neko_log%message(log_buf)
391
392 call neko_scratch_registry%relinquish_field(temp_indices)
393
394 end subroutine pde_filter_backward_mapping
395
396 subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype)
397
398 implicit none
399 class(pc_t), allocatable, target, intent(inout) :: pc
400 class(ksp_t), target, intent(inout) :: ksp
401 type(coef_t), target, intent(in) :: coef
402 type(dofmap_t), target, intent(in) :: dof
403 type(gs_t), target, intent(inout) :: gs
404 type(bc_list_t), target, intent(inout) :: bclst
405 character(len=*) :: pctype
406
407 call precon_factory(pc, pctype)
408
409 select type (pcp => pc)
410 type is (jacobi_t)
411 call pcp%init(coef, dof, gs)
412 type is (sx_jacobi_t)
413 call pcp%init(coef, dof, gs)
414 type is (device_jacobi_t)
415 call pcp%init(coef, dof, gs)
416 end select
417
418 call ksp%set_pc(pc)
419
420 end subroutine filter_precon_factory
421
422end module pde_filter
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:46
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .