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