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