Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_fluid_pnpn.f90
Go to the documentation of this file.
1
34!
37 use, intrinsic :: iso_fortran_env, only: error_unit
38 use coefs, only: coef_t
39 use symmetry, only: symmetry_t
40 use registry, only: neko_registry
41 use logger, only: neko_log, log_size, neko_log_debug
42 use num_types, only: rp
43 use krylov, only: ksp_monitor_t
44 use adjoint_pnpn_residual, only: adjoint_pnpn_prs_res_t, &
45 adjoint_pnpn_vel_res_t, adjoint_pnpn_prs_res_factory, &
46 adjoint_pnpn_vel_res_factory
47 use rhs_maker, only: rhs_maker_sumab_t, rhs_maker_bdf_t, rhs_maker_ext_t, &
48 rhs_maker_oifs_t, rhs_maker_sumab_fctry, rhs_maker_bdf_fctry, &
49 rhs_maker_ext_fctry, rhs_maker_oifs_fctry
52 use device_mathops, only: device_opadd2cm
53 use fluid_aux, only: fluid_step_info
54 use scratch_registry, only: neko_scratch_registry
55 use projection, only: projection_t
56 use projection_vel, only: projection_vel_t
57 use device, only: device_memcpy, host_to_device, device_event_sync, &
58 glb_cmd_event
59 use advection_adjoint, only: advection_adjoint_t, advection_adjoint_factory
60 use profiler, only: profiler_start_region, profiler_end_region
61 use json_module, only: json_file, json_core, json_value
62 use json_utils, only: json_get, json_get_or_default, json_extract_item
63 use ax_product, only: ax_t, ax_helm_factory
64 use field, only: field_t
65 use dirichlet, only: dirichlet_t
66 use shear_stress, only: shear_stress_t
67 use wall_model_bc, only: wall_model_bc_t
68 use facet_normal, only: facet_normal_t
69 use non_normal, only: non_normal_t
70 use checkpoint, only: chkp_t
71 use mesh, only: mesh_t
72 use user_intf, only: user_t
73 use time_step_controller, only: time_step_controller_t
74 use gs_ops, only: gs_op_add
75 use neko_config, only: neko_bcknd_device
76 use mathops, only: opadd2cm
77 use bc_list, only: bc_list_t
78 use zero_dirichlet, only: zero_dirichlet_t
79 use utils, only: neko_error
80 use field_math, only: field_add2, field_copy, &
81 field_add2s2
82 use bc, only: bc_t
83 use file, only: file_t
84 use operators, only: ortho
85 use opr_device, only: device_ortho
86 use inflow, only: inflow_t
87 use field_dirichlet, only: field_dirichlet_t
88 use blasius, only: blasius_t
89 use field_dirichlet_vector, only: field_dirichlet_vector_t
90 use dong_outflow, only: dong_outflow_t
91 use time_state, only: time_state_t
92 use vector, only: vector_t
93 use device_math, only: device_vlsc3, device_cmult, device_col2
94 use math, only: vlsc3, cmult, col2
95 use, intrinsic :: iso_c_binding, only: c_ptr, c_null_ptr, c_associated
96 use comm, only: neko_comm, mpi_real_precision
97 use mpi_f08, only: mpi_sum, mpi_max, mpi_allreduce, mpi_comm_world, &
98 mpi_integer, mpi_logical, mpi_lor
99 use operators, only : opgrad, curl, grad
100 use normal_vec_bcs, only: normal_vec_bcs_t
101
102 implicit none
103 private
104
105 type, public, extends(adjoint_fluid_scheme_incompressible_t) :: &
107
109 type(field_t) :: p_res, u_res, v_res, w_res
110
113 type(field_t) :: dp, du, dv, dw
114
115 !
116 ! Implicit operators, i.e. the left-hand-side of the Helmholz problem.
117 !
118
119 ! Coupled Helmholz operator for velocity
120 class(ax_t), allocatable :: ax_vel
121 ! Helmholz operator for pressure
122 class(ax_t), allocatable :: ax_prs
123
124 !
125 ! Projections for solver speed-up
126 !
127
129 type(projection_t) :: proj_prs
130 type(projection_vel_t) :: proj_vel
131
132 !
133 ! Special Karniadakis scheme boundary conditions in the pressure equation
134 !
135
137 type(facet_normal_t) :: bc_prs_surface
138
140 type(facet_normal_t) :: bc_sym_surface
141
143 type(normal_vec_bcs_t) :: bc_curl_curl
144
145 !
146 ! Boundary conditions and lists for residuals and solution increments
147 !
148
150 type(zero_dirichlet_t) :: bc_vel_res
152 type(zero_dirichlet_t) :: bc_du
154 type(zero_dirichlet_t) :: bc_dv
156 type(zero_dirichlet_t) :: bc_dw
158 type(zero_dirichlet_t) :: bc_dp
159
161 type(bc_list_t) :: bclst_vel_res
162 type(bc_list_t) :: bclst_du
163 type(bc_list_t) :: bclst_dv
164 type(bc_list_t) :: bclst_dw
165 type(bc_list_t) :: bclst_dp
166
167
168 ! Checker for wether we have a strong pressure bc. If not, the pressure
169 ! is demeaned at every time step.
170 logical :: prs_dirichlet = .false.
171
172
173 ! The advection operator.
174 class(advection_adjoint_t), allocatable :: adv
175
176 ! Time OIFS interpolation scheme for advection.
177 logical :: oifs
178
179 ! Time variables
180 type(field_t) :: abx1, aby1, abz1
181 type(field_t) :: abx2, aby2, abz2
182
183 ! Advection terms for the oifs method
184 type(field_t) :: advx, advy, advz
185
187 class(adjoint_pnpn_prs_res_t), allocatable :: prs_res
188
190 class(adjoint_pnpn_vel_res_t), allocatable :: vel_res
191
193 class(rhs_maker_sumab_t), allocatable :: sumab
194
196 class(rhs_maker_ext_t), allocatable :: makeabf
197
199 class(rhs_maker_bdf_t), allocatable :: makebdf
200
202 class(rhs_maker_oifs_t), allocatable :: makeoifs
203
204 ! !> Adjust flow volume
205 ! type(fluid_volflow_t) :: vol_flow
206
208 logical :: full_stress_formulation = .false.
209
210 ! ======================================================================= !
211 ! Addressable attributes
212
213 real(kind=rp) :: norm_scaling
214 real(kind=rp) :: norm_target
215 real(kind=rp) :: norm_tolerance
216
217 ! ======================================================================= !
218 ! Definition of shorthands and local variables
219
221 real(kind=rp) :: norm_l2_base
222
224 real(kind=rp) :: norm_l2_upper
226 real(kind=rp) :: norm_l2_lower
227
229 type(file_t) :: file_output
230
231 contains
233 procedure, pass(this) :: init => adjoint_fluid_pnpn_init
235 procedure, pass(this) :: free => adjoint_fluid_pnpn_free
237 procedure, pass(this) :: step => adjoint_fluid_pnpn_step
239 procedure, pass(this) :: restart => adjoint_fluid_pnpn_restart
241 procedure, pass(this) :: setup_bcs => adjoint_fluid_pnpn_setup_bcs
243 procedure, pass(this) :: write_boundary_conditions => &
244 adjoint_fluid_pnpn_write_boundary_conditions
245
247 procedure, public, pass(this) :: pw_compute_ => power_iterations_compute
248
249 end type adjoint_fluid_pnpn_t
250
259 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
260 class(bc_t), pointer, intent(inout) :: object
261 type(adjoint_fluid_pnpn_t), intent(in) :: scheme
262 type(json_file), intent(inout) :: json
263 type(coef_t), intent(in) :: coef
264 type(user_t), intent(in) :: user
265 end subroutine pressure_bc_factory
266 end interface pressure_bc_factory
267
275 interface velocity_bc_factory
276 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
277 class(bc_t), pointer, intent(inout) :: object
278 type(adjoint_fluid_pnpn_t), intent(in) :: scheme
279 type(json_file), intent(inout) :: json
280 type(coef_t), intent(in) :: coef
281 type(user_t), intent(in) :: user
282 end subroutine velocity_bc_factory
283 end interface velocity_bc_factory
284
285contains
286
287 subroutine adjoint_fluid_pnpn_init(this, msh, lx, params, user, chkp)
288 class(adjoint_fluid_pnpn_t), target, intent(inout) :: this
289 type(mesh_t), target, intent(inout) :: msh
290 integer, intent(in) :: lx
291 type(json_file), target, intent(inout) :: params
292 type(user_t), target, intent(in) :: user
293 type(chkp_t), target, intent(inout) :: chkp
294 character(len=15), parameter :: scheme = 'Adjoint (Pn/Pn)'
295 real(kind=rp) :: abs_tol
296 character(len=LOG_SIZE) :: log_buf
297 integer :: integer_val, solver_maxiter
298 character(len=:), allocatable :: solver_type, precon_type
299 logical :: monitor, found
300 logical :: advection
301 type(json_file) :: numerics_params, precon_params
302
303 ! Temporary field pointers
304 character(len=:), allocatable :: file_name
305 character(len=256) :: header_line
306
307 call this%free()
308
309 ! Initialize base class
310 call this%init_base(msh, lx, params, scheme, user, .true.)
311
312 ! Add pressure field to the registry. For this scheme it is in the same
313 ! Xh as the velocity
314 call neko_registry%add_field(this%dm_Xh, 'p_adj')
315 this%p_adj => neko_registry%get_field('p_adj')
316
317 !
318 ! Select governing equations via associated residual and Ax types
319 !
320
321 call json_get(params, 'case.numerics.time_order', integer_val)
322 allocate(this%ext_bdf)
323 call this%ext_bdf%init(integer_val)
324
325 call json_get_or_default(params, "case.fluid.full_stress_formulation", &
326 this%full_stress_formulation, .false.)
327
328 if (this%full_stress_formulation .eqv. .true.) then
329 call neko_error( &
330 "Full stress formulation is not supported in the adjoint module.")
331 ! ! Setup backend dependent Ax routines
332 ! call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
333
334 ! ! Setup backend dependent prs residual routines
335 ! call pnpn_prs_res_stress_factory(this%prs_res)
336
337 ! ! Setup backend dependent vel residual routines
338 ! call pnpn_vel_res_stress_factory(this%vel_res)
339 else
340 ! Setup backend dependent Ax routines
341 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
342
343 ! Setup backend dependent prs residual routines
344 call adjoint_pnpn_prs_res_factory(this%prs_res)
345
346 ! Setup backend dependent vel residual routines
347 call adjoint_pnpn_vel_res_factory(this%vel_res)
348 end if
349
350 if (params%valid_path('case.fluid.nut_field')) then
351 if (this%full_stress_formulation .eqv. .false.) then
352 call neko_error("You need to set full_stress_formulation to " // &
353 "true for the fluid to have a spatially varying " // &
354 "viscocity field.")
355 end if
356 call json_get(params, 'case.fluid.nut_field', this%nut_field_name)
357 else
358 this%nut_field_name = ""
359 end if
360
361 ! Setup Ax for the pressure
362 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
363
364
365 ! Setup backend dependent summation of AB/BDF
366 call rhs_maker_sumab_fctry(this%sumab)
367
368 ! Setup backend dependent summation of extrapolation scheme
369 call rhs_maker_ext_fctry(this%makeabf)
370
371 ! Setup backend depenent contributions to F from lagged BD terms
372 call rhs_maker_bdf_fctry(this%makebdf)
373
374 ! Setup backend dependent summations of the OIFS method
375 call rhs_maker_oifs_fctry(this%makeoifs)
376
377 ! Initialize variables specific to this plan
378 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
379 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
380
381 call this%p_res%init(dm_xh, "p_res")
382 call this%u_res%init(dm_xh, "u_res")
383 call this%v_res%init(dm_xh, "v_res")
384 call this%w_res%init(dm_xh, "w_res")
385 call this%abx1%init(dm_xh, "abx1")
386 call this%aby1%init(dm_xh, "aby1")
387 call this%abz1%init(dm_xh, "abz1")
388 call this%abx2%init(dm_xh, "abx2")
389 call this%aby2%init(dm_xh, "aby2")
390 call this%abz2%init(dm_xh, "abz2")
391 call this%advx%init(dm_xh, "advx")
392 call this%advy%init(dm_xh, "advy")
393 call this%advz%init(dm_xh, "advz")
394 end associate
395
396 call this%du%init(this%dm_Xh, 'du')
397 call this%dv%init(this%dm_Xh, 'dv')
398 call this%dw%init(this%dm_Xh, 'dw')
399 call this%dp%init(this%dm_Xh, 'dp')
400
401 ! Set up boundary conditions
402 call this%setup_bcs(user, params)
403
404 ! Check if we need to output boundaries
405 call json_get_or_default(params, 'case.output_boundary', found, .false.)
406 if (found) call this%write_boundary_conditions()
407
408 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
409 this%pr_projection_activ_step)
410
411 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
412 this%vel_projection_activ_step)
413
414 ! Determine the time-interpolation scheme
415 call json_get_or_default(params, 'case.numerics.oifs', this%oifs, .false.)
416 if (params%valid_path('case.fluid.flow_rate_force')) then
417 call neko_error("Flow rate forcing not available for adjoint_fluid_pnpn")
418 end if
419
420 ! Setup pressure solver
421 call neko_log%section("Pressure solver")
422
423 call json_get_or_default(params, &
424 'case.fluid.pressure_solver.max_iterations', &
425 solver_maxiter, 800)
426 call json_get(params, 'case.fluid.pressure_solver.type', solver_type)
427 call json_get(params, 'case.fluid.pressure_solver.preconditioner.type', &
428 precon_type)
429 call json_get(params, &
430 'case.fluid.pressure_solver.preconditioner', precon_params)
431 call json_get(params, 'case.fluid.pressure_solver.absolute_tolerance', &
432 abs_tol)
433 call json_get_or_default(params, 'case.fluid.pressure_solver.monitor', &
434 monitor, .false.)
435 call neko_log%message('Type : ('// trim(solver_type) // &
436 ', ' // trim(precon_type) // ')')
437 write(log_buf, '(A,ES13.6)') 'Abs tol :', abs_tol
438 call neko_log%message(log_buf)
439
440 call this%solver_factory(this%ksp_prs, this%dm_Xh%size(), &
441 solver_type, solver_maxiter, abs_tol, monitor)
442 call this%precon_factory_(this%pc_prs, this%ksp_prs, &
443 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, &
444 precon_type, precon_params)
445 call neko_log%end_section()
446
447 ! Initialize the advection factory
448 call neko_log%section("Advection factory")
449 call json_get_or_default(params, 'case.fluid.advection', advection, .true.)
450 call json_get(params, 'case.numerics', numerics_params)
451 call advection_adjoint_factory(this%adv, numerics_params, this%c_Xh, &
452 this%ulag, this%vlag, this%wlag, &
453 chkp%dtlag, chkp%tlag, this%ext_bdf, &
454 .not. advection)
455 ! Should be in init_base maybe?
456 this%chkp => chkp
457 ! This is probably scheme specific
458 ! Should not be init really, but more like, add fluid or something...
459 call this%chkp%add_fluid(this%u_adj, this%v_adj, this%w_adj, this%p_adj)
460
461 this%chkp%abx1 => this%abx1
462 this%chkp%abx2 => this%abx2
463 this%chkp%aby1 => this%aby1
464 this%chkp%aby2 => this%aby2
465 this%chkp%abz1 => this%abz1
466 this%chkp%abz2 => this%abz2
467 call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
468
469 call neko_log%end_section()
470
471 ! ------------------------------------------------------------------------ !
472 ! Handling the rescaling and baseflow
473
474 ! Read the norm scaling from the json file
475 call json_get_or_default(params, 'norm_scaling', &
476 this%norm_scaling, 0.5_rp)
477
478 ! The baseflow is the solution to the forward.
479 ! Userdefined baseflows can be invoked via setting initial conditions
480 ! call neko_registry%add_field(this%dm_Xh, 'u')
481 ! call neko_registry%add_field(this%dm_Xh, 'v')
482 ! call neko_registry%add_field(this%dm_Xh, 'w')
483 ! call neko_registry%add_field(this%dm_Xh, 'p')
484 this%u_b => neko_registry%get_field('u')
485 this%v_b => neko_registry%get_field('v')
486 this%w_b => neko_registry%get_field('w')
487 this%p_b => neko_registry%get_field('p')
488
489 ! Read the json file
490 call json_get_or_default(params, 'norm_target', &
491 this%norm_target, -1.0_rp)
492 call json_get_or_default(params, 'norm_tolerance', &
493 this%norm_tolerance, 10.0_rp)
494
495 ! Build the header
496 call json_get_or_default(params, 'output_file', &
497 file_name, 'power_iterations.csv')
498 call this%file_output%init(trim(file_name))
499 write(header_line, '(A)') 'Time, Norm, Scaling'
500 call this%file_output%set_header(header_line)
501
502 end subroutine adjoint_fluid_pnpn_init
503
504 subroutine adjoint_fluid_pnpn_restart(this, chkp)
505 class(adjoint_fluid_pnpn_t), target, intent(inout) :: this
506 type(chkp_t), intent(inout) :: chkp
507 real(kind=rp) :: dtlag(10), tlag(10)
508 integer :: i, n
509
510 dtlag = chkp%dtlag
511 tlag = chkp%tlag
512
513 n = this%u_adj%dof%size()
514 if (allocated(this%chkp%previous_mesh%elements) .or. &
515 chkp%previous_Xh%lx .ne. this%Xh%lx) then
516 associate(u => this%u_adj, v => this%v_adj, w => this%w_adj, &
517 p => this%p_adj, c_xh => this%c_Xh, ulag => this%ulag, &
518 vlag => this%vlag, wlag => this%wlag)
519 call col2(u%x, c_xh%mult, n)
520 call col2(v%x, c_xh%mult, n)
521 call col2(w%x, c_xh%mult, n)
522 call col2(p%x, c_xh%mult, n)
523 do i = 1, this%ulag%size()
524 call col2(ulag%lf(i)%x, c_xh%mult, n)
525 call col2(vlag%lf(i)%x, c_xh%mult, n)
526 call col2(wlag%lf(i)%x, c_xh%mult, n)
527 end do
528 end associate
529 end if
530
531 if (neko_bcknd_device .eq. 1) then
532 associate(u => this%u_adj, v => this%v_adj, w => this%w_adj, &
533 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
534 p => this%p_adj)
535 call device_memcpy(u%x, u%x_d, u%dof%size(), &
536 host_to_device, sync = .false.)
537 call device_memcpy(v%x, v%x_d, v%dof%size(), &
538 host_to_device, sync = .false.)
539 call device_memcpy(w%x, w%x_d, w%dof%size(), &
540 host_to_device, sync = .false.)
541 call device_memcpy(p%x, p%x_d, p%dof%size(), &
542 host_to_device, sync = .false.)
543 call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
544 u%dof%size(), host_to_device, sync = .false.)
545 call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
546 u%dof%size(), host_to_device, sync = .false.)
547
548 call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
549 v%dof%size(), host_to_device, sync = .false.)
550 call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
551 v%dof%size(), host_to_device, sync = .false.)
552
553 call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
554 w%dof%size(), host_to_device, sync = .false.)
555 call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
556 w%dof%size(), host_to_device, sync = .false.)
557 call device_memcpy(this%abx1%x, this%abx1%x_d, &
558 w%dof%size(), host_to_device, sync = .false.)
559 call device_memcpy(this%abx2%x, this%abx2%x_d, &
560 w%dof%size(), host_to_device, sync = .false.)
561 call device_memcpy(this%aby1%x, this%aby1%x_d, &
562 w%dof%size(), host_to_device, sync = .false.)
563 call device_memcpy(this%aby2%x, this%aby2%x_d, &
564 w%dof%size(), host_to_device, sync = .false.)
565 call device_memcpy(this%abz1%x, this%abz1%x_d, &
566 w%dof%size(), host_to_device, sync = .false.)
567 call device_memcpy(this%abz2%x, this%abz2%x_d, &
568 w%dof%size(), host_to_device, sync = .false.)
569 call device_memcpy(this%advx%x, this%advx%x_d, &
570 w%dof%size(), host_to_device, sync = .false.)
571 call device_memcpy(this%advy%x, this%advy%x_d, &
572 w%dof%size(), host_to_device, sync = .false.)
573 call device_memcpy(this%advz%x, this%advz%x_d, &
574 w%dof%size(), host_to_device, sync = .false.)
575 end associate
576 end if
577 ! Make sure that continuity is maintained (important for interpolation)
578 ! Do not do this for lagged rhs
579 ! (derivatives are not necessairly coninous across elements)
580
581 if (allocated(this%chkp%previous_mesh%elements) &
582 .or. this%chkp%previous_Xh%lx .ne. this%Xh%lx) then
583 call this%gs_Xh%op(this%u_adj, gs_op_add)
584 call this%gs_Xh%op(this%v_adj, gs_op_add)
585 call this%gs_Xh%op(this%w_adj, gs_op_add)
586 call this%gs_Xh%op(this%p_adj, gs_op_add)
587
588 do i = 1, this%ulag%size()
589 call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
590 call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
591 call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
592 end do
593 end if
594
595 end subroutine adjoint_fluid_pnpn_restart
596
597 subroutine adjoint_fluid_pnpn_free(this)
598 class(adjoint_fluid_pnpn_t), intent(inout) :: this
599
600 !Deallocate velocity and pressure fields
601 call this%scheme_free()
602
603 call this%bc_prs_surface%free()
604 call this%bc_sym_surface%free()
605 call this%bc_curl_curl%free()
606 call this%bclst_vel_res%free()
607 call this%bclst_dp%free()
608 call this%proj_prs%free()
609 call this%proj_vel%free()
610
611 call this%p_res%free()
612 call this%u_res%free()
613 call this%v_res%free()
614 call this%w_res%free()
615
616 call this%du%free()
617 call this%dv%free()
618 call this%dw%free()
619 call this%dp%free()
620
621 call this%abx1%free()
622 call this%aby1%free()
623 call this%abz1%free()
624
625 call this%abx2%free()
626 call this%aby2%free()
627 call this%abz2%free()
628
629 call this%advx%free()
630 call this%advy%free()
631 call this%advz%free()
632
633 if (allocated(this%Ax_vel)) then
634 deallocate(this%Ax_vel)
635 end if
636
637 if (allocated(this%Ax_prs)) then
638 deallocate(this%Ax_prs)
639 end if
640
641 if (allocated(this%prs_res)) then
642 deallocate(this%prs_res)
643 end if
644
645 if (allocated(this%vel_res)) then
646 deallocate(this%vel_res)
647 end if
648
649 if (allocated(this%sumab)) then
650 deallocate(this%sumab)
651 end if
652
653 if (allocated(this%makeabf)) then
654 deallocate(this%makeabf)
655 end if
656
657 if (allocated(this%makebdf)) then
658 deallocate(this%makebdf)
659 end if
660
661 if (allocated(this%makeoifs)) then
662 deallocate(this%makeoifs)
663 end if
664
665 if (allocated(this%ext_bdf)) then
666 deallocate(this%ext_bdf)
667 end if
668
669 ! call this%vol_flow%free()
670
671 end subroutine adjoint_fluid_pnpn_free
672
677 subroutine adjoint_fluid_pnpn_step(this, time, dt_controller)
678 class(adjoint_fluid_pnpn_t), target, intent(inout) :: this
679 type(time_state_t), intent(in) :: time
680 type(time_step_controller_t), intent(in) :: dt_controller
681 ! number of degrees of freedom
682 integer :: n
683 ! Solver results monitors (pressure + 3 velocity)
684 type(ksp_monitor_t) :: ksp_results(4)
685 type(field_t), pointer :: dx_p_adj, dy_p_adj, dz_p_adj, nx1, nx2, nx3, &
686 work1, work2
687 integer :: temp_indices(3)
688 integer :: cc_indices(8)
689 real(kind=rp) :: rho_val, mu_val
690
691 if (this%freeze) return
692
693 n = this%dm_Xh%size()
694
695 call profiler_start_region('Adjoint')
696 associate(u => this%u_adj, v => this%v_adj, w => this%w_adj, &
697 p => this%p_adj, &
698 u_e => this%u_adj_e, v_e => this%v_adj_e, w_e => this%w_adj_e, &
699 du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
700 u_b => this%u_b, v_b => this%v_b, w_b => this%w_b, &
701 u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
702 p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
703 xh => this%Xh, &
704 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
705 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
706 msh => this%msh, prs_res => this%prs_res, &
707 source_term => this%source_term, vel_res => this%vel_res, &
708 sumab => this%sumab, &
709 makeabf => this%makeabf, makebdf => this%makebdf, &
710 vel_projection_dim => this%vel_projection_dim, &
711 pr_projection_dim => this%pr_projection_dim, &
712 oifs => this%oifs, &
713 rho => this%rho, mu => this%mu, &
714 f_x => this%f_adj_x, f_y => this%f_adj_y, f_z => this%f_adj_z, &
715 t => time%t, tstep => time%tstep, dt => time%dt, &
716 ext_bdf => this%ext_bdf, event => glb_cmd_event)
717
718 ! Extrapolate the velocity if it's not done in nut_field estimation
719 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
720 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
721
722 ! Compute the source terms
723 call this%source_term%compute(time)
724
725 ! Add Neumann bc contributions to the RHS
726 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
727 this%dm_Xh%size(), time, strong = .false.)
728
729 if (oifs) then
730 call neko_error("OIFS not implemented for adjoint")
731
732 else
733 ! Add the advection operators to the right-hand-side.
734 call this%adv%compute_adjoint(u, v, w, u_b, v_b, w_b, &
735 f_x, f_y, f_z, &
736 xh, this%c_Xh, dm_xh%size())
737
738 ! At this point the RHS contains the sum of the advection operator and
739 ! additional source terms, evaluated using the velocity field from the
740 ! previous time-step. Now, this value is used in the explicit time
741 ! scheme to advance both terms in time.
742 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
743 this%abx2, this%aby2, this%abz2, &
744 f_x%x, f_y%x, f_z%x, &
745 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
746
747 ! Add the RHS contributions coming from the BDF scheme.
748 call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
749 u, v, w, c_xh%B, rho%x(1,1,1,1), dt, &
750 ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
751 end if
752
753 call ulag%update()
754 call vlag%update()
755 call wlag%update()
756
757 call this%bc_apply_vel(time, strong = .true.)
758 call this%bc_apply_prs(time)
759
760 ! Now we need the surface contribution of the curl curl BC.(explicit in p)
761 call neko_scratch_registry%request_field(dx_p_adj, cc_indices(1), .false.)
762 call neko_scratch_registry%request_field(dy_p_adj, cc_indices(2), .false.)
763 call neko_scratch_registry%request_field(dz_p_adj, cc_indices(3), .false.)
764
765 ! Note: zero interior
766 call neko_scratch_registry%request_field(nx1, cc_indices(4), .true.)
767 call neko_scratch_registry%request_field(nx2, cc_indices(5), .true.)
768 call neko_scratch_registry%request_field(nx3, cc_indices(6), .true.)
769
770 call neko_scratch_registry%request_field(work1, cc_indices(7), .false.)
771 call neko_scratch_registry%request_field(work2, cc_indices(8), .false.)
772
773 ! gradient of adjoint pressure (explicit)
774 call grad(dx_p_adj%x, dy_p_adj%x, dz_p_adj%x, this%p_adj%x, c_xh)
775
776 ! Now we compute the n x grad(p) (this include 2D weights)
777 call this%bc_curl_curl%apply_n_cross(nx1%x, nx2%x, nx3%x, dx_p_adj%x, &
778 dy_p_adj%x, dz_p_adj%x, dx_p_adj%size())
779
780 ! Now we need curl on the test function, note that transpose of curl is
781 ! negative curl
782 ! reuse dx_p_adj etc as fx, fy, fz etc
783 call curl(dx_p_adj, dy_p_adj, dz_p_adj, nx1, nx2, nx3, work1, work2, c_xh)
784
785 ! Forward does gsop on the residual (which has the pressure gradient)
786 call gs_xh%op(dx_p_adj, gs_op_add, event)
787 call device_event_sync(event)
788 call gs_xh%op(dy_p_adj, gs_op_add, event)
789 call device_event_sync(event)
790 call gs_xh%op(dz_p_adj, gs_op_add, event)
791 call device_event_sync(event)
792
793 ! multiplcity
794 if (neko_bcknd_device .eq. 1) then
795 call device_col2(dx_p_adj%x_d, c_xh%mult_d, dx_p_adj%size())
796 call device_col2(dy_p_adj%x_d, c_xh%mult_d, dx_p_adj%size())
797 call device_col2(dz_p_adj%x_d, c_xh%mult_d, dx_p_adj%size())
798 else
799 call col2(dx_p_adj%x, c_xh%mult, dx_p_adj%size())
800 call col2(dy_p_adj%x, c_xh%mult, dx_p_adj%size())
801 call col2(dz_p_adj%x, c_xh%mult, dx_p_adj%size())
802 end if
803
804 rho_val = rho%x(1,1,1,1)
805 mu_val = mu%x(1,1,1,1)
806 call field_add2s2(f_x, dx_p_adj, -mu_val / rho_val)
807 call field_add2s2(f_y, dy_p_adj, -mu_val / rho_val)
808 call field_add2s2(f_z, dz_p_adj, -mu_val / rho_val)
809
810 call neko_scratch_registry%relinquish_field(cc_indices)
811
812 ! Update material properties if necessary
813 call this%update_material_properties(time)
814
815 ! Compute intermediate velocity residual.
816
817 call profiler_start_region('Adjoint_velocity_residual')
818
819 call vel_res%compute(ax_vel, u, v, w, &
820 u_res, v_res, w_res, &
821 p, &
822 f_x, f_y, f_z, &
823 c_xh, msh, xh, &
824 mu, rho, ext_bdf%diffusion_coeffs(1), &
825 dt, dm_xh%size())
826
827 call gs_xh%op(u_res, gs_op_add, event)
828 call device_event_sync(event)
829 call gs_xh%op(v_res, gs_op_add, event)
830 call device_event_sync(event)
831 call gs_xh%op(w_res, gs_op_add, event)
832 call device_event_sync(event)
833
834 ! Set residual to zero at strong velocity boundaries.
835 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
836
837 call profiler_end_region('Adjoint_velocity_residual')
838
839 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
840 tstep, c_xh, n, dt_controller, 'Velocity')
841
842 call this%pc_vel%update()
843
844 call profiler_start_region("Adjoint_velocity_solve")
845 ksp_results(1:3) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
846 u_res%x, v_res%x, w_res%x, n, c_xh, &
847 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
848 this%ksp_vel%max_iter)
849
850 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
851 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
852 dt_controller)
853
854 if (neko_bcknd_device .eq. 1) then
855 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
856 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
857 else
858 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
859 end if
860 call profiler_end_region("Adjoint_velocity_solve")
861
862 !------------------------------------------------------------------------!
863 ! now the RHS of our pressure eqn is the new adjoint velocity
864 ! be careful with the order of the gsops here, we will handle this in the
865 ! residual calculation. So we enter WITHOUT a mass matrix.
866 call field_copy(f_x, u)
867 call field_copy(f_y, v)
868 call field_copy(f_z, w)
869 !------------------------------------------------------------------------!
870 call profiler_start_region('Adjoint_pressure_residual')
871
872 call prs_res%compute(p, p_res, &
873 u, v, w, &
874 f_x, f_y, f_z, &
875 c_xh, gs_xh, &
876 this%bc_prs_surface, this%bc_sym_surface, &
877 ax_prs, ext_bdf%diffusion_coeffs(1), dt, &
878 mu, rho, event)
879
880 ! De-mean the pressure residual when no strong pressure boundaries present
881 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1) then
882 call device_ortho(p_res%x_d, this%glb_n_points, n)
883 else if (.not. this%prs_dirichlet) then
884 call ortho(p_res%x, this%glb_n_points, n)
885 end if
886
887 call gs_xh%op(p_res, gs_op_add, event)
888 call device_event_sync(event)
889
890 ! Set the residual to zero at strong pressure boundaries.
891 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
892
893
894 call profiler_end_region('Adjoint_pressure_residual')
895
896
897 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
898 'Pressure')
899
900 call this%pc_prs%update()
901
902 call profiler_start_region('Adjoint_pressure_solve')
903
904 ! Solve for the pressure increment.
905 ksp_results(4) = &
906 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
907 this%bclst_dp, gs_xh)
908
909
910 call profiler_end_region('Adjoint_pressure_solve')
911
912 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
913 this%bclst_dp, gs_xh, n, tstep, dt_controller)
914
915 ! Update the pressure with the increment. Demean if necessary.
916 call field_add2(p, dp, n)
917 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1) then
918 call device_ortho(p%x_d, this%glb_n_points, n)
919 else if (.not. this%prs_dirichlet) then
920 call ortho(p%x, this%glb_n_points, n)
921 end if
922
923 ksp_results(4)%name = 'Adjoint Pressure'
924 ksp_results(1)%name = 'Adjoint Velocity U'
925 ksp_results(2)%name = 'Adjoint Velocity V'
926 ksp_results(3)%name = 'Adjoint Velocity W'
927
928 if (this%forced_flow_rate) then
929 call neko_error('Forced flow rate is not implemented for the adjoint')
930
931 end if
932
933 !------------------------------------------------------------------------!
934 ! correct the velocity with the pressure
935 call neko_scratch_registry%request_field(dx_p_adj, temp_indices(1), &
936 .false.)
937 call neko_scratch_registry%request_field(dy_p_adj, temp_indices(2), &
938 .false.)
939 call neko_scratch_registry%request_field(dz_p_adj, temp_indices(3), &
940 .false.)
941
942 ! gradient of adjoint pressure (explicit)
943 call opgrad(dx_p_adj%x, dy_p_adj%x, dz_p_adj%x, this%p_adj%x, c_xh)
944
945 ! they gsop the residual (which has the pressure gradient)
946 call gs_xh%op(dx_p_adj, gs_op_add, event)
947 call device_event_sync(event)
948 call gs_xh%op(dy_p_adj, gs_op_add, event)
949 call device_event_sync(event)
950 call gs_xh%op(dz_p_adj, gs_op_add, event)
951 call device_event_sync(event)
952
953 ! divide by mass matrix
954 if (neko_bcknd_device .eq. 1) then
955 call device_col2(dx_p_adj%x_d, c_xh%Binv_d, dx_p_adj%size())
956 call device_col2(dy_p_adj%x_d, c_xh%Binv_d, dx_p_adj%size())
957 call device_col2(dz_p_adj%x_d, c_xh%Binv_d, dx_p_adj%size())
958 else
959 ! NOTE. This term comes from the handling of the pressure RHS, which
960 ! DOES include the multiplicity in the op.
961 call col2(dx_p_adj%x, c_xh%Binv, dx_p_adj%size())
962 call col2(dy_p_adj%x, c_xh%Binv, dx_p_adj%size())
963 call col2(dz_p_adj%x, c_xh%Binv, dx_p_adj%size())
964 end if
965
966 if (neko_bcknd_device .eq. 1) then
967 call device_opadd2cm(u%x_d, v%x_d, w%x_d, dx_p_adj%x_d, &
968 dy_p_adj%x_d, dz_p_adj%x_d, -1.0_rp, n, msh%gdim)
969 else
970 call opadd2cm(u%x, v%x, w%x, dx_p_adj%x, dy_p_adj%x, dz_p_adj%x, &
971 -1.0_rp, n, msh%gdim)
972 end if
973
974 call neko_scratch_registry%relinquish_field(temp_indices)
975 !------------------------------------------------------------------------!
976
977 call fluid_step_info(time, ksp_results, &
978 this%full_stress_formulation, this%strict_convergence)
979
980 end associate
981 call profiler_end_region('Adjoint')
982
983 end subroutine adjoint_fluid_pnpn_step
984
989 subroutine adjoint_fluid_pnpn_setup_bcs(this, user, params)
990 use mpi_f08, only: mpi_in_place
991 class(adjoint_fluid_pnpn_t), intent(inout) :: this
992 type(user_t), target, intent(in) :: user
993 type(json_file), intent(inout) :: params
994 integer :: i, n_bcs, j, zone_size, global_zone_size, ierr
995 class(bc_t), pointer :: bc_i
996 type(json_core) :: core
997 type(json_value), pointer :: bc_object
998 type(json_file) :: bc_subdict
999 logical :: found
1000 ! Monitor which boundary zones have been marked
1001 logical, allocatable :: marked_zones(:)
1002 integer, allocatable :: zone_indices(:)
1003 character(len=:), allocatable :: json_key
1004
1005 ! Lists for the residuals and solution increments
1006 call this%bclst_vel_res%init()
1007 call this%bclst_du%init()
1008 call this%bclst_dv%init()
1009 call this%bclst_dw%init()
1010 call this%bclst_dp%init()
1011
1012 call this%bc_vel_res%init_from_components(this%c_Xh)
1013 call this%bc_du%init_from_components(this%c_Xh)
1014 call this%bc_dv%init_from_components(this%c_Xh)
1015 call this%bc_dw%init_from_components(this%c_Xh)
1016 call this%bc_dp%init_from_components(this%c_Xh)
1017
1018 ! Special PnPn boundary conditions for pressure
1019 call this%bc_prs_surface%init_from_components(this%c_Xh)
1020 call this%bc_sym_surface%init_from_components(this%c_Xh)
1021 call this%bc_curl_curl%init_from_components(this%c_Xh)
1022
1023 json_key = 'case.adjoint_fluid.boundary_conditions'
1024
1025 ! Populate bcs_vel and bcs_prs based on the case file
1026 if (params%valid_path(json_key)) then
1027 call params%info(json_key, n_children = n_bcs)
1028 call params%get_core(core)
1029 call params%get(json_key, bc_object, found)
1030
1031 !
1032 ! Velocity bcs
1033 !
1034 call this%bcs_vel%init(n_bcs)
1035
1036 allocate(marked_zones(size(this%msh%labeled_zones)))
1037 marked_zones = .false.
1038
1039 do i = 1, n_bcs
1040 ! Create a new json containing just the subdict for this bc
1041 call json_extract_item(core, bc_object, i, bc_subdict)
1042
1043 call json_get(bc_subdict, "zone_indices", zone_indices)
1044
1045 ! Check that we are not trying to assing a bc to zone, for which one
1046 ! has already been assigned and that the zone has more than 0 size
1047 ! in the mesh.
1048 do j = 1, size(zone_indices)
1049 zone_size = this%msh%labeled_zones(zone_indices(j))%size
1050 call mpi_allreduce(zone_size, global_zone_size, 1, &
1051 mpi_integer, mpi_max, neko_comm, ierr)
1052
1053 if (global_zone_size .eq. 0) then
1054 write(error_unit, '(A, A, I0, A, A, I0, A)') "*** ERROR ***: ",&
1055 "Zone index ", zone_indices(j), &
1056 " is invalid as this zone has 0 size, meaning it does ", &
1057 "not in the mesh. Check adjoint boundary condition ", &
1058 i, "."
1059 error stop
1060 end if
1061
1062 if (marked_zones(zone_indices(j)) .eqv. .true.) then
1063 write(error_unit, '(A, A, I0, A, A, A, A)') "*** ERROR ***: ", &
1064 "Zone with index ", zone_indices(j), &
1065 " has already been assigned a boundary condition. ", &
1066 "Please check your boundary_conditions entry for the ", &
1067 "adjoint and make sure that each zone index appears ", &
1068 "only in a single boundary condition."
1069 error stop
1070 else
1071 marked_zones(zone_indices(j)) = .true.
1072 end if
1073 end do
1074
1075 bc_i => null()
1076 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
1077
1078 ! Not all bcs require an allocation for velocity in particular,
1079 ! so we check.
1080 if (associated(bc_i)) then
1081
1082 ! We need to treat mixed bcs separately because they are by
1083 ! convention marked weak and currently contain nested
1084 ! bcs, some of which are strong.
1085 select type (bc_i)
1086 type is (symmetry_t)
1087 ! Symmetry has 3 internal bcs, but only one actually contains
1088 ! markings.
1089 ! Symmetry's apply_scalar doesn't do anything, so we need to
1090 ! mark individual nested bcs to the du,dv,dw, whereas the
1091 ! vel_res can just get symmetry as a whole, because on this
1092 ! list we call apply_vector.
1093 ! Additionally we have to mark the special surface bc for p.
1094 call this%bclst_vel_res%append(bc_i)
1095 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1096 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1097 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1098
1099 call this%bcs_vel%append(bc_i)
1100
1101 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
1102 type is (non_normal_t)
1103 ! This is a bc for the residuals and increments, not the
1104 ! velocity itself. So, don't append to bcs_vel
1105 call this%bclst_vel_res%append(bc_i)
1106 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1107 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1108 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1109 type is (shear_stress_t)
1110 ! Same as symmetry
1111 call this%bclst_vel_res%append(bc_i%symmetry)
1112 call this%bclst_du%append(bc_i%symmetry%bc_x)
1113 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1114 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1115
1116 call this%bcs_vel%append(bc_i)
1117 type is (wall_model_bc_t)
1118 ! Same as symmetry
1119 call this%bclst_vel_res%append(bc_i%symmetry)
1120 call this%bclst_du%append(bc_i%symmetry%bc_x)
1121 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1122 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1123
1124 call this%bcs_vel%append(bc_i)
1125 class default
1126
1127 ! For the default case we use our dummy zero_dirichlet bcs to
1128 ! mark the same faces as in ordinary velocity dirichlet
1129 ! conditions.
1130 ! Additionally we mark the special PnPn pressure bc.
1131 if (bc_i%strong .eqv. .true.) then
1132 call this%bc_vel_res%mark_facets(bc_i%marked_facet)
1133 call this%bc_du%mark_facets(bc_i%marked_facet)
1134 call this%bc_dv%mark_facets(bc_i%marked_facet)
1135 call this%bc_dw%mark_facets(bc_i%marked_facet)
1136
1137 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
1138 end if
1139
1140 ! add all BCs to curl curl
1141 call this%bc_curl_curl%mark_facets(bc_i%marked_facet)
1142
1143 call this%bcs_vel%append(bc_i)
1144 end select
1145 end if
1146 end do
1147
1148 ! Make sure all labeled zones with non-zero size have been marked
1149 do i = 1, size(this%msh%labeled_zones)
1150 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1151 (marked_zones(i) .eqv. .false.)) then
1152 write(error_unit, '(A, A, I0)') "*** ERROR ***: ", &
1153 "No adjoint boundary condition assigned to zone ", i
1154 error stop
1155 end if
1156 end do
1157
1158 !
1159 ! Pressure bcs
1160 !
1161 call this%bcs_prs%init(n_bcs)
1162
1163 do i = 1, n_bcs
1164 ! Create a new json containing just the subdict for this bc
1165 call json_extract_item(core, bc_object, i, bc_subdict)
1166 bc_i => null()
1167 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
1168
1169 ! Not all bcs require an allocation for pressure in particular,
1170 ! so we check.
1171 if (associated(bc_i)) then
1172 call this%bcs_prs%append(bc_i)
1173
1174 ! Mark strong bcs in the dummy dp bc to force zero change.
1175 if (bc_i%strong .eqv. .true.) then
1176 call this%bc_dp%mark_facets(bc_i%marked_facet)
1177 end if
1178
1179 end if
1180
1181 end do
1182 else
1183 ! Check that there are no labeled zones, i.e. all are periodic.
1184 do i = 1, size(this%msh%labeled_zones)
1185 if (this%msh%labeled_zones(i)%size .gt. 0) then
1186 call neko_error("No boundary_conditions entry in the case file!")
1187 end if
1188 end do
1189
1190 end if
1191
1192 call this%bc_prs_surface%finalize()
1193 call this%bc_sym_surface%finalize()
1194 call this%bc_curl_curl%finalize()
1195
1196 call this%bc_vel_res%finalize()
1197 call this%bc_du%finalize()
1198 call this%bc_dv%finalize()
1199 call this%bc_dw%finalize()
1200 call this%bc_dp%finalize()
1201
1202 call this%bclst_vel_res%append(this%bc_vel_res)
1203 call this%bclst_du%append(this%bc_du)
1204 call this%bclst_dv%append(this%bc_dv)
1205 call this%bclst_dw%append(this%bc_dw)
1206 call this%bclst_dp%append(this%bc_dp)
1207
1208 ! If we have no strong pressure bcs, we will demean the pressure
1209 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1210 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1211 mpi_logical, mpi_lor, neko_comm)
1212
1213 end subroutine adjoint_fluid_pnpn_setup_bcs
1214
1216 subroutine adjoint_fluid_pnpn_write_boundary_conditions(this)
1217 class(adjoint_fluid_pnpn_t), target, intent(inout) :: this
1218 type(dirichlet_t) :: bdry_mask
1219 type(field_t), pointer :: bdry_field
1220 type(file_t) :: bdry_file
1221 integer :: temp_index, i
1222 class(bc_t), pointer :: bci
1223 character(len=LOG_SIZE) :: log_buf
1224
1225 call neko_log%section("Adjoint boundary conditions")
1226 write(log_buf, '(A)') &
1227 'Marking using integer keys in boundary_adjoint0.f00000'
1228 call neko_log%message(log_buf)
1229 write(log_buf, '(A)') 'Condition-value pairs: '
1230 call neko_log%message(log_buf)
1231 write(log_buf, '(A)') ' no_slip = 1'
1232 call neko_log%message(log_buf)
1233 write(log_buf, '(A)') ' velocity_value = 2'
1234 call neko_log%message(log_buf)
1235 write(log_buf, '(A)') ' outflow, normal_outflow (+dong) = 3'
1236 call neko_log%message(log_buf)
1237 write(log_buf, '(A)') ' symmetry = 4'
1238 call neko_log%message(log_buf)
1239 write(log_buf, '(A)') ' user_velocity_pointwise = 5'
1240 call neko_log%message(log_buf)
1241 write(log_buf, '(A)') ' periodic = 6'
1242 call neko_log%message(log_buf)
1243 write(log_buf, '(A)') ' user_velocity = 7'
1244 call neko_log%message(log_buf)
1245 write(log_buf, '(A)') ' user_pressure = 8'
1246 call neko_log%message(log_buf)
1247 write(log_buf, '(A)') ' shear_stress = 9'
1248 call neko_log%message(log_buf)
1249 write(log_buf, '(A)') ' wall_modelling = 10'
1250 call neko_log%message(log_buf)
1251 write(log_buf, '(A)') ' blasius_profile = 11'
1252 call neko_log%message(log_buf)
1253 call neko_log%end_section()
1254
1255 call neko_scratch_registry%request_field(bdry_field, temp_index, .true.)
1256
1257
1258
1259 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1260 call bdry_mask%mark_zone(this%msh%periodic)
1261 call bdry_mask%finalize()
1262 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1263 call bdry_mask%free()
1264
1265 do i = 1, this%bcs_prs%size()
1266 bci => this%bcs_prs%get(i)
1267 select type (bc => bci)
1268 type is (zero_dirichlet_t)
1269 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1270 call bdry_mask%mark_facets(bci%marked_facet)
1271 call bdry_mask%finalize()
1272 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1273 call bdry_mask%free()
1274 type is (dong_outflow_t)
1275 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1276 call bdry_mask%mark_facets(bci%marked_facet)
1277 call bdry_mask%finalize()
1278 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1279 call bdry_mask%free()
1280 type is (field_dirichlet_t)
1281 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1282 call bdry_mask%mark_facets(bci%marked_facet)
1283 call bdry_mask%finalize()
1284 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1285 call bdry_mask%free()
1286 end select
1287 end do
1288
1289 do i = 1, this%bcs_vel%size()
1290 bci => this%bcs_vel%get(i)
1291 select type (bc => bci)
1292 type is (zero_dirichlet_t)
1293 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1294 call bdry_mask%mark_facets(bci%marked_facet)
1295 call bdry_mask%finalize()
1296 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1297 call bdry_mask%free()
1298 type is (inflow_t)
1299 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1300 call bdry_mask%mark_facets(bci%marked_facet)
1301 call bdry_mask%finalize()
1302 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1303 call bdry_mask%free()
1304 type is (symmetry_t)
1305 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1306 call bdry_mask%mark_facets(bci%marked_facet)
1307 call bdry_mask%finalize()
1308 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1309 call bdry_mask%free()
1310 type is (field_dirichlet_vector_t)
1311 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1312 call bdry_mask%mark_facets(bci%marked_facet)
1313 call bdry_mask%finalize()
1314 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1315 call bdry_mask%free()
1316 type is (shear_stress_t)
1317 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1318 call bdry_mask%mark_facets(bci%marked_facet)
1319 call bdry_mask%finalize()
1320 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1321 call bdry_mask%free()
1322 type is (wall_model_bc_t)
1323 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1324 call bdry_mask%mark_facets(bci%marked_facet)
1325 call bdry_mask%finalize()
1326 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1327 call bdry_mask%free()
1328 type is (blasius_t)
1329 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1330 call bdry_mask%mark_facets(bci%marked_facet)
1331 call bdry_mask%finalize()
1332 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1333 call bdry_mask%free()
1334 end select
1335 end do
1336
1337
1338 call bdry_file%init('boundary_adjoint.fld')
1339 call bdry_file%write(bdry_field)
1340
1341 call neko_scratch_registry%relinquish_field(temp_index)
1342 end subroutine adjoint_fluid_pnpn_write_boundary_conditions
1343
1344 ! End of section to verify
1345 ! ========================================================================== !
1346
1347 subroutine rescale_fluid(fluid_data, scale)
1348
1349 class(adjoint_fluid_pnpn_t), intent(inout) :: fluid_data
1351 real(kind=rp), intent(in) :: scale
1352
1353 ! Local variables
1354 integer :: i
1355
1356 ! Scale the velocity fields
1357 if (neko_bcknd_device .eq. 1) then
1358 call device_cmult(fluid_data%u_adj%x_d, scale, fluid_data%u_adj%size())
1359 call device_cmult(fluid_data%v_adj%x_d, scale, fluid_data%v_adj%size())
1360 call device_cmult(fluid_data%w_adj%x_d, scale, fluid_data%w_adj%size())
1361 else
1362 call cmult(fluid_data%u_adj%x, scale, fluid_data%u_adj%size())
1363 call cmult(fluid_data%v_adj%x, scale, fluid_data%v_adj%size())
1364 call cmult(fluid_data%w_adj%x, scale, fluid_data%w_adj%size())
1365 end if
1366
1367 ! Scale the right hand sides
1368 if (neko_bcknd_device .eq. 1) then
1369 call device_cmult(fluid_data%f_adj_x%x_d, scale, &
1370 fluid_data%f_adj_x%size())
1371 call device_cmult(fluid_data%f_adj_y%x_d, scale, &
1372 fluid_data%f_adj_y%size())
1373 call device_cmult(fluid_data%f_adj_z%x_d, scale, &
1374 fluid_data%f_adj_z%size())
1375 ! HARRY
1376 ! maybe the abx's too
1377 call device_cmult(fluid_data%abx1%x_d, scale, fluid_data%abx1%size())
1378 call device_cmult(fluid_data%aby1%x_d, scale, fluid_data%aby1%size())
1379 call device_cmult(fluid_data%abz1%x_d, scale, fluid_data%abz1%size())
1380 call device_cmult(fluid_data%abx2%x_d, scale, fluid_data%abx2%size())
1381 call device_cmult(fluid_data%aby2%x_d, scale, fluid_data%aby2%size())
1382 call device_cmult(fluid_data%abz2%x_d, scale, fluid_data%abz2%size())
1383
1384 else
1385 call cmult(fluid_data%f_adj_x%x, scale, fluid_data%f_adj_x%size())
1386 call cmult(fluid_data%f_adj_y%x, scale, fluid_data%f_adj_y%size())
1387 call cmult(fluid_data%f_adj_z%x, scale, fluid_data%f_adj_z%size())
1388
1389 call cmult(fluid_data%abx1%x, scale, fluid_data%abx1%size())
1390 call cmult(fluid_data%aby1%x, scale, fluid_data%aby1%size())
1391 call cmult(fluid_data%abz1%x, scale, fluid_data%abz1%size())
1392
1393 call cmult(fluid_data%abx2%x, scale, fluid_data%abx2%size())
1394 call cmult(fluid_data%aby2%x, scale, fluid_data%aby2%size())
1395 call cmult(fluid_data%abz2%x, scale, fluid_data%abz2%size())
1396 end if
1397
1398 ! Scale the lag terms
1399 if (neko_bcknd_device .eq. 1) then
1400 do i = 1, fluid_data%ulag%size()
1401 call device_cmult(fluid_data%ulag%lf(i)%x_d, &
1402 scale, fluid_data%ulag%lf(i)%size())
1403 end do
1404
1405 do i = 1, fluid_data%vlag%size()
1406 call device_cmult(fluid_data%vlag%lf(i)%x_d, &
1407 scale, fluid_data%vlag%lf(i)%size())
1408 end do
1409
1410 do i = 1, fluid_data%wlag%size()
1411 call device_cmult(fluid_data%wlag%lf(i)%x_d, &
1412 scale, fluid_data%wlag%lf(i)%size())
1413 end do
1414 else
1415 do i = 1, fluid_data%ulag%size()
1416 call cmult(fluid_data%ulag%lf(i)%x, &
1417 scale, fluid_data%ulag%lf(i)%size())
1418 end do
1419
1420 do i = 1, fluid_data%vlag%size()
1421 call cmult(fluid_data%vlag%lf(i)%x, &
1422 scale, fluid_data%vlag%lf(i)%size())
1423 end do
1424
1425 do i = 1, fluid_data%wlag%size()
1426 call cmult(fluid_data%wlag%lf(i)%x, &
1427 scale, fluid_data%wlag%lf(i)%size())
1428 end do
1429 end if
1430
1431 end subroutine rescale_fluid
1432
1433 function norm(x, y, z, B, volume, n)
1434 use mpi_f08, only: mpi_in_place
1435
1436 integer, intent(in) :: n
1437 real(kind=rp), dimension(n), intent(in) :: x, y, z
1438 real(kind=rp), dimension(n), intent(in) :: b
1439 real(kind=rp), intent(in) :: volume
1440
1441 real(kind=rp) :: norm
1442
1443 norm = vlsc3(x, x, b, n) + vlsc3(y, y, b, n) + vlsc3(z, z, b, n)
1444
1445 call mpi_allreduce(mpi_in_place, norm, 1, &
1446 mpi_real_precision, mpi_sum, mpi_comm_world)
1447
1448 norm = sqrt(norm / volume)
1449 end function norm
1450
1451 function device_norm(x_d, y_d, z_d, B_d, volume, n)
1452 use mpi_f08, only: mpi_in_place
1453
1454 type(c_ptr), intent(in) :: x_d, y_d, z_d
1455 type(c_ptr), intent(in) :: B_d
1456 real(kind=rp), intent(in) :: volume
1457 integer, intent(in) :: n
1458
1459 real(kind=rp) :: device_norm
1460
1461 device_norm = device_vlsc3(x_d, x_d, b_d, n) + &
1462 device_vlsc3(y_d, y_d, b_d, n) + &
1463 device_vlsc3(z_d, z_d, b_d, n)
1464
1465 call mpi_allreduce(mpi_in_place, device_norm, 1, &
1466 mpi_real_precision, mpi_sum, mpi_comm_world)
1467
1468 device_norm = sqrt(device_norm / volume)
1469
1470 end function device_norm
1471
1476 subroutine power_iterations_compute(this, t, tstep)
1477 class(adjoint_fluid_pnpn_t), target, intent(inout) :: this
1478 real(kind=rp), intent(in) :: t
1479 integer, intent(in) :: tstep
1480
1481 ! Local variables
1482 real(kind=rp) :: scaling_factor
1483 real(kind=rp) :: norm_l2, norm_l2_base
1484 character(len=256) :: log_message
1485 type(vector_t) :: data_line
1486 integer :: n
1487
1488 n = this%c_Xh%dof%size()
1489 if (tstep .eq. 1) then
1490 if (neko_bcknd_device .eq. 1) then
1491 norm_l2_base = device_norm(this%u_adj%x_d, this%v_adj%x_d, &
1492 this%w_adj%x_d, &
1493 this%c_Xh%B_d, this%c_Xh%volume, n)
1494 else
1495 norm_l2_base = this%norm_scaling * norm(this%u_adj%x, this%v_adj%x, &
1496 this%w_adj%x, &
1497 this%c_Xh%B, this%c_Xh%volume, n)
1498 end if
1499 if (this%norm_target .lt. 0.0_rp) then
1500 this%norm_target = norm_l2_base
1501 end if
1502
1503 this%norm_l2_upper = this%norm_tolerance * this%norm_target
1504 this%norm_l2_lower = this%norm_target / this%norm_tolerance
1505
1506 end if
1507
1508 ! Compute the norm of the velocity field and eigenvalue estimate
1509 if (neko_bcknd_device .eq. 1) then
1510 norm_l2 = device_norm(this%u_adj%x_d, this%v_adj%x_d, this%w_adj%x_d, &
1511 this%c_Xh%B_d, this%c_Xh%volume, n)
1512 else
1513 norm_l2 = norm(this%u_adj%x, this%v_adj%x, this%w_adj%x, &
1514 this%c_Xh%B, this%c_Xh%volume, n)
1515 end if
1516 norm_l2 = sqrt(this%norm_scaling) * norm_l2
1517 scaling_factor = 1.0_rp
1518
1519 ! Rescale the flow if necessary
1520 if (norm_l2 .gt. this%norm_l2_upper &
1521 .or. norm_l2 .lt. this%norm_l2_lower) then
1522 scaling_factor = this%norm_target / norm_l2
1523 call rescale_fluid(this, scaling_factor)
1524 norm_l2 = this%norm_target
1525
1526 if (tstep .eq. 1) then
1527 scaling_factor = 1.0_rp
1528 end if
1529 end if
1530
1531 ! Log the results
1532 !call neko_log%section('Power Iterations', lvl = NEKO_LOG_DEBUG)
1533 call neko_log%section('Power Iterations')
1534
1535 write (log_message, '(A7,E20.14)') 'Norm: ', norm_l2
1536 call neko_log%message(log_message, lvl = neko_log_debug)
1537 write (log_message, '(A7,E20.14)') 'Scaling: ', scaling_factor
1538 call neko_log%message(log_message, lvl = neko_log_debug)
1539
1540 ! Save to file
1541 call data_line%init(2)
1542 data_line%x = [norm_l2, scaling_factor]
1543 call this%file_output%write(data_line, t)
1544
1545 !call neko_log%end_section('Power Iterations', lvl = NEKO_LOG_DEBUG)
1546 call neko_log%end_section('Power Iterations')
1547 end subroutine power_iterations_compute
1548
1549end module adjoint_fluid_pnpn
Boundary condition factory for pressure.
Adjoint Pn/Pn formulation.
Subroutines to add advection terms to the RHS of a transport equation.
Abstract type to compute pressure residual.
Abstract type to compute velocity residual.
Base abstract type for computing the advection operator.
Dirichlet condition in facet normal direction.