250 subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
251 class(fluid_pnpn_t),
target,
intent(inout) :: this
252 type(
mesh_t),
target,
intent(inout) :: msh
253 integer,
intent(in) :: lx
254 type(json_file),
target,
intent(inout) :: params
255 type(
user_t),
target,
intent(in) :: user
256 type(
chkp_t),
target,
intent(inout) :: chkp
257 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
259 class(
bc_t),
pointer :: bc_i, vel_bc
260 real(kind=
rp) :: abs_tol
261 character(len=LOG_SIZE) :: log_buf
262 integer :: ierr, integer_val, solver_maxiter
263 character(len=:),
allocatable :: solver_type, precon_type
264 logical :: monitor, found
266 type(json_file) :: numerics_params, precon_params
271 call this%init_base(msh, lx, params, scheme,
user, .true.)
283 allocate(this%ext_bdf)
284 call this%ext_bdf%init(integer_val)
287 this%full_stress_formulation, .false.)
291 call this%c_Xh%generate_cyclic_bc()
293 if (this%full_stress_formulation)
then
295 call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
298 call pnpn_prs_res_stress_factory(this%prs_res)
301 call pnpn_vel_res_stress_factory(this%vel_res)
304 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
307 call pnpn_prs_res_factory(this%prs_res)
310 call pnpn_vel_res_factory(this%vel_res)
315 if (params%valid_path(
'case.fluid.nut_field'))
then
316 if (.not. this%full_stress_formulation)
then
317 call neko_error(
"You need to set full_stress_formulation to " // &
318 "true for the fluid to have a spatially varying " // &
321 call json_get(params,
'case.fluid.nut_field', this%nut_field_name)
323 this%nut_field_name =
""
327 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
331 call rhs_maker_sumab_fctry(this%sumab)
334 call rhs_maker_ext_fctry(this%makeabf)
337 call rhs_maker_bdf_fctry(this%makebdf)
340 call rhs_maker_oifs_fctry(this%makeoifs)
343 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
344 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
346 call this%p_res%init(dm_xh,
"p_res")
347 call this%u_res%init(dm_xh,
"u_res")
348 call this%v_res%init(dm_xh,
"v_res")
349 call this%w_res%init(dm_xh,
"w_res")
350 call this%abx1%init(dm_xh,
"abx1")
351 call this%aby1%init(dm_xh,
"aby1")
352 call this%abz1%init(dm_xh,
"abz1")
353 call this%abx2%init(dm_xh,
"abx2")
354 call this%aby2%init(dm_xh,
"aby2")
355 call this%abz2%init(dm_xh,
"abz2")
356 call this%advx%init(dm_xh,
"advx")
357 call this%advy%init(dm_xh,
"advy")
358 call this%advz%init(dm_xh,
"advz")
361 call this%du%init(this%dm_Xh,
'du')
362 call this%dv%init(this%dm_Xh,
'dv')
363 call this%dw%init(this%dm_Xh,
'dw')
364 call this%dp%init(this%dm_Xh,
'dp')
366 call this%ale%init(this%c_Xh, params,
user)
368 call neko_log%section(
"Fluid boundary conditions")
370 call this%setup_bcs(
user, params)
374 if (found)
call this%write_boundary_conditions()
377 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
378 this%pr_projection_activ_step, &
379 this%pr_projection_reorthogonalize_basis)
381 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
382 this%vel_projection_activ_step)
388 if (params%valid_path(
'case.fluid.flow_rate_force'))
then
389 call this%vol_flow%init(this%dm_Xh, params)
393 call neko_log%section(
"Pressure solver")
396 'case.fluid.pressure_solver.max_iterations', &
398 call json_get(params,
'case.fluid.pressure_solver.type', solver_type)
399 call json_get(params,
'case.fluid.pressure_solver.preconditioner.type', &
402 'case.fluid.pressure_solver.preconditioner', precon_params)
404 'case.fluid.pressure_solver.absolute_tolerance', &
408 call neko_log%message(
'Type : ('// trim(solver_type) // &
409 ', ' // trim(precon_type) //
')')
410 write(log_buf,
'(A,ES13.6)')
'Abs tol :', abs_tol
413 call this%solver_factory(this%ksp_prs, this%dm_Xh%size(), &
414 solver_type, solver_maxiter, abs_tol, monitor)
415 call this%precon_factory_(this%pc_prs, this%ksp_prs, &
416 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, &
417 precon_type, precon_params)
422 call json_get(params,
'case.numerics', numerics_params)
423 call advection_factory(this%adv, numerics_params, this%c_Xh, &
424 this%ulag, this%vlag, this%wlag, &
425 chkp%dtlag, chkp%tlag, this%ext_bdf, &
430 call this%chkp%add_fluid(this%u, this%v, this%w, this%p)
432 this%chkp%abx1 => this%abx1
433 this%chkp%abx2 => this%abx2
434 this%chkp%aby1 => this%aby1
435 this%chkp%aby2 => this%aby2
436 this%chkp%abz1 => this%abz1
437 this%chkp%abz2 => this%abz2
438 call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
441 if (this%ale%active)
then
442 call this%chkp%add_ale(this%c_Xh%dof%x, this%c_Xh%dof%y, &
444 this%c_Xh%Blag, this%c_Xh%Blaglag, &
445 this%ale%wm_x, this%ale%wm_y, this%ale%wm_z, &
446 this%ale%wm_x_lag, this%ale%wm_y_lag, &
448 this%ale%global_pivot_pos, &
449 this%ale%global_pivot_vel_lag, &
450 this%ale%global_basis_pos, &
451 this%ale%global_basis_vel_lag)
458 subroutine fluid_pnpn_restart(this, chkp)
459 class(fluid_pnpn_t),
target,
intent(inout) :: this
460 type(chkp_t),
intent(inout) :: chkp
461 real(kind=rp) :: dtlag(10), tlag(10)
462 type(field_t) :: u_temp, v_temp, w_temp
468 n = this%u%dof%size()
469 if (
allocated(chkp%previous_mesh%elements) .or. &
470 chkp%previous_Xh%lx .ne. this%Xh%lx)
then
471 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
472 c_xh => this%c_Xh, ulag => this%ulag, vlag => this%vlag, &
474 do concurrent(j = 1:n)
475 u%x(j,1,1,1) = u%x(j,1,1,1) * c_xh%mult(j,1,1,1)
476 v%x(j,1,1,1) = v%x(j,1,1,1) * c_xh%mult(j,1,1,1)
477 w%x(j,1,1,1) = w%x(j,1,1,1) * c_xh%mult(j,1,1,1)
478 p%x(j,1,1,1) = p%x(j,1,1,1) * c_xh%mult(j,1,1,1)
480 do i = 1, this%ulag%size()
481 do concurrent(j = 1:n)
482 ulag%lf(i)%x(j,1,1,1) = ulag%lf(i)%x(j,1,1,1) &
484 vlag%lf(i)%x(j,1,1,1) = vlag%lf(i)%x(j,1,1,1) &
486 wlag%lf(i)%x(j,1,1,1) = wlag%lf(i)%x(j,1,1,1) &
493 if (neko_bcknd_device .eq. 1)
then
494 associate(u => this%u, v => this%v, w => this%w, &
495 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
497 call device_memcpy(u%x, u%x_d, u%dof%size(), &
498 host_to_device, sync = .false.)
499 call device_memcpy(v%x, v%x_d, v%dof%size(), &
500 host_to_device, sync = .false.)
501 call device_memcpy(w%x, w%x_d, w%dof%size(), &
502 host_to_device, sync = .false.)
503 call device_memcpy(p%x, p%x_d, p%dof%size(), &
504 host_to_device, sync = .false.)
505 call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
506 u%dof%size(), host_to_device, sync = .false.)
507 call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
508 u%dof%size(), host_to_device, sync = .false.)
510 call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
511 v%dof%size(), host_to_device, sync = .false.)
512 call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
513 v%dof%size(), host_to_device, sync = .false.)
515 call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
516 w%dof%size(), host_to_device, sync = .false.)
517 call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
518 w%dof%size(), host_to_device, sync = .false.)
519 call device_memcpy(this%abx1%x, this%abx1%x_d, &
520 w%dof%size(), host_to_device, sync = .false.)
521 call device_memcpy(this%abx2%x, this%abx2%x_d, &
522 w%dof%size(), host_to_device, sync = .false.)
523 call device_memcpy(this%aby1%x, this%aby1%x_d, &
524 w%dof%size(), host_to_device, sync = .false.)
525 call device_memcpy(this%aby2%x, this%aby2%x_d, &
526 w%dof%size(), host_to_device, sync = .false.)
527 call device_memcpy(this%abz1%x, this%abz1%x_d, &
528 w%dof%size(), host_to_device, sync = .false.)
529 call device_memcpy(this%abz2%x, this%abz2%x_d, &
530 w%dof%size(), host_to_device, sync = .false.)
531 call device_memcpy(this%advx%x, this%advx%x_d, &
532 w%dof%size(), host_to_device, sync = .false.)
533 call device_memcpy(this%advy%x, this%advy%x_d, &
534 w%dof%size(), host_to_device, sync = .false.)
535 call device_memcpy(this%advz%x, this%advz%x_d, &
536 w%dof%size(), host_to_device, sync = .false.)
543 if (
allocated(chkp%previous_mesh%elements) &
544 .or. chkp%previous_Xh%lx .ne. this%Xh%lx)
then
546 call rotate_cyc(this%u%x, this%v%x, this%w%x, 1, this%c_Xh)
547 call this%gs_Xh%op(this%u, gs_op_add)
548 call this%gs_Xh%op(this%v, gs_op_add)
549 call this%gs_Xh%op(this%w, gs_op_add)
550 call this%gs_Xh%op(this%p, gs_op_add)
551 call rotate_cyc(this%u%x, this%v%x, this%w%x, 0, this%c_Xh)
553 do i = 1, this%ulag%size()
554 call rotate_cyc(this%ulag%lf(i)%x, this%vlag%lf(i)%x, &
555 this%wlag%lf(i)%x, 1, this%c_Xh)
556 call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
557 call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
558 call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
559 call rotate_cyc(this%ulag%lf(i)%x, this%vlag%lf(i)%x, &
560 this%wlag%lf(i)%x, 0, this%c_Xh)
564 call this%ale%set_coef_restart(this%c_Xh, this%adv, chkp%t)
671 subroutine fluid_pnpn_step(this, time, dt_controller)
672 class(fluid_pnpn_t),
target,
intent(inout) :: this
673 type(time_state_t),
intent(in) :: time
674 type(time_step_controller_t),
intent(in) :: dt_controller
678 type(ksp_monitor_t) :: ksp_results(4)
680 type(file_t) :: dump_file
681 class(bc_t),
pointer :: bc_i
682 type(non_normal_t),
pointer :: bc_j
684 if (this%freeze)
return
686 n = this%dm_Xh%size()
688 call profiler_start_region(
'Fluid', 1)
689 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
690 u_e => this%u_e, v_e => this%v_e, w_e => this%w_e, &
691 du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
692 u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
693 p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
695 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
696 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
697 msh => this%msh, prs_res => this%prs_res, &
698 source_term => this%source_term, vel_res => this%vel_res, &
699 sumab => this%sumab, makeoifs => this%makeoifs, &
700 makeabf => this%makeabf, makebdf => this%makebdf, &
701 vel_projection_dim => this%vel_projection_dim, &
702 pr_projection_dim => this%pr_projection_dim, &
704 rho => this%rho, mu_tot => this%mu_tot, &
705 f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
706 t => time%t, tstep => time%tstep, dt => time%dt, &
707 ext_bdf => this%ext_bdf, event => glb_cmd_event, &
711 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
712 ulag, vlag, wlag, ext_bdf%advection_coeffs%x, ext_bdf%nadv)
715 call this%source_term%compute(time)
718 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
719 this%dm_Xh%size(), time, strong = .false.)
721 if (this%ale%active)
then
723 call neko_error(
"ALE is not yet supported " // &
724 "with OIFS time integration.")
727 call this%adv%compute_ale(u, v, w, &
728 ale%wm_x, ale%wm_y, ale%wm_z, &
730 xh, c_xh, dm_xh%size())
736 call this%adv%compute(u, v, w, &
737 this%advx, this%advy, this%advz, &
738 xh, this%c_Xh, dm_xh%size(), dt)
745 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
746 this%abx2, this%aby2, this%abz2, &
747 f_x%x, f_y%x, f_z%x, &
748 rho%x(1,1,1,1), ext_bdf%advection_coeffs%x, n)
751 call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
752 f_x%x, f_y%x, f_z%x, &
753 rho%x(1,1,1,1), dt, n)
756 call this%adv%compute(u, v, w, &
758 xh, this%c_Xh, dm_xh%size())
765 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
766 this%abx2, this%aby2, this%abz2, &
767 f_x%x, f_y%x, f_z%x, &
768 rho%x(1,1,1,1), ext_bdf%advection_coeffs%x, n)
774 call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
775 u, v, w, c_xh%B, rho%x(1,1,1,1), dt, &
776 ext_bdf%diffusion_coeffs%x, ext_bdf%ndiff, n, &
777 c_xh%Blag, c_xh%Blaglag)
781 if (this%ale%active)
then
783 call this%ale%advance_mesh(c_xh, time, ext_bdf%nadv)
785 call profiler_start_region(
'ALE recompute metrics')
787 call c_xh%recompute_metrics()
790 call this%adv%recompute_metrics(c_xh, .true.)
791 call profiler_end_region(
'ALE recompute metrics')
798 call this%bc_apply_vel(time, strong = .true.)
799 call this%bc_apply_prs(time)
802 call this%update_material_properties(time)
805 call profiler_start_region(
'Pressure_residual', 18)
806 call prs_res%compute(p, p_res,&
811 this%bc_prs_surface, this%bc_sym_surface,&
812 ax_prs, ext_bdf%diffusion_coeffs%x(1), dt, &
817 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
818 call device_ortho(p_res%x_d, this%glb_n_points, n)
819 else if (.not. this%prs_dirichlet)
then
820 call ortho(p_res%x, this%glb_n_points, n)
823 call gs_xh%op(p_res, gs_op_add, event)
824 call device_event_sync(event)
827 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
830 call profiler_end_region(
'Pressure_residual', 18)
832 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
833 ax = ax_prs, gs_h = gs_xh, bclst = this%bclst_dp, &
836 call this%pc_prs%update()
838 call profiler_start_region(
'Pressure_solve', 3)
842 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
843 this%bclst_dp, gs_xh)
844 ksp_results(1)%name =
'Pressure'
847 call profiler_end_region(
'Pressure_solve', 3)
849 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
850 this%bclst_dp, gs_xh, n, tstep, dt_controller)
853 call field_add2(p, dp, n)
854 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
855 call device_ortho(p%x_d, this%glb_n_points, n)
856 else if (.not. this%prs_dirichlet)
then
857 call ortho(p%x, this%glb_n_points, n)
861 call profiler_start_region(
'Velocity_residual', 19)
862 call vel_res%compute(ax_vel, u, v, w, &
863 u_res, v_res, w_res, &
867 mu_tot, rho, ext_bdf%diffusion_coeffs%x(1), &
870 call rotate_cyc(u_res%x, v_res%x, w_res%x, 1, c_xh)
871 call gs_xh%op(u_res, gs_op_add, event)
872 call device_event_sync(event)
873 call gs_xh%op(v_res, gs_op_add, event)
874 call device_event_sync(event)
875 call gs_xh%op(w_res, gs_op_add, event)
876 call device_event_sync(event)
877 call rotate_cyc(u_res%x, v_res%x, w_res%x, 0, c_xh)
880 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
883 call profiler_end_region(
'Velocity_residual', 19)
885 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
886 tstep, c_xh, n, dt_controller,
'Velocity')
888 call this%pc_vel%update()
890 call profiler_start_region(
"Velocity_solve", 4)
891 ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
892 u_res%x, v_res%x, w_res%x, n, c_xh, &
893 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
894 this%ksp_vel%max_iter)
895 call profiler_end_region(
"Velocity_solve", 4)
896 if (this%full_stress_formulation)
then
897 ksp_results(2)%name =
'Momentum'
899 ksp_results(2)%name =
'X-Velocity'
900 ksp_results(3)%name =
'Y-Velocity'
901 ksp_results(4)%name =
'Z-Velocity'
904 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
905 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
908 if (neko_bcknd_device .eq. 1)
then
909 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
910 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
912 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
915 if (this%forced_flow_rate)
then
917 call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
918 c_xh, gs_xh, ext_bdf, rho%x(1,1,1,1), mu_tot, &
919 dt, time, this%bclst_dp, this%bclst_du, this%bclst_dv, &
920 this%bclst_dw, this%bclst_vel_res, ax_vel, ax_prs, this%ksp_prs, &
921 this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
922 this%ksp_vel%max_iter)
928 call this%ale%update_mesh_velocity(c_xh, time)
930 call fluid_step_info(time, ksp_results, &
931 this%full_stress_formulation, this%strict_convergence, &
932 this%allow_stabilization)
935 call profiler_end_region(
'Fluid', 1)
940 subroutine fluid_pnpn_setup_bcs(this, user, params)
941 class(fluid_pnpn_t),
target,
intent(inout) :: this
942 type(user_t),
target,
intent(in) :: user
943 type(json_file),
intent(inout) :: params
944 integer :: i, n_bcs, zone_index, j, zone_size, global_zone_size, ierr
945 class(bc_t),
pointer :: bc_i
946 type(json_core) :: core
947 type(json_value),
pointer :: bc_object
948 type(json_file) :: bc_subdict
949 logical :: ale_active_local, any_moving_wall, moving_
952 logical,
allocatable :: marked_zones(:)
953 integer,
allocatable :: zone_indices(:)
956 character(len=:),
allocatable :: bc_type_str
957 this%ale%has_moving_boundary = .false.
958 any_moving_wall = .false.
959 ale_active_local = .false.
960 call json_get_or_default(params, &
961 'case.fluid.ale.enabled', &
962 ale_active_local, .false.)
965 call this%bclst_vel_res%init()
966 call this%bclst_du%init()
967 call this%bclst_dv%init()
968 call this%bclst_dw%init()
969 call this%bclst_dp%init()
971 call this%bc_vel_res%init_from_components(this%c_Xh)
972 call this%bc_du%init_from_components(this%c_Xh)
973 call this%bc_dv%init_from_components(this%c_Xh)
974 call this%bc_dw%init_from_components(this%c_Xh)
975 call this%bc_dp%init_from_components(this%c_Xh)
978 call this%bc_prs_surface%init_from_components(this%c_Xh)
979 call this%bc_sym_surface%init_from_components(this%c_Xh)
982 if (params%valid_path(
'case.fluid.boundary_conditions'))
then
983 call params%info(
'case.fluid.boundary_conditions', n_children = n_bcs)
984 call params%get_core(core)
985 call params%get(
'case.fluid.boundary_conditions', bc_object, found)
990 call this%bcs_vel%init(n_bcs)
992 allocate(marked_zones(
size(this%msh%labeled_zones)))
993 marked_zones = .false.
997 call json_extract_item(core, bc_object, i, bc_subdict)
999 call json_get_or_lookup(bc_subdict,
"zone_indices", zone_indices)
1002 call json_get(bc_subdict,
"type", bc_type_str)
1004 if (trim(bc_type_str) .eq.
"no_slip")
then
1005 call json_get_or_default(bc_subdict,
"moving", moving_, .false.)
1008 this%ale%has_moving_boundary = .true.
1014 do j = 1,
size(zone_indices)
1015 zone_size = this%msh%labeled_zones(zone_indices(j))%size
1016 call mpi_allreduce(zone_size, global_zone_size, 1, &
1017 mpi_integer, mpi_max, neko_comm, ierr)
1019 if (global_zone_size .eq. 0)
then
1020 write(error_unit,
'(A, A, I0, A, A, I0, A)')
"*** ERROR ***: ",&
1021 "Zone index ", zone_indices(j), &
1022 " is invalid as this zone has 0 size, meaning it ", &
1023 "is not in the mesh. Check fluid boundary condition ", &
1028 if (marked_zones(zone_indices(j)))
then
1029 write(error_unit,
'(A, A, I0, A, A, A, A)')
"*** ERROR ***: ", &
1030 "Zone with index ", zone_indices(j), &
1031 " has already been assigned a boundary condition. ", &
1032 "Please check your boundary_conditions entry for the ", &
1033 "fluid and make sure that each zone index appears only ", &
1034 "in a single boundary condition."
1037 marked_zones(zone_indices(j)) = .true.
1042 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
1046 if (
associated(bc_i))
then
1052 type is (symmetry_t)
1060 call this%bclst_vel_res%append(bc_i)
1061 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1062 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1063 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1065 call this%bcs_vel%append(bc_i)
1067 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
1068 type is (non_normal_t)
1071 call this%bclst_vel_res%append(bc_i)
1072 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1073 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1074 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1075 type is (shear_stress_t)
1077 call this%bclst_vel_res%append(bc_i%symmetry)
1078 call this%bclst_du%append(bc_i%symmetry%bc_x)
1079 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1080 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1082 call this%bcs_vel%append(bc_i)
1083 type is (wall_model_bc_t)
1085 call this%bclst_vel_res%append(bc_i%symmetry)
1086 call this%bclst_du%append(bc_i%symmetry%bc_x)
1087 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1088 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1090 call this%bcs_vel%append(bc_i)
1097 if (bc_i%strong)
then
1098 call this%bc_vel_res%mark_facets(bc_i%marked_facet)
1099 call this%bc_du%mark_facets(bc_i%marked_facet)
1100 call this%bc_dv%mark_facets(bc_i%marked_facet)
1101 call this%bc_dw%mark_facets(bc_i%marked_facet)
1103 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
1106 call this%bcs_vel%append(bc_i)
1111 if (this%ale%active .and. (.not. this%ale%has_moving_boundary))
then
1112 call neko_error(
"Case file error: ALE is active, " // &
1113 "but no moving wall was found. " // &
1114 "Use type = 'no_slip' with 'moving': true in case file.")
1118 do i = 1,
size(this%msh%labeled_zones)
1119 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1120 (.not. marked_zones(i)))
then
1121 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
1122 "No fluid boundary condition assigned to zone ", i
1130 call this%bcs_prs%init(n_bcs)
1134 call json_extract_item(core, bc_object, i, bc_subdict)
1136 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
1140 if (
associated(bc_i))
then
1141 call this%bcs_prs%append(bc_i)
1144 if (bc_i%strong)
then
1145 call this%bc_dp%mark_facets(bc_i%marked_facet)
1153 do i = 1,
size(this%msh%labeled_zones)
1154 if (this%msh%labeled_zones(i)%size .gt. 0)
then
1155 call neko_error(
"No boundary_conditions entry in the case file!")
1161 call this%bcs_vel%init()
1162 call this%bcs_prs%init()
1166 call this%bc_prs_surface%finalize()
1167 call this%bc_sym_surface%finalize()
1169 call this%bc_vel_res%finalize()
1170 call this%bc_du%finalize()
1171 call this%bc_dv%finalize()
1172 call this%bc_dw%finalize()
1173 call this%bc_dp%finalize()
1175 call this%bclst_vel_res%append(this%bc_vel_res)
1176 call this%bclst_du%append(this%bc_du)
1177 call this%bclst_dv%append(this%bc_dv)
1178 call this%bclst_dw%append(this%bc_dw)
1179 call this%bclst_dp%append(this%bc_dp)
1182 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1183 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1184 mpi_logical, mpi_lor, neko_comm)
1187 if (
allocated(marked_zones))
then
1188 deallocate(marked_zones)
1191 if (
allocated(zone_indices))
then
1192 deallocate(zone_indices)
1198 subroutine fluid_pnpn_write_boundary_conditions(this)
1205 class(fluid_pnpn_t),
target,
intent(inout) :: this
1206 type(dirichlet_t) :: bdry_mask
1207 type(field_t),
pointer :: bdry_field
1208 type(file_t) :: bdry_file
1209 integer :: temp_index, i
1210 class(bc_t),
pointer :: bci
1211 character(len=LOG_SIZE) :: log_buf
1213 write(log_buf,
'(A)')
'Marking using integer keys in bdry0.f00000'
1214 call neko_log%message(log_buf)
1215 write(log_buf,
'(A)')
'Condition-value pairs: '
1216 call neko_log%message(log_buf)
1217 write(log_buf,
'(A)')
' no_slip (stationary wall) = 1'
1218 call neko_log%message(log_buf)
1219 write(log_buf,
'(A)')
' velocity_value = 2'
1220 call neko_log%message(log_buf)
1221 write(log_buf,
'(A)')
' outflow, normal_outflow (+dong) = 3'
1222 call neko_log%message(log_buf)
1223 write(log_buf,
'(A)')
' symmetry = 4'
1224 call neko_log%message(log_buf)
1225 write(log_buf,
'(A)')
' periodic = 6'
1226 call neko_log%message(log_buf)
1227 write(log_buf,
'(A)')
' user_velocity = 7'
1228 call neko_log%message(log_buf)
1229 write(log_buf,
'(A)')
' user_pressure = 8'
1230 call neko_log%message(log_buf)
1231 write(log_buf,
'(A)')
' shear_stress = 9'
1232 call neko_log%message(log_buf)
1233 write(log_buf,
'(A)')
' wall_modelling = 10'
1234 call neko_log%message(log_buf)
1235 write(log_buf,
'(A)')
' blasius_profile = 11'
1236 call neko_log%message(log_buf)
1237 write(log_buf,
'(A)')
' no_slip (moving wall) = 12'
1238 call neko_log%message(log_buf)
1240 call neko_scratch_registry%request_field(bdry_field, temp_index, .true.)
1244 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1245 call bdry_mask%mark_zone(this%msh%periodic)
1246 call bdry_mask%finalize()
1247 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1248 call bdry_mask%free()
1250 do i = 1, this%bcs_prs%size()
1251 bci => this%bcs_prs%get(i)
1252 select type (
bc => bci)
1253 type is (zero_dirichlet_t)
1254 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1255 call bdry_mask%mark_facets(bci%marked_facet)
1256 call bdry_mask%finalize()
1257 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1258 call bdry_mask%free()
1260 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1261 call bdry_mask%mark_facets(bci%marked_facet)
1262 call bdry_mask%finalize()
1263 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1264 call bdry_mask%free()
1266 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1267 call bdry_mask%mark_facets(bci%marked_facet)
1268 call bdry_mask%finalize()
1269 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1270 call bdry_mask%free()
1274 do i = 1, this%bcs_vel%size()
1275 bci => this%bcs_vel%get(i)
1276 select type (
bc => bci)
1278 if (
bc%is_moving)
then
1280 call bdry_mask%init_from_components(this%c_Xh, 12.0_rp)
1283 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1285 call bdry_mask%mark_facets(bci%marked_facet)
1286 call bdry_mask%finalize()
1287 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1288 call bdry_mask%free()
1290 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1291 call bdry_mask%mark_facets(bci%marked_facet)
1292 call bdry_mask%finalize()
1293 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1294 call bdry_mask%free()
1295 type is (symmetry_t)
1296 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1297 call bdry_mask%mark_facets(bci%marked_facet)
1298 call bdry_mask%finalize()
1299 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1300 call bdry_mask%free()
1302 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1303 call bdry_mask%mark_facets(bci%marked_facet)
1304 call bdry_mask%finalize()
1305 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1306 call bdry_mask%free()
1307 type is (shear_stress_t)
1308 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1309 call bdry_mask%mark_facets(bci%marked_facet)
1310 call bdry_mask%finalize()
1311 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1312 call bdry_mask%free()
1313 type is (wall_model_bc_t)
1314 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1315 call bdry_mask%mark_facets(bci%marked_facet)
1316 call bdry_mask%finalize()
1317 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1318 call bdry_mask%free()
1320 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1321 call bdry_mask%mark_facets(bci%marked_facet)
1322 call bdry_mask%finalize()
1323 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1324 call bdry_mask%free()
1329 call bdry_file%init(
'bdry.fld')
1330 call bdry_file%write(bdry_field)
1332 call neko_scratch_registry%relinquish_field(temp_index)