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