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
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, &
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, &
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_integer, &
99 use operators,
only : opgrad, curl, grad
109 type(field_t) :: p_res, u_res, v_res, w_res
113 type(field_t) :: dp, du, dv, dw
120 class(ax_t),
allocatable :: ax_vel
122 class(ax_t),
allocatable :: ax_prs
129 type(projection_t) :: proj_prs
130 type(projection_vel_t) :: proj_vel
137 type(facet_normal_t) :: bc_prs_surface
140 type(facet_normal_t) :: bc_sym_surface
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
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
170 logical :: prs_dirichlet = .false.
180 type(field_t) :: abx1, aby1, abz1
181 type(field_t) :: abx2, aby2, abz2
184 type(field_t) :: advx, advy, advz
193 class(rhs_maker_sumab_t),
allocatable :: sumab
196 class(rhs_maker_ext_t),
allocatable :: makeabf
199 class(rhs_maker_bdf_t),
allocatable :: makebdf
202 class(rhs_maker_oifs_t),
allocatable :: makeoifs
208 logical :: full_stress_formulation = .false.
213 real(kind=rp) :: norm_scaling
214 real(kind=rp) :: norm_target
215 real(kind=rp) :: norm_tolerance
221 real(kind=rp) :: norm_l2_base
224 real(kind=rp) :: norm_l2_upper
226 real(kind=rp) :: norm_l2_lower
229 type(file_t) :: file_output
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
247 procedure,
public, pass(this) :: pw_compute_ => power_iterations_compute
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
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
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
301 type(json_file) :: numerics_params, precon_params
304 character(len=:),
allocatable :: file_name
305 character(len=256) :: header_line
310 call this%init_base(msh, lx, params, scheme, user, .true.)
314 call neko_registry%add_field(this%dm_Xh,
'p_adj')
315 this%p_adj => neko_registry%get_field(
'p_adj')
321 call json_get(params,
'case.numerics.time_order', integer_val)
322 allocate(this%ext_bdf)
323 call this%ext_bdf%init(integer_val)
325 call json_get_or_default(params,
"case.fluid.full_stress_formulation", &
326 this%full_stress_formulation, .false.)
328 if (this%full_stress_formulation .eqv. .true.)
then
330 "Full stress formulation is not supported in the adjoint module.")
341 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
344 call adjoint_pnpn_prs_res_factory(this%prs_res)
347 call adjoint_pnpn_vel_res_factory(this%vel_res)
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 " // &
356 call json_get(params,
'case.fluid.nut_field', this%nut_field_name)
358 this%nut_field_name =
""
362 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
366 call rhs_maker_sumab_fctry(this%sumab)
369 call rhs_maker_ext_fctry(this%makeabf)
372 call rhs_maker_bdf_fctry(this%makebdf)
375 call rhs_maker_oifs_fctry(this%makeoifs)
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)
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")
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')
402 call this%setup_bcs(user, params)
405 call json_get_or_default(params,
'case.output_boundary', found, .false.)
406 if (found)
call this%write_boundary_conditions()
408 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
409 this%pr_projection_activ_step)
411 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
412 this%vel_projection_activ_step)
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")
421 call neko_log%section(
"Pressure solver")
423 call json_get_or_default(params, &
424 'case.fluid.pressure_solver.max_iterations', &
426 call json_get(params,
'case.fluid.pressure_solver.type', solver_type)
427 call json_get(params,
'case.fluid.pressure_solver.preconditioner.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', &
433 call json_get_or_default(params,
'case.fluid.pressure_solver.monitor', &
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)
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()
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, &
459 call this%chkp%add_fluid(this%u_adj, this%v_adj, this%w_adj, this%p_adj)
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)
469 call neko_log%end_section()
475 call json_get_or_default(params,
'norm_scaling', &
476 this%norm_scaling, 0.5_rp)
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')
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)
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)
502 end subroutine adjoint_fluid_pnpn_init
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)
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)
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,&
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.)
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.)
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.)
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)
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)
595 end subroutine adjoint_fluid_pnpn_restart
597 subroutine adjoint_fluid_pnpn_free(this)
598 class(adjoint_fluid_pnpn_t),
intent(inout) :: this
601 call this%scheme_free()
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()
611 call this%p_res%free()
612 call this%u_res%free()
613 call this%v_res%free()
614 call this%w_res%free()
621 call this%abx1%free()
622 call this%aby1%free()
623 call this%abz1%free()
625 call this%abx2%free()
626 call this%aby2%free()
627 call this%abz2%free()
629 call this%advx%free()
630 call this%advy%free()
631 call this%advz%free()
633 if (
allocated(this%Ax_vel))
then
634 deallocate(this%Ax_vel)
637 if (
allocated(this%Ax_prs))
then
638 deallocate(this%Ax_prs)
641 if (
allocated(this%prs_res))
then
642 deallocate(this%prs_res)
645 if (
allocated(this%vel_res))
then
646 deallocate(this%vel_res)
649 if (
allocated(this%sumab))
then
650 deallocate(this%sumab)
653 if (
allocated(this%makeabf))
then
654 deallocate(this%makeabf)
657 if (
allocated(this%makebdf))
then
658 deallocate(this%makebdf)
661 if (
allocated(this%makeoifs))
then
662 deallocate(this%makeoifs)
665 if (
allocated(this%ext_bdf))
then
666 deallocate(this%ext_bdf)
671 end subroutine adjoint_fluid_pnpn_free
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
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, &
687 integer :: temp_indices(3)
688 integer :: cc_indices(8)
689 real(kind=rp) :: rho_val, mu_val
691 if (this%freeze)
return
693 n = this%dm_Xh%size()
695 call profiler_start_region(
'Adjoint')
696 associate(u => this%u_adj, v => this%v_adj, w => this%w_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, &
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, &
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)
719 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
720 ulag, vlag, wlag, ext_bdf%advection_coeffs%x, ext_bdf%nadv)
723 call this%source_term%compute(time)
726 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
727 this%dm_Xh%size(), time, strong = .false.)
730 call neko_error(
"OIFS not implemented for adjoint")
734 call this%adv%compute_adjoint(u, v, w, u_b, v_b, w_b, &
736 xh, this%c_Xh, dm_xh%size())
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%x, n)
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%x, ext_bdf%ndiff, n, &
751 c_xh%Blag, c_xh%Blaglag)
758 call this%bc_apply_vel(time, strong = .true.)
759 call this%bc_apply_prs(time)
762 call neko_scratch_registry%request_field(dx_p_adj, cc_indices(1), .false.)
763 call neko_scratch_registry%request_field(dy_p_adj, cc_indices(2), .false.)
764 call neko_scratch_registry%request_field(dz_p_adj, cc_indices(3), .false.)
767 call neko_scratch_registry%request_field(nx1, cc_indices(4), .true.)
768 call neko_scratch_registry%request_field(nx2, cc_indices(5), .true.)
769 call neko_scratch_registry%request_field(nx3, cc_indices(6), .true.)
771 call neko_scratch_registry%request_field(work1, cc_indices(7), .false.)
772 call neko_scratch_registry%request_field(work2, cc_indices(8), .false.)
775 call grad(dx_p_adj%x, dy_p_adj%x, dz_p_adj%x, this%p_adj%x, c_xh)
778 call this%bc_curl_curl%apply_n_cross(nx1%x, nx2%x, nx3%x, dx_p_adj%x, &
779 dy_p_adj%x, dz_p_adj%x, dx_p_adj%size())
784 call curl(dx_p_adj, dy_p_adj, dz_p_adj, nx1, nx2, nx3, work1, work2, c_xh)
787 call gs_xh%op(dx_p_adj, gs_op_add, event)
788 call device_event_sync(event)
789 call gs_xh%op(dy_p_adj, gs_op_add, event)
790 call device_event_sync(event)
791 call gs_xh%op(dz_p_adj, gs_op_add, event)
792 call device_event_sync(event)
795 if (neko_bcknd_device .eq. 1)
then
796 call device_col2(dx_p_adj%x_d, c_xh%mult_d, dx_p_adj%size())
797 call device_col2(dy_p_adj%x_d, c_xh%mult_d, dx_p_adj%size())
798 call device_col2(dz_p_adj%x_d, c_xh%mult_d, dx_p_adj%size())
800 call col2(dx_p_adj%x, c_xh%mult, dx_p_adj%size())
801 call col2(dy_p_adj%x, c_xh%mult, dx_p_adj%size())
802 call col2(dz_p_adj%x, c_xh%mult, dx_p_adj%size())
805 rho_val = rho%x(1,1,1,1)
806 mu_val = mu%x(1,1,1,1)
807 call field_add2s2(f_x, dx_p_adj, -mu_val / rho_val)
808 call field_add2s2(f_y, dy_p_adj, -mu_val / rho_val)
809 call field_add2s2(f_z, dz_p_adj, -mu_val / rho_val)
811 call neko_scratch_registry%relinquish_field(cc_indices)
814 call this%update_material_properties(time)
818 call profiler_start_region(
'Adjoint_velocity_residual')
820 call vel_res%compute(ax_vel, u, v, w, &
821 u_res, v_res, w_res, &
825 mu, rho, ext_bdf%diffusion_coeffs%x(1), &
828 call gs_xh%op(u_res, gs_op_add, event)
829 call device_event_sync(event)
830 call gs_xh%op(v_res, gs_op_add, event)
831 call device_event_sync(event)
832 call gs_xh%op(w_res, gs_op_add, event)
833 call device_event_sync(event)
836 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
838 call profiler_end_region(
'Adjoint_velocity_residual')
840 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
841 tstep, c_xh, n, dt_controller,
'Velocity')
843 call this%pc_vel%update()
845 call profiler_start_region(
"Adjoint_velocity_solve")
846 ksp_results(1:3) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
847 u_res%x, v_res%x, w_res%x, n, c_xh, &
848 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
849 this%ksp_vel%max_iter)
851 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
852 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
855 if (neko_bcknd_device .eq. 1)
then
856 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
857 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
859 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
861 call profiler_end_region(
"Adjoint_velocity_solve")
867 call field_copy(f_x, u)
868 call field_copy(f_y, v)
869 call field_copy(f_z, w)
871 call profiler_start_region(
'Adjoint_pressure_residual')
873 call prs_res%compute(p, p_res, &
877 this%bc_prs_surface, this%bc_sym_surface, &
878 ax_prs, ext_bdf%diffusion_coeffs%x(1), dt, &
882 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
883 call device_ortho(p_res%x_d, this%glb_n_points, n)
884 else if (.not. this%prs_dirichlet)
then
885 call ortho(p_res%x, this%glb_n_points, n)
888 call gs_xh%op(p_res, gs_op_add, event)
889 call device_event_sync(event)
892 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
895 call profiler_end_region(
'Adjoint_pressure_residual')
898 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
901 call this%pc_prs%update()
903 call profiler_start_region(
'Adjoint_pressure_solve')
907 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
908 this%bclst_dp, gs_xh)
911 call profiler_end_region(
'Adjoint_pressure_solve')
913 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
914 this%bclst_dp, gs_xh, n, tstep, dt_controller)
917 call field_add2(p, dp, n)
918 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
919 call device_ortho(p%x_d, this%glb_n_points, n)
920 else if (.not. this%prs_dirichlet)
then
921 call ortho(p%x, this%glb_n_points, n)
924 ksp_results(4)%name =
'Adjoint Pressure'
925 ksp_results(1)%name =
'Adjoint Velocity U'
926 ksp_results(2)%name =
'Adjoint Velocity V'
927 ksp_results(3)%name =
'Adjoint Velocity W'
929 if (this%forced_flow_rate)
then
930 call neko_error(
'Forced flow rate is not implemented for the adjoint')
936 call neko_scratch_registry%request_field(dx_p_adj, temp_indices(1), &
938 call neko_scratch_registry%request_field(dy_p_adj, temp_indices(2), &
940 call neko_scratch_registry%request_field(dz_p_adj, temp_indices(3), &
944 call opgrad(dx_p_adj%x, dy_p_adj%x, dz_p_adj%x, this%p_adj%x, c_xh)
947 call gs_xh%op(dx_p_adj, gs_op_add, event)
948 call device_event_sync(event)
949 call gs_xh%op(dy_p_adj, gs_op_add, event)
950 call device_event_sync(event)
951 call gs_xh%op(dz_p_adj, gs_op_add, event)
952 call device_event_sync(event)
955 if (neko_bcknd_device .eq. 1)
then
956 call device_col2(dx_p_adj%x_d, c_xh%Binv_d, dx_p_adj%size())
957 call device_col2(dy_p_adj%x_d, c_xh%Binv_d, dx_p_adj%size())
958 call device_col2(dz_p_adj%x_d, c_xh%Binv_d, dx_p_adj%size())
962 call col2(dx_p_adj%x, c_xh%Binv, dx_p_adj%size())
963 call col2(dy_p_adj%x, c_xh%Binv, dx_p_adj%size())
964 call col2(dz_p_adj%x, c_xh%Binv, dx_p_adj%size())
967 if (neko_bcknd_device .eq. 1)
then
968 call device_opadd2cm(u%x_d, v%x_d, w%x_d, dx_p_adj%x_d, &
969 dy_p_adj%x_d, dz_p_adj%x_d, -1.0_rp, n, msh%gdim)
971 call opadd2cm(u%x, v%x, w%x, dx_p_adj%x, dy_p_adj%x, dz_p_adj%x, &
972 -1.0_rp, n, msh%gdim)
975 call neko_scratch_registry%relinquish_field(temp_indices)
978 call fluid_step_info(time, ksp_results, &
979 this%full_stress_formulation, this%strict_convergence)
982 call profiler_end_region(
'Adjoint')
984 end subroutine adjoint_fluid_pnpn_step
990 subroutine adjoint_fluid_pnpn_setup_bcs(this, user, params)
991 use mpi_f08,
only: mpi_in_place
992 class(adjoint_fluid_pnpn_t),
intent(inout) :: this
993 type(user_t),
target,
intent(in) :: user
994 type(json_file),
intent(inout) :: params
995 integer :: i, n_bcs, j, zone_size, global_zone_size, ierr
996 class(bc_t),
pointer :: bc_i
997 type(json_core) :: core
998 type(json_value),
pointer :: bc_object
999 type(json_file) :: bc_subdict
1002 logical,
allocatable :: marked_zones(:)
1003 integer,
allocatable :: zone_indices(:)
1004 character(len=:),
allocatable :: json_key
1007 call this%bclst_vel_res%init()
1008 call this%bclst_du%init()
1009 call this%bclst_dv%init()
1010 call this%bclst_dw%init()
1011 call this%bclst_dp%init()
1013 call this%bc_vel_res%init_from_components(this%c_Xh)
1014 call this%bc_du%init_from_components(this%c_Xh)
1015 call this%bc_dv%init_from_components(this%c_Xh)
1016 call this%bc_dw%init_from_components(this%c_Xh)
1017 call this%bc_dp%init_from_components(this%c_Xh)
1020 call this%bc_prs_surface%init_from_components(this%c_Xh)
1021 call this%bc_sym_surface%init_from_components(this%c_Xh)
1022 call this%bc_curl_curl%init_from_components(this%c_Xh)
1024 json_key =
'case.adjoint_fluid.boundary_conditions'
1027 if (params%valid_path(json_key))
then
1028 call params%info(json_key, n_children = n_bcs)
1029 call params%get_core(core)
1030 call params%get(json_key, bc_object, found)
1035 call this%bcs_vel%init(n_bcs)
1037 allocate(marked_zones(
size(this%msh%labeled_zones)))
1038 marked_zones = .false.
1042 call json_extract_item(core, bc_object, i, bc_subdict)
1044 call json_get(bc_subdict,
"zone_indices", zone_indices)
1049 do j = 1,
size(zone_indices)
1050 zone_size = this%msh%labeled_zones(zone_indices(j))%size
1051 call mpi_allreduce(zone_size, global_zone_size, 1, &
1052 mpi_integer, mpi_max, neko_comm, ierr)
1054 if (global_zone_size .eq. 0)
then
1055 write(error_unit,
'(A, A, I0, A, A, I0, A)')
"*** ERROR ***: ",&
1056 "Zone index ", zone_indices(j), &
1057 " is invalid as this zone has 0 size, meaning it does ", &
1058 "not in the mesh. Check adjoint boundary condition ", &
1063 if (marked_zones(zone_indices(j)) .eqv. .true.)
then
1064 write(error_unit,
'(A, A, I0, A, A, A, A)')
"*** ERROR ***: ", &
1065 "Zone with index ", zone_indices(j), &
1066 " has already been assigned a boundary condition. ", &
1067 "Please check your boundary_conditions entry for the ", &
1068 "adjoint and make sure that each zone index appears ", &
1069 "only in a single boundary condition."
1072 marked_zones(zone_indices(j)) = .true.
1077 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
1081 if (
associated(bc_i))
then
1087 type is (symmetry_t)
1095 call this%bclst_vel_res%append(bc_i)
1096 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1097 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1098 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1100 call this%bcs_vel%append(bc_i)
1102 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
1103 type is (non_normal_t)
1106 call this%bclst_vel_res%append(bc_i)
1107 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1108 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1109 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1110 type is (shear_stress_t)
1112 call this%bclst_vel_res%append(bc_i%symmetry)
1113 call this%bclst_du%append(bc_i%symmetry%bc_x)
1114 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1115 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1117 call this%bcs_vel%append(bc_i)
1118 type is (wall_model_bc_t)
1120 call this%bclst_vel_res%append(bc_i%symmetry)
1121 call this%bclst_du%append(bc_i%symmetry%bc_x)
1122 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1123 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1125 call this%bcs_vel%append(bc_i)
1132 if (bc_i%strong .eqv. .true.)
then
1133 call this%bc_vel_res%mark_facets(bc_i%marked_facet)
1134 call this%bc_du%mark_facets(bc_i%marked_facet)
1135 call this%bc_dv%mark_facets(bc_i%marked_facet)
1136 call this%bc_dw%mark_facets(bc_i%marked_facet)
1138 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
1142 call this%bc_curl_curl%mark_facets(bc_i%marked_facet)
1144 call this%bcs_vel%append(bc_i)
1150 do i = 1,
size(this%msh%labeled_zones)
1151 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1152 (marked_zones(i) .eqv. .false.))
then
1153 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
1154 "No adjoint boundary condition assigned to zone ", i
1162 call this%bcs_prs%init(n_bcs)
1166 call json_extract_item(core, bc_object, i, bc_subdict)
1168 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
1172 if (
associated(bc_i))
then
1173 call this%bcs_prs%append(bc_i)
1176 if (bc_i%strong .eqv. .true.)
then
1177 call this%bc_dp%mark_facets(bc_i%marked_facet)
1185 do i = 1,
size(this%msh%labeled_zones)
1186 if (this%msh%labeled_zones(i)%size .gt. 0)
then
1187 call neko_error(
"No boundary_conditions entry in the case file!")
1193 call this%bc_prs_surface%finalize()
1194 call this%bc_sym_surface%finalize()
1195 call this%bc_curl_curl%finalize()
1197 call this%bc_vel_res%finalize()
1198 call this%bc_du%finalize()
1199 call this%bc_dv%finalize()
1200 call this%bc_dw%finalize()
1201 call this%bc_dp%finalize()
1203 call this%bclst_vel_res%append(this%bc_vel_res)
1204 call this%bclst_du%append(this%bc_du)
1205 call this%bclst_dv%append(this%bc_dv)
1206 call this%bclst_dw%append(this%bc_dw)
1207 call this%bclst_dp%append(this%bc_dp)
1210 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1211 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1212 mpi_logical, mpi_lor, neko_comm)
1214 end subroutine adjoint_fluid_pnpn_setup_bcs
1217 subroutine adjoint_fluid_pnpn_write_boundary_conditions(this)
1218 class(adjoint_fluid_pnpn_t),
target,
intent(inout) :: this
1219 type(dirichlet_t) :: bdry_mask
1220 type(field_t),
pointer :: bdry_field
1221 type(file_t) :: bdry_file
1222 integer :: temp_index, i
1223 class(bc_t),
pointer :: bci
1224 character(len=LOG_SIZE) :: log_buf
1226 call neko_log%section(
"Adjoint boundary conditions")
1227 write(log_buf,
'(A)') &
1228 'Marking using integer keys in boundary_adjoint0.f00000'
1229 call neko_log%message(log_buf)
1230 write(log_buf,
'(A)')
'Condition-value pairs: '
1231 call neko_log%message(log_buf)
1232 write(log_buf,
'(A)')
' no_slip = 1'
1233 call neko_log%message(log_buf)
1234 write(log_buf,
'(A)')
' velocity_value = 2'
1235 call neko_log%message(log_buf)
1236 write(log_buf,
'(A)')
' outflow, normal_outflow (+dong) = 3'
1237 call neko_log%message(log_buf)
1238 write(log_buf,
'(A)')
' symmetry = 4'
1239 call neko_log%message(log_buf)
1240 write(log_buf,
'(A)')
' user_velocity_pointwise = 5'
1241 call neko_log%message(log_buf)
1242 write(log_buf,
'(A)')
' periodic = 6'
1243 call neko_log%message(log_buf)
1244 write(log_buf,
'(A)')
' user_velocity = 7'
1245 call neko_log%message(log_buf)
1246 write(log_buf,
'(A)')
' user_pressure = 8'
1247 call neko_log%message(log_buf)
1248 write(log_buf,
'(A)')
' shear_stress = 9'
1249 call neko_log%message(log_buf)
1250 write(log_buf,
'(A)')
' wall_modelling = 10'
1251 call neko_log%message(log_buf)
1252 write(log_buf,
'(A)')
' blasius_profile = 11'
1253 call neko_log%message(log_buf)
1254 call neko_log%end_section()
1256 call neko_scratch_registry%request_field(bdry_field, temp_index, .true.)
1260 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1261 call bdry_mask%mark_zone(this%msh%periodic)
1262 call bdry_mask%finalize()
1263 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1264 call bdry_mask%free()
1266 do i = 1, this%bcs_prs%size()
1267 bci => this%bcs_prs%get(i)
1268 select type (bc => bci)
1269 type is (zero_dirichlet_t)
1270 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1271 call bdry_mask%mark_facets(bci%marked_facet)
1272 call bdry_mask%finalize()
1273 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1274 call bdry_mask%free()
1275 type is (dong_outflow_t)
1276 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1277 call bdry_mask%mark_facets(bci%marked_facet)
1278 call bdry_mask%finalize()
1279 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1280 call bdry_mask%free()
1281 type is (field_dirichlet_t)
1282 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1283 call bdry_mask%mark_facets(bci%marked_facet)
1284 call bdry_mask%finalize()
1285 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1286 call bdry_mask%free()
1290 do i = 1, this%bcs_vel%size()
1291 bci => this%bcs_vel%get(i)
1292 select type (bc => bci)
1293 type is (zero_dirichlet_t)
1294 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1295 call bdry_mask%mark_facets(bci%marked_facet)
1296 call bdry_mask%finalize()
1297 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1298 call bdry_mask%free()
1300 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1301 call bdry_mask%mark_facets(bci%marked_facet)
1302 call bdry_mask%finalize()
1303 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1304 call bdry_mask%free()
1305 type is (symmetry_t)
1306 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1307 call bdry_mask%mark_facets(bci%marked_facet)
1308 call bdry_mask%finalize()
1309 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1310 call bdry_mask%free()
1311 type is (field_dirichlet_vector_t)
1312 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1313 call bdry_mask%mark_facets(bci%marked_facet)
1314 call bdry_mask%finalize()
1315 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1316 call bdry_mask%free()
1317 type is (shear_stress_t)
1318 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1319 call bdry_mask%mark_facets(bci%marked_facet)
1320 call bdry_mask%finalize()
1321 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1322 call bdry_mask%free()
1323 type is (wall_model_bc_t)
1324 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1325 call bdry_mask%mark_facets(bci%marked_facet)
1326 call bdry_mask%finalize()
1327 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1328 call bdry_mask%free()
1330 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1331 call bdry_mask%mark_facets(bci%marked_facet)
1332 call bdry_mask%finalize()
1333 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1334 call bdry_mask%free()
1339 call bdry_file%init(
'boundary_adjoint.fld')
1340 call bdry_file%write(bdry_field)
1342 call neko_scratch_registry%relinquish_field(temp_index)
1343 end subroutine adjoint_fluid_pnpn_write_boundary_conditions
1348 subroutine rescale_fluid(fluid_data, scale)
1350 class(adjoint_fluid_pnpn_t),
intent(inout) :: fluid_data
1352 real(kind=rp),
intent(in) :: scale
1358 if (neko_bcknd_device .eq. 1)
then
1359 call device_cmult(fluid_data%u_adj%x_d, scale, fluid_data%u_adj%size())
1360 call device_cmult(fluid_data%v_adj%x_d, scale, fluid_data%v_adj%size())
1361 call device_cmult(fluid_data%w_adj%x_d, scale, fluid_data%w_adj%size())
1363 call cmult(fluid_data%u_adj%x, scale, fluid_data%u_adj%size())
1364 call cmult(fluid_data%v_adj%x, scale, fluid_data%v_adj%size())
1365 call cmult(fluid_data%w_adj%x, scale, fluid_data%w_adj%size())
1369 if (neko_bcknd_device .eq. 1)
then
1370 call device_cmult(fluid_data%f_adj_x%x_d, scale, &
1371 fluid_data%f_adj_x%size())
1372 call device_cmult(fluid_data%f_adj_y%x_d, scale, &
1373 fluid_data%f_adj_y%size())
1374 call device_cmult(fluid_data%f_adj_z%x_d, scale, &
1375 fluid_data%f_adj_z%size())
1378 call device_cmult(fluid_data%abx1%x_d, scale, fluid_data%abx1%size())
1379 call device_cmult(fluid_data%aby1%x_d, scale, fluid_data%aby1%size())
1380 call device_cmult(fluid_data%abz1%x_d, scale, fluid_data%abz1%size())
1381 call device_cmult(fluid_data%abx2%x_d, scale, fluid_data%abx2%size())
1382 call device_cmult(fluid_data%aby2%x_d, scale, fluid_data%aby2%size())
1383 call device_cmult(fluid_data%abz2%x_d, scale, fluid_data%abz2%size())
1386 call cmult(fluid_data%f_adj_x%x, scale, fluid_data%f_adj_x%size())
1387 call cmult(fluid_data%f_adj_y%x, scale, fluid_data%f_adj_y%size())
1388 call cmult(fluid_data%f_adj_z%x, scale, fluid_data%f_adj_z%size())
1390 call cmult(fluid_data%abx1%x, scale, fluid_data%abx1%size())
1391 call cmult(fluid_data%aby1%x, scale, fluid_data%aby1%size())
1392 call cmult(fluid_data%abz1%x, scale, fluid_data%abz1%size())
1394 call cmult(fluid_data%abx2%x, scale, fluid_data%abx2%size())
1395 call cmult(fluid_data%aby2%x, scale, fluid_data%aby2%size())
1396 call cmult(fluid_data%abz2%x, scale, fluid_data%abz2%size())
1400 if (neko_bcknd_device .eq. 1)
then
1401 do i = 1, fluid_data%ulag%size()
1402 call device_cmult(fluid_data%ulag%lf(i)%x_d, &
1403 scale, fluid_data%ulag%lf(i)%size())
1406 do i = 1, fluid_data%vlag%size()
1407 call device_cmult(fluid_data%vlag%lf(i)%x_d, &
1408 scale, fluid_data%vlag%lf(i)%size())
1411 do i = 1, fluid_data%wlag%size()
1412 call device_cmult(fluid_data%wlag%lf(i)%x_d, &
1413 scale, fluid_data%wlag%lf(i)%size())
1416 do i = 1, fluid_data%ulag%size()
1417 call cmult(fluid_data%ulag%lf(i)%x, &
1418 scale, fluid_data%ulag%lf(i)%size())
1421 do i = 1, fluid_data%vlag%size()
1422 call cmult(fluid_data%vlag%lf(i)%x, &
1423 scale, fluid_data%vlag%lf(i)%size())
1426 do i = 1, fluid_data%wlag%size()
1427 call cmult(fluid_data%wlag%lf(i)%x, &
1428 scale, fluid_data%wlag%lf(i)%size())
1432 end subroutine rescale_fluid
1434 function norm(x, y, z, B, volume, n)
1435 use mpi_f08,
only: mpi_in_place
1437 integer,
intent(in) :: n
1438 real(kind=rp),
dimension(n),
intent(in) :: x, y, z
1439 real(kind=rp),
dimension(n),
intent(in) :: b
1440 real(kind=rp),
intent(in) :: volume
1442 real(kind=rp) :: norm
1444 norm = vlsc3(x, x, b, n) + vlsc3(y, y, b, n) + vlsc3(z, z, b, n)
1446 call mpi_allreduce(mpi_in_place, norm, 1, &
1447 mpi_real_precision, mpi_sum, neko_comm)
1449 norm = sqrt(norm / volume)
1452 function device_norm(x_d, y_d, z_d, B_d, volume, n)
1453 use mpi_f08,
only: mpi_in_place
1455 type(c_ptr),
intent(in) :: x_d, y_d, z_d
1456 type(c_ptr),
intent(in) :: B_d
1457 real(kind=rp),
intent(in) :: volume
1458 integer,
intent(in) :: n
1460 real(kind=rp) :: device_norm
1462 device_norm = device_vlsc3(x_d, x_d, b_d, n) + &
1463 device_vlsc3(y_d, y_d, b_d, n) + &
1464 device_vlsc3(z_d, z_d, b_d, n)
1466 call mpi_allreduce(mpi_in_place, device_norm, 1, &
1467 mpi_real_precision, mpi_sum, neko_comm)
1469 device_norm = sqrt(device_norm / volume)
1471 end function device_norm
1477 subroutine power_iterations_compute(this, t, tstep)
1478 class(adjoint_fluid_pnpn_t),
target,
intent(inout) :: this
1479 real(kind=rp),
intent(in) :: t
1480 integer,
intent(in) :: tstep
1483 real(kind=rp) :: scaling_factor
1484 real(kind=rp) :: norm_l2, norm_l2_base
1485 character(len=256) :: log_message
1486 type(vector_t) :: data_line
1489 n = this%c_Xh%dof%size()
1490 if (tstep .eq. 1)
then
1491 if (neko_bcknd_device .eq. 1)
then
1492 norm_l2_base = device_norm(this%u_adj%x_d, this%v_adj%x_d, &
1494 this%c_Xh%B_d, this%c_Xh%volume, n)
1496 norm_l2_base = this%norm_scaling * norm(this%u_adj%x, this%v_adj%x, &
1498 this%c_Xh%B, this%c_Xh%volume, n)
1500 if (this%norm_target .lt. 0.0_rp)
then
1501 this%norm_target = norm_l2_base
1504 this%norm_l2_upper = this%norm_tolerance * this%norm_target
1505 this%norm_l2_lower = this%norm_target / this%norm_tolerance
1510 if (neko_bcknd_device .eq. 1)
then
1511 norm_l2 = device_norm(this%u_adj%x_d, this%v_adj%x_d, this%w_adj%x_d, &
1512 this%c_Xh%B_d, this%c_Xh%volume, n)
1514 norm_l2 = norm(this%u_adj%x, this%v_adj%x, this%w_adj%x, &
1515 this%c_Xh%B, this%c_Xh%volume, n)
1517 norm_l2 = sqrt(this%norm_scaling) * norm_l2
1518 scaling_factor = 1.0_rp
1521 if (norm_l2 .gt. this%norm_l2_upper &
1522 .or. norm_l2 .lt. this%norm_l2_lower)
then
1523 scaling_factor = this%norm_target / norm_l2
1524 call rescale_fluid(this, scaling_factor)
1525 norm_l2 = this%norm_target
1527 if (tstep .eq. 1)
then
1528 scaling_factor = 1.0_rp
1534 call neko_log%section(
'Power Iterations')
1536 write (log_message,
'(A7,E20.14)')
'Norm: ', norm_l2
1537 call neko_log%message(log_message, lvl = neko_log_debug)
1538 write (log_message,
'(A7,E20.14)')
'Scaling: ', scaling_factor
1539 call neko_log%message(log_message, lvl = neko_log_debug)
1542 call data_line%init(2)
1543 data_line%x = [norm_l2, scaling_factor]
1544 call this%file_output%write(data_line, t)
1547 call neko_log%end_section(
'Power Iterations')
1548 end subroutine power_iterations_compute
Boundary condition factory for pressure.
Adjoint Pn/Pn formulation.
Subroutines to add advection terms to the RHS of a transport equation.
Base type of all fluid formulations.
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.