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! Copyright (c) 2023, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
33!
36 use num_types, only: rp
37 use json_module, only: json_file
38 use field_registry, only: neko_field_registry
39 use field, only: field_t
40 use coefs, only: coef_t
41 use ax_product, only: ax_t, ax_helm_factory
42 use krylov, only: ksp_t, ksp_monitor_t, krylov_solver_factory, &
43 krylov_solver_destroy
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 field_registry, only: neko_field_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 hsmg, only: hsmg_t
63 use utils, only: neko_error
64 use device_math, only: device_cfill, device_subcol3, device_cmult
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 => &
106 procedure, pass(this) :: free => pde_filter_free
108 procedure, pass(this) :: apply_forward => pde_filter_apply
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) :: apply_backward => pde_filter_apply_backward
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
129 ! TODO
130 ! I'll do the json stuff later...
131 this%r = 0.01
132 this%abstol_filt = 0.0000000001_rp
133 this%ksp_max_iter = 200
134 this%ksp_solver = "cg"
135 this%precon_type_filt = "jacobi"
136
137 call this%init_base(json, coef)
138 call pde_filter_init_from_attributes(this, coef)
139
140 end subroutine pde_filter_init_from_json
141
143 subroutine pde_filter_init_from_attributes(this, coef)
144 class(pde_filter_t), intent(inout) :: this
145 type(coef_t), intent(inout) :: coef
146 integer :: n
147
148 ! set the number of dofs
149 n = this%coef%dof%size()
150
151 ! init the bc list (all Neuman BCs, will remain empty)
152 call this%bclst_filt%init()
153
154 ! Setup backend dependent Ax routines
155 call ax_helm_factory(this%Ax, full_formulation = .false.)
156
157 ! set up krylov solver
158 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
159 this%ksp_max_iter, this%abstol_filt)
160
161 ! set up preconditioner
162 call filter_precon_factory(this%pc_filt, this%ksp_filt, &
163 this%coef, this%coef%dof, this%coef%gs_h, this%bclst_filt, &
164 this%precon_type_filt)
165
167
169 subroutine pde_filter_free(this)
170 class(pde_filter_t), intent(inout) :: this
171
172 if (allocated(this%Ax)) then
173 deallocate(this%Ax)
174 end if
175
176 if (allocated(this%ksp_filt)) then
177 call krylov_solver_destroy(this%ksp_filt)
178 deallocate(this%ksp_filt)
179 end if
180
181 if (allocated(this%pc_filt)) then
182 call precon_destroy(this%pc_filt)
183 deallocate(this%pc_filt)
184 end if
185
186 call this%bclst_filt%free()
187
188 call this%free_base()
189
190 end subroutine pde_filter_free
191
195 subroutine pde_filter_apply(this, X_out, X_in)
196 class(pde_filter_t), intent(inout) :: this
197 type(field_t), intent(in) :: X_in
198 type(field_t), intent(inout) :: X_out
199 integer :: n, i
200 type(field_t), pointer :: RHS, d_X_out
201 character(len=LOG_SIZE) :: log_buf
202 integer :: temp_indices(2)
203
204 n = this%coef%dof%size()
205 call neko_scratch_registry%request_field(rhs, temp_indices(1))
206 call neko_scratch_registry%request_field(d_x_out, temp_indices(2))
207 ! in a similar fasion to pressure/velocity, we will solve for d_X_out.
208
209 ! to improve convergence, we use X_in as an initial guess for X_out.
210 ! so X_out = X_in + d_X_in.
211
212 ! Defining the operator A = -r^2 \nabla^2 + I
213 ! the system changes from:
214 ! A (X_out) = X_in
215 ! to
216 ! A (d_X_out) = X_in - A(X_in)
217
218 ! set up Helmholtz operators and RHS
219 if (neko_bcknd_device .eq. 1) then
220 ! TODO
221 ! I think this is correct but I've never tested it
222 call device_cfill(this%coef%h1_d, this%r**2, n)
223 call device_cfill(this%coef%h2_d, 1.0_rp, n)
224 else
225 ! h1 is already negative in its definition
226 this%coef%h1 = this%r**2
227 ! ax_helm includes the mass matrix in h2
228 this%coef%h2 = 1.0_rp
229 end if
230 this%coef%ifh2 = .true.
231
232 ! compute the A(X_in) component of the RHS
233 ! (note, to be safe with the inout intent we first copy X_in to the
234 ! temporary d_X_out)
235 call field_copy(d_x_out, x_in)
236 call this%Ax%compute(rhs%x, d_x_out%x, this%coef, this%coef%msh, &
237 this%coef%Xh)
238
239 if (neko_bcknd_device .eq. 1) then
240 call device_subcol3(rhs%x_d, x_in%x_d, this%coef%B_d, n)
241 call device_cmult(rhs%x_d, -1.0_rp, n)
242 else
243 do i = 1, n
244 ! mass matrix should be included here
245 rhs%x(i,1,1,1) = x_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
246 - rhs%x(i,1,1,1)
247 end do
248 end if
249
250 ! gather scatter
251 call this%coef%gs_h%op(rhs, gs_op_add)
252
253 ! set BCs
254 call this%bclst_filt%apply_scalar(rhs%x, n)
255
256 ! Solve Helmholtz equation
257 call profiler_start_region('filter solve')
258 this%ksp_results(1) = &
259 this%ksp_filt%solve(this%Ax, d_x_out, rhs%x, n, this%coef, &
260 this%bclst_filt, this%coef%gs_h)
261
262 call profiler_end_region
263
264 ! add result
265 call field_add3(x_out, x_in, d_x_out)
266 ! update preconditioner (needed?)
267 call this%pc_filt%update()
268
269 ! write it all out
270 call neko_log%message('Filter')
271
272 write(log_buf, '(A,A,A)') 'Iterations: ',&
273 'Start residual: ', 'Final residual:'
274 call neko_log%message(log_buf)
275 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
276 this%ksp_results%res_start, this%ksp_results%res_final
277 call neko_log%message(log_buf)
278
279 call neko_scratch_registry%relinquish_field(temp_indices)
280
281
282
283 end subroutine pde_filter_apply
284
289 ! TODO
290 ! this really confuses me!
291 ! it's not really a chain rule back, it's just a filtering of the sensitivity
292 ! ?
293 !
294 ! Update:
295 ! After an email exchange with Niels:
296 ! We DON'T want to be filtering the sensitivity, this IS just the chain rule.
297 ! So to the best of my knowledge, we're just applying the same filter
298 ! on the sensitivity field.
299 !
300 ! Niels did mention the order of the RHS assembly should be reversed however.
301 ! I'm not exactly sure how this applies to us, but it should be brought up
302 ! in the next group meeting.
303 subroutine pde_filter_apply_backward(this, dF_dX_in, dF_dX_out, X_in)
304 class(pde_filter_t), intent(inout) :: this
305 type(field_t), intent(in) :: X_in
306 type(field_t), intent(in) :: df_dX_out
307 type(field_t), intent(inout) :: df_dX_in
308 integer :: n, i
309 type(field_t), pointer :: RHS, delta ! I'm so sorry for this notation..
310 integer :: temp_indices(2)
311 character(len=LOG_SIZE) :: log_buf
312
313 n = this%coef%dof%size()
314
315 call neko_scratch_registry%request_field(rhs, temp_indices(1))
316 call neko_scratch_registry%request_field(delta, temp_indices(2))
317
318 ! set up Helmholtz operators and RHS
319 if (neko_bcknd_device .eq. 1) then
320 ! TODO
321 ! I think this is correct but I've never tested it
322 call device_cfill(this%coef%h1_d, this%r**2, n)
323 call device_cfill(this%coef%h2_d, 1.0_rp, n)
324 else
325 ! h1 is already negative in its definition
326 this%coef%h1 = this%r**2
327 ! ax_helm includes the mass matrix in h2
328 this%coef%h2 = 1.0_rp
329 end if
330 this%coef%ifh2 = .true.
331
332 ! compute the A(dF_dX_out) component of the RHS
333 ! (note, to be safe with the inout intent we first copy dF_dX_out to the
334 ! temporary delta)
335 call field_copy(delta, df_dx_out)
336 call this%Ax%compute(rhs%x, delta%x, this%coef, this%coef%msh, &
337 this%coef%Xh)
338
339 if (neko_bcknd_device .eq. 1) then
340 call device_subcol3(rhs%x_d, df_dx_out%x_d, this%coef%B_d, n)
341 call device_cmult(rhs%x_d, -1.0_rp, n)
342 else
343 do i = 1, n
344 ! mass matrix should be included here
345 rhs%x(i,1,1,1) = df_dx_out%x(i,1,1,1) * this%coef%B(i,1,1,1) &
346 - rhs%x(i,1,1,1)
347 end do
348 end if
349
350 ! gather scatter
351 call this%coef%gs_h%op(rhs, gs_op_add)
352
353 ! set BCs
354 call this%bclst_filt%apply_scalar(rhs%x, n)
355
356 ! Solve Helmholtz equation
357 call profiler_start_region('filter solve')
358 this%ksp_results(1) = &
359 this%ksp_filt%solve(this%Ax, delta, rhs%x, n, this%coef, &
360 this%bclst_filt, this%coef%gs_h)
361
362 ! add result
363 call field_add3(df_dx_in, df_dx_out, delta)
364
365 call profiler_end_region
366
367 ! update preconditioner (needed?)
368 call this%pc_filt%update()
369
370 ! write it all out
371 call neko_log%message('Filter')
372
373 write(log_buf, '(A,A,A)') 'Iterations: ',&
374 'Start residual: ', 'Final residual:'
375 call neko_log%message(log_buf)
376 write(log_buf, '(I11,3x, E15.7,5x, E15.7)') this%ksp_results%iter, &
377 this%ksp_results%res_start, this%ksp_results%res_final
378 call neko_log%message(log_buf)
379
380 call neko_scratch_registry%relinquish_field(temp_indices)
381
382 end subroutine pde_filter_apply_backward
383
384 subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype)
385
386 implicit none
387 class(pc_t), allocatable, target, intent(inout) :: pc
388 class(ksp_t), target, intent(inout) :: ksp
389 type(coef_t), target, intent(in) :: coef
390 type(dofmap_t), target, intent(in) :: dof
391 type(gs_t), target, intent(inout) :: gs
392 type(bc_list_t), target, intent(inout) :: bclst
393 character(len=*) :: pctype
394
395 call precon_factory(pc, pctype)
396
397 select type (pcp => pc)
398 type is (jacobi_t)
399 call pcp%init(coef, dof, gs)
400 type is (sx_jacobi_t)
401 call pcp%init(coef, dof, gs)
402 type is (device_jacobi_t)
403 call pcp%init(coef, dof, gs)
404 type is (hsmg_t)
405 if (len_trim(pctype) .gt. 4) then
406 if (index(pctype, '+') .eq. 5) then
407 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, &
408 bclst, trim(pctype(6:)))
409 else
410 call neko_error('Unknown coarse grid solver')
411 end if
412 else
413 call pcp%init(dof%msh, dof%Xh, coef, dof, gs, bclst)
414 end if
415 end select
416
417 call ksp%set_pc(pc)
418
419 end subroutine filter_precon_factory
420
421end module pde_filter
Mappings to be applied to a scalar field.
Definition mapping.f90:35
A PDE based filter.
subroutine pde_filter_free(this)
Destructor.
subroutine pde_filter_init_from_json(this, json, coef)
Constructor from json.
subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype)
subroutine pde_filter_apply_backward(this, df_dx_in, df_dx_out, x_in)
Apply the adjoint filter.
subroutine pde_filter_apply(this, x_out, x_in)
Apply the filter.
subroutine pde_filter_init_from_attributes(this, coef)
Actual constructor.
Base abstract class for mapping.
Definition mapping.f90:44
A PDE based filter mapping $\rho \mapsto \tilde{\rho}$, see Lazarov & O. Sigmund 2010,...