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
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
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
89 real(kind=rp) :: abstol_filt
91 integer :: ksp_max_iter
93 character(len=:),
allocatable :: ksp_solver
95 character(len=:),
allocatable :: precon_type_filt
96 integer :: ksp_n, n, i
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
119 procedure, pass(this) :: backward_mapping => pde_filter_backward_mapping
127 type(json_file),
intent(inout) :: json
128 type(coef_t),
intent(inout) :: coef
129 real(kind=rp) :: r, tol
131 character(len=:),
allocatable :: ksp_solver, precon_type
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")
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)
147 subroutine pde_filter_init_from_attributes(this, coef, r, tol, max_iter, &
148 ksp_solver, precon_type)
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
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
163 n = this%coef%dof%size()
166 call this%bclst_filt%init()
169 call ax_helm_factory(this%Ax, full_formulation = .false.)
172 call krylov_solver_factory(this%ksp_filt, n, this%ksp_solver, &
173 this%ksp_max_iter, this%abstol_filt)
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)
180 end subroutine pde_filter_init_from_attributes
183 subroutine pde_filter_free(this)
186 if (
allocated(this%Ax))
then
190 if (
allocated(this%ksp_filt))
then
191 call this%ksp_filt%free()
192 deallocate(this%ksp_filt)
195 if (
allocated(this%pc_filt))
then
196 call precon_destroy(this%pc_filt)
197 deallocate(this%pc_filt)
200 call this%bclst_filt%free()
202 call this%free_base()
204 end subroutine pde_filter_free
210 subroutine pde_filter_forward_mapping(this, X_out, X_in)
212 type(field_t),
intent(in) :: X_in
213 type(field_t),
intent(inout) :: X_out
215 type(field_t),
pointer :: RHS, d_X_out
216 character(len=LOG_SIZE) :: log_buf
217 integer :: temp_indices(2)
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.)
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)
239 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
241 this%coef%h2 = 1.0_rp
243 this%coef%ifh2 = .true.
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, &
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)
258 rhs%x(i,1,1,1) = x_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
264 call this%coef%gs_h%op(rhs, gs_op_add)
267 call this%bclst_filt%apply_scalar(rhs%x, n)
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)
275 call profiler_end_region
278 call field_add3(x_out, x_in, d_x_out)
280 call this%pc_filt%update()
283 call neko_log%message(
'Filter')
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)
292 call neko_scratch_registry%relinquish_field(temp_indices)
296 end subroutine pde_filter_forward_mapping
317 subroutine pde_filter_backward_mapping(this, sens_out, sens_in, X_in)
319 type(field_t),
intent(in) :: X_in
320 type(field_t),
intent(in) :: sens_in
321 type(field_t),
intent(inout) :: sens_out
323 type(field_t),
pointer :: RHS, delta
324 integer :: temp_indices(2)
325 character(len=LOG_SIZE) :: log_buf
327 n = this%coef%dof%size()
329 call neko_scratch_registry%request_field(rhs, temp_indices(1), .false.)
330 call neko_scratch_registry%request_field(delta, temp_indices(2), .false.)
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)
338 this%coef%h1 = (this%r / (2.0_rp * sqrt(3.0_rp)))**2
340 this%coef%h2 = 1.0_rp
342 this%coef%ifh2 = .true.
347 call field_copy(delta, sens_in)
348 call this%Ax%compute(rhs%x, delta%x, this%coef, this%coef%msh, &
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)
357 rhs%x(i,1,1,1) = sens_in%x(i,1,1,1) * this%coef%B(i,1,1,1) &
363 call this%coef%gs_h%op(rhs, gs_op_add)
366 call this%bclst_filt%apply_scalar(rhs%x, n)
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)
375 call field_add3(sens_out, sens_in, delta)
377 call profiler_end_region
380 call this%pc_filt%update()
383 call neko_log%message(
'Filter')
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)
392 call neko_scratch_registry%relinquish_field(temp_indices)
394 end subroutine pde_filter_backward_mapping
396 subroutine filter_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype)
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
407 call precon_factory(pc, pctype)
409 select type (pcp => pc)
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)
420 end subroutine filter_precon_factory
Continuation scheduler for the optimization loop.
Mappings to be applied to a scalar field.
subroutine pde_filter_init_from_json(this, json, coef)
Constructor from json.
Base abstract class for mapping.
A PDE based filter mapping , see Lazarov & O. Sigmund 2010, by solving an equation of the form .