254 subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
255 class(fluid_pnpn_t),
target,
intent(inout) :: this
256 type(
mesh_t),
target,
intent(inout) :: msh
257 integer,
intent(in) :: lx
258 type(json_file),
target,
intent(inout) :: params
259 type(
user_t),
target,
intent(in) :: user
260 type(
chkp_t),
target,
intent(inout) :: chkp
261 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
263 class(
bc_t),
pointer :: bc_i, vel_bc
264 real(kind=
rp) :: abs_tol
265 character(len=LOG_SIZE) :: log_buf
266 integer :: ierr, integer_val, solver_maxiter
267 character(len=:),
allocatable :: solver_type, precon_type
268 logical :: monitor, found
270 type(json_file) :: numerics_params, precon_params
275 call this%init_base(msh, lx, params, scheme,
user, .true.)
287 allocate(this%ext_bdf)
288 call this%ext_bdf%init(integer_val)
291 this%full_stress_formulation, .false.)
295 call this%c_Xh%generate_cyclic_bc()
297 if (this%full_stress_formulation)
then
299 call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
302 call pnpn_prs_res_stress_factory(this%prs_res)
305 call pnpn_vel_res_stress_factory(this%vel_res)
308 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
311 call pnpn_prs_res_factory(this%prs_res)
314 call pnpn_vel_res_factory(this%vel_res)
319 if (params%valid_path(
'case.fluid.nut_field'))
then
320 if (.not. this%full_stress_formulation)
then
321 call neko_error(
"You need to set full_stress_formulation to " // &
322 "true for the fluid to have a spatially varying " // &
325 call json_get(params,
'case.fluid.nut_field', this%nut_field_name)
327 this%nut_field_name =
""
331 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
335 call rhs_maker_sumab_fctry(this%sumab)
338 call rhs_maker_ext_fctry(this%makeabf)
341 call rhs_maker_bdf_fctry(this%makebdf)
344 call rhs_maker_oifs_fctry(this%makeoifs)
347 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
348 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
350 call this%p_res%init(dm_xh,
"p_res")
351 call this%u_res%init(dm_xh,
"u_res")
352 call this%v_res%init(dm_xh,
"v_res")
353 call this%w_res%init(dm_xh,
"w_res")
354 call this%abx1%init(dm_xh,
"abx1")
355 call this%aby1%init(dm_xh,
"aby1")
356 call this%abz1%init(dm_xh,
"abz1")
357 call this%abx2%init(dm_xh,
"abx2")
358 call this%aby2%init(dm_xh,
"aby2")
359 call this%abz2%init(dm_xh,
"abz2")
360 call this%advx%init(dm_xh,
"advx")
361 call this%advy%init(dm_xh,
"advy")
362 call this%advz%init(dm_xh,
"advz")
365 call this%du%init(this%dm_Xh,
'du')
366 call this%dv%init(this%dm_Xh,
'dv')
367 call this%dw%init(this%dm_Xh,
'dw')
368 call this%dp%init(this%dm_Xh,
'dp')
370 call this%ale%init(this%c_Xh, params,
user)
372 call neko_log%section(
"Fluid boundary conditions")
374 call this%setup_bcs(
user, params)
378 if (found)
call this%write_boundary_conditions()
381 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
382 this%pr_projection_activ_step, &
383 this%pr_projection_reorthogonalize_basis)
385 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
386 this%vel_projection_activ_step)
392 if (params%valid_path(
'case.fluid.flow_rate_force'))
then
393 call this%vol_flow%init(this%dm_Xh, params)
397 call neko_log%section(
"Pressure solver")
400 'case.fluid.pressure_solver.max_iterations', &
402 call json_get(params,
'case.fluid.pressure_solver.type', solver_type)
403 call json_get(params,
'case.fluid.pressure_solver.preconditioner.type', &
406 'case.fluid.pressure_solver.preconditioner', precon_params)
408 'case.fluid.pressure_solver.absolute_tolerance', &
412 call neko_log%message(
'Type : ('// trim(solver_type) // &
413 ', ' // trim(precon_type) //
')')
414 write(log_buf,
'(A,ES13.6)')
'Abs tol :', abs_tol
417 call this%solver_factory(this%ksp_prs, this%dm_Xh%size(), &
418 solver_type, solver_maxiter, abs_tol, monitor)
419 call this%precon_factory_(this%pc_prs, this%ksp_prs, &
420 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, &
421 precon_type, precon_params)
426 call json_get(params,
'case.numerics', numerics_params)
427 call advection_factory(this%adv, numerics_params, this%c_Xh, &
428 this%ulag, this%vlag, this%wlag, &
429 chkp%dtlag, chkp%tlag, this%ext_bdf, &
434 call this%chkp%add_fluid(this%u, this%v, this%w, this%p)
436 this%chkp%abx1 => this%abx1
437 this%chkp%abx2 => this%abx2
438 this%chkp%aby1 => this%aby1
439 this%chkp%aby2 => this%aby2
440 this%chkp%abz1 => this%abz1
441 this%chkp%abz2 => this%abz2
442 call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
445 if (this%ale%active)
then
446 call this%chkp%add_ale(this%c_Xh%dof%x, this%c_Xh%dof%y, &
448 this%c_Xh%Blag, this%c_Xh%Blaglag, &
449 this%ale%wm_x, this%ale%wm_y, this%ale%wm_z, &
450 this%ale%wm_x_lag, this%ale%wm_y_lag, &
452 this%ale%global_pivot_pos, &
453 this%ale%global_pivot_vel_lag, &
454 this%ale%global_basis_pos, &
455 this%ale%global_basis_vel_lag)
460 this%schwarz_iterations, 0)
466 subroutine fluid_pnpn_restart(this, chkp)
467 class(fluid_pnpn_t),
target,
intent(inout) :: this
468 type(chkp_t),
intent(inout) :: chkp
469 real(kind=rp) :: dtlag(10), tlag(10)
470 type(field_t) :: u_temp, v_temp, w_temp
476 n = this%u%dof%size()
477 if (
allocated(chkp%previous_mesh%elements) .or. &
478 chkp%previous_Xh%lx .ne. this%Xh%lx)
then
479 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
480 c_xh => this%c_Xh, ulag => this%ulag, vlag => this%vlag, &
482 do concurrent(j = 1:n)
483 u%x(j,1,1,1) = u%x(j,1,1,1) * c_xh%mult(j,1,1,1)
484 v%x(j,1,1,1) = v%x(j,1,1,1) * c_xh%mult(j,1,1,1)
485 w%x(j,1,1,1) = w%x(j,1,1,1) * c_xh%mult(j,1,1,1)
486 p%x(j,1,1,1) = p%x(j,1,1,1) * c_xh%mult(j,1,1,1)
488 do i = 1, this%ulag%size()
489 do concurrent(j = 1:n)
490 ulag%lf(i)%x(j,1,1,1) = ulag%lf(i)%x(j,1,1,1) &
492 vlag%lf(i)%x(j,1,1,1) = vlag%lf(i)%x(j,1,1,1) &
494 wlag%lf(i)%x(j,1,1,1) = wlag%lf(i)%x(j,1,1,1) &
501 if (neko_bcknd_device .eq. 1)
then
502 associate(u => this%u, v => this%v, w => this%w, &
503 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
505 call device_memcpy(u%x, u%x_d, u%dof%size(), &
506 host_to_device, sync = .false.)
507 call device_memcpy(v%x, v%x_d, v%dof%size(), &
508 host_to_device, sync = .false.)
509 call device_memcpy(w%x, w%x_d, w%dof%size(), &
510 host_to_device, sync = .false.)
511 call device_memcpy(p%x, p%x_d, p%dof%size(), &
512 host_to_device, sync = .false.)
513 call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
514 u%dof%size(), host_to_device, sync = .false.)
515 call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
516 u%dof%size(), host_to_device, sync = .false.)
518 call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
519 v%dof%size(), host_to_device, sync = .false.)
520 call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
521 v%dof%size(), host_to_device, sync = .false.)
523 call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
524 w%dof%size(), host_to_device, sync = .false.)
525 call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
526 w%dof%size(), host_to_device, sync = .false.)
527 call device_memcpy(this%abx1%x, this%abx1%x_d, &
528 w%dof%size(), host_to_device, sync = .false.)
529 call device_memcpy(this%abx2%x, this%abx2%x_d, &
530 w%dof%size(), host_to_device, sync = .false.)
531 call device_memcpy(this%aby1%x, this%aby1%x_d, &
532 w%dof%size(), host_to_device, sync = .false.)
533 call device_memcpy(this%aby2%x, this%aby2%x_d, &
534 w%dof%size(), host_to_device, sync = .false.)
535 call device_memcpy(this%abz1%x, this%abz1%x_d, &
536 w%dof%size(), host_to_device, sync = .false.)
537 call device_memcpy(this%abz2%x, this%abz2%x_d, &
538 w%dof%size(), host_to_device, sync = .false.)
539 call device_memcpy(this%advx%x, this%advx%x_d, &
540 w%dof%size(), host_to_device, sync = .false.)
541 call device_memcpy(this%advy%x, this%advy%x_d, &
542 w%dof%size(), host_to_device, sync = .false.)
543 call device_memcpy(this%advz%x, this%advz%x_d, &
544 w%dof%size(), host_to_device, sync = .false.)
551 if (
allocated(chkp%previous_mesh%elements) &
552 .or. chkp%previous_Xh%lx .ne. this%Xh%lx)
then
554 call rotate_cyc(this%u, this%v, this%w, 1, this%c_Xh)
555 call this%gs_Xh%op(this%u, gs_op_add)
556 call this%gs_Xh%op(this%v, gs_op_add)
557 call this%gs_Xh%op(this%w, gs_op_add)
558 call this%gs_Xh%op(this%p, gs_op_add)
559 call rotate_cyc(this%u, this%v, this%w, 0, this%c_Xh)
561 do i = 1, this%ulag%size()
562 call rotate_cyc(this%ulag%lf(i), this%vlag%lf(i), &
563 this%wlag%lf(i), 1, this%c_Xh)
564 call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
565 call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
566 call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
567 call rotate_cyc(this%ulag%lf(i), this%vlag%lf(i), &
568 this%wlag%lf(i), 0, this%c_Xh)
572 call this%ale%set_coef_restart(this%c_Xh, this%adv, chkp%t)
679 subroutine fluid_pnpn_step(this, time, dt_controller)
680 class(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
686 type(ksp_monitor_t) :: ksp_results(4)
689 type(file_t) :: dump_file
690 class(bc_t),
pointer :: bc_i
691 type(non_normal_t),
pointer :: bc_j
693 if (this%freeze)
return
695 n = this%dm_Xh%size()
697 call profiler_start_region(
'Fluid', 1)
698 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
699 u_e => this%u_e, v_e => this%v_e, w_e => this%w_e, &
700 du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
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, makeoifs => this%makeoifs, &
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_tot => this%mu_tot, &
714 f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
715 t => time%t, tstep => time%tstep, dt => time%dt, &
716 ext_bdf => this%ext_bdf, event => glb_cmd_event, &
720 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
721 ulag, vlag, wlag, ext_bdf%advection_coeffs%x, ext_bdf%nadv)
724 call this%source_term%compute(time)
727 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
728 this%dm_Xh%size(), time, strong = .false.)
730 if (this%ale%active)
then
732 call neko_error(
"ALE is not yet supported " // &
733 "with OIFS time integration.")
736 call this%adv%compute_ale(u, v, w, &
737 ale%wm_x, ale%wm_y, ale%wm_z, &
739 xh, c_xh, dm_xh%size())
745 call this%adv%compute(u, v, w, &
746 this%advx, this%advy, this%advz, &
747 xh, this%c_Xh, dm_xh%size(), dt)
754 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
755 this%abx2, this%aby2, this%abz2, &
756 f_x%x, f_y%x, f_z%x, &
757 rho%x(1,1,1,1), ext_bdf%advection_coeffs%x, n)
760 call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
761 f_x%x, f_y%x, f_z%x, &
762 rho%x(1,1,1,1), dt, n)
765 call this%adv%compute(u, v, w, &
767 xh, this%c_Xh, dm_xh%size())
774 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
775 this%abx2, this%aby2, this%abz2, &
776 f_x%x, f_y%x, f_z%x, &
777 rho%x(1,1,1,1), ext_bdf%advection_coeffs%x, n)
783 call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
784 u, v, w, c_xh%B, rho%x(1,1,1,1), dt, &
785 ext_bdf%diffusion_coeffs%x, ext_bdf%ndiff, n, &
786 c_xh%Blag, c_xh%Blaglag)
790 if (this%ale%active)
then
792 call this%ale%advance_mesh(c_xh, time, ext_bdf%nadv)
794 call profiler_start_region(
'ALE recompute metrics')
796 call c_xh%recompute_metrics()
799 call this%adv%recompute_metrics(c_xh, .true.)
800 call profiler_end_region(
'ALE recompute metrics')
808 call this%update_material_properties(time)
810 do iter = 1, 1 + this%schwarz_iterations
812 call this%bc_apply_vel(time, strong = .true.)
813 call this%bc_apply_prs(time)
816 call profiler_start_region(
'Pressure_residual', 18)
817 call prs_res%compute(p, p_res,&
822 this%bc_prs_surface, this%bc_sym_surface,&
823 ax_prs, ext_bdf%diffusion_coeffs%x(1), dt, &
828 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
829 call device_ortho(p_res%x_d, this%glb_n_points, n)
830 else if (.not. this%prs_dirichlet)
then
831 call ortho(p_res%x, this%glb_n_points, n)
834 call gs_xh%op(p_res, gs_op_add, event)
835 call device_event_sync(event)
838 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
841 call profiler_end_region(
'Pressure_residual', 18)
845 if (iter .eq. 1)
then
846 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
847 ax = ax_prs, gs_h = gs_xh, bclst = this%bclst_dp, &
851 call this%pc_prs%update()
853 call profiler_start_region(
'Pressure_solve', 3)
857 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
858 this%bclst_dp, gs_xh)
859 ksp_results(1)%name =
'Pressure'
862 call profiler_end_region(
'Pressure_solve', 3)
864 if (iter .eq. 1)
then
865 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
866 this%bclst_dp, gs_xh, n, tstep, dt_controller)
870 call field_add2(p, dp, n)
871 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
872 call device_ortho(p%x_d, this%glb_n_points, n)
873 else if (.not. this%prs_dirichlet)
then
874 call ortho(p%x, this%glb_n_points, n)
878 call profiler_start_region(
'Velocity_residual', 19)
879 call vel_res%compute(ax_vel, u, v, w, &
880 u_res, v_res, w_res, &
884 mu_tot, rho, ext_bdf%diffusion_coeffs%x(1), &
887 call rotate_cyc(u_res, v_res, w_res, 1, c_xh)
888 call gs_xh%op(u_res, gs_op_add, event)
889 call device_event_sync(event)
890 call gs_xh%op(v_res, gs_op_add, event)
891 call device_event_sync(event)
892 call gs_xh%op(w_res, gs_op_add, event)
893 call device_event_sync(event)
894 call rotate_cyc(u_res, v_res, w_res, 0, c_xh)
897 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
900 call profiler_end_region(
'Velocity_residual', 19)
902 if (iter .eq. 1)
then
903 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
904 tstep, c_xh, n, dt_controller,
'Velocity')
907 call this%pc_vel%update()
909 call profiler_start_region(
"Velocity_solve", 4)
910 ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
911 u_res%x, v_res%x, w_res%x, n, c_xh, &
912 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
913 this%ksp_vel%max_iter)
914 call profiler_end_region(
"Velocity_solve", 4)
915 if (this%full_stress_formulation)
then
916 ksp_results(2)%name =
'Momentum'
918 ksp_results(2)%name =
'X-Velocity'
919 ksp_results(3)%name =
'Y-Velocity'
920 ksp_results(4)%name =
'Z-Velocity'
923 if (iter .eq. 1)
then
924 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
925 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
929 if (neko_bcknd_device .eq. 1)
then
930 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
931 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
933 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
936 call fluid_step_info(time, ksp_results, &
937 this%full_stress_formulation, this%strict_convergence, &
938 this%allow_stabilization, iter)
942 if (this%forced_flow_rate)
then
944 call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
945 c_xh, gs_xh, ext_bdf, rho%x(1,1,1,1), mu_tot, &
946 dt, time, this%bclst_dp, this%bclst_du, this%bclst_dv, &
947 this%bclst_dw, this%bclst_vel_res, ax_vel, ax_prs, this%ksp_prs, &
948 this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
949 this%ksp_vel%max_iter)
955 call this%ale%update_mesh_velocity(c_xh, time)
958 call profiler_end_region(
'Fluid', 1)
963 subroutine fluid_pnpn_setup_bcs(this, user, params)
964 class(fluid_pnpn_t),
target,
intent(inout) :: this
965 type(user_t),
target,
intent(in) :: user
966 type(json_file),
intent(inout) :: params
967 integer :: i, n_bcs, zone_index, j, zone_size, global_zone_size, ierr
968 class(bc_t),
pointer :: bc_i
969 type(json_core) :: core
970 type(json_value),
pointer :: bc_object
971 type(json_file) :: bc_subdict
972 logical :: ale_active_local, any_moving_wall, moving_
975 logical,
allocatable :: marked_zones(:)
976 integer,
allocatable :: zone_indices(:)
979 character(len=:),
allocatable :: bc_type_str
980 this%ale%has_moving_boundary = .false.
981 any_moving_wall = .false.
982 ale_active_local = .false.
983 call json_get_or_default(params, &
984 'case.fluid.ale.enabled', &
985 ale_active_local, .false.)
988 call this%bclst_vel_res%init()
989 call this%bclst_du%init()
990 call this%bclst_dv%init()
991 call this%bclst_dw%init()
992 call this%bclst_dp%init()
994 call this%bc_vel_res%init_from_components(this%c_Xh)
995 call this%bc_du%init_from_components(this%c_Xh)
996 call this%bc_dv%init_from_components(this%c_Xh)
997 call this%bc_dw%init_from_components(this%c_Xh)
998 call this%bc_dp%init_from_components(this%c_Xh)
1001 call this%bc_prs_surface%init_from_components(this%c_Xh)
1002 call this%bc_sym_surface%init_from_components(this%c_Xh)
1005 if (params%valid_path(
'case.fluid.boundary_conditions'))
then
1006 call params%info(
'case.fluid.boundary_conditions', n_children = n_bcs)
1007 call params%get_core(core)
1008 call params%get(
'case.fluid.boundary_conditions', bc_object, found)
1013 call this%bcs_vel%init(n_bcs)
1015 allocate(marked_zones(
size(this%msh%labeled_zones)))
1016 marked_zones = .false.
1020 call json_extract_item(core, bc_object, i, bc_subdict)
1022 call json_get_or_lookup(bc_subdict,
"zone_indices", zone_indices)
1025 call json_get(bc_subdict,
"type", bc_type_str)
1027 if (trim(bc_type_str) .eq.
"no_slip")
then
1028 call json_get_or_default(bc_subdict,
"moving", moving_, .false.)
1031 this%ale%has_moving_boundary = .true.
1037 do j = 1,
size(zone_indices)
1038 zone_size = this%msh%labeled_zones(zone_indices(j))%size
1039 call mpi_allreduce(zone_size, global_zone_size, 1, &
1040 mpi_integer, mpi_max, neko_comm, ierr)
1042 if (global_zone_size .eq. 0)
then
1043 write(error_unit,
'(A, A, I0, A, A, I0, A)')
"*** ERROR ***: ",&
1044 "Zone index ", zone_indices(j), &
1045 " is invalid as this zone has 0 size, meaning it ", &
1046 "is not in the mesh. Check fluid boundary condition ", &
1051 if (marked_zones(zone_indices(j)))
then
1052 write(error_unit,
'(A, A, I0, A, A, A, A)')
"*** ERROR ***: ", &
1053 "Zone with index ", zone_indices(j), &
1054 " has already been assigned a boundary condition. ", &
1055 "Please check your boundary_conditions entry for the ", &
1056 "fluid and make sure that each zone index appears only ", &
1057 "in a single boundary condition."
1060 marked_zones(zone_indices(j)) = .true.
1065 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
1069 if (
associated(bc_i))
then
1075 type is (symmetry_t)
1083 call this%bclst_vel_res%append(bc_i)
1084 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1085 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1086 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1088 call this%bcs_vel%append(bc_i)
1090 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
1091 type is (non_normal_t)
1094 call this%bclst_vel_res%append(bc_i)
1095 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1096 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1097 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1098 type is (shear_stress_t)
1100 call this%bclst_vel_res%append(bc_i%symmetry)
1101 call this%bclst_du%append(bc_i%symmetry%bc_x)
1102 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1103 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1105 call this%bcs_vel%append(bc_i)
1106 type is (wall_model_bc_t)
1108 call this%bclst_vel_res%append(bc_i%symmetry)
1109 call this%bclst_du%append(bc_i%symmetry%bc_x)
1110 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1111 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1113 call this%bcs_vel%append(bc_i)
1120 if (bc_i%strong)
then
1121 call this%bc_vel_res%mark_facets(bc_i%marked_facet)
1122 call this%bc_du%mark_facets(bc_i%marked_facet)
1123 call this%bc_dv%mark_facets(bc_i%marked_facet)
1124 call this%bc_dw%mark_facets(bc_i%marked_facet)
1126 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
1129 call this%bcs_vel%append(bc_i)
1134 if (this%ale%active .and. (.not. this%ale%has_moving_boundary))
then
1135 call neko_error(
"Case file error: ALE is active, " // &
1136 "but no moving wall was found. " // &
1137 "Use type = 'no_slip' with 'moving': true in case file.")
1141 do i = 1,
size(this%msh%labeled_zones)
1142 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1143 (.not. marked_zones(i)))
then
1144 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
1145 "No fluid boundary condition assigned to zone ", i
1153 call this%bcs_prs%init(n_bcs)
1157 call json_extract_item(core, bc_object, i, bc_subdict)
1159 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
1163 if (
associated(bc_i))
then
1164 call this%bcs_prs%append(bc_i)
1167 if (bc_i%strong)
then
1168 call this%bc_dp%mark_facets(bc_i%marked_facet)
1176 do i = 1,
size(this%msh%labeled_zones)
1177 if (this%msh%labeled_zones(i)%size .gt. 0)
then
1178 call neko_error(
"No boundary_conditions entry in the case file!")
1184 call this%bcs_vel%init()
1185 call this%bcs_prs%init()
1189 call this%bc_prs_surface%finalize()
1190 call this%bc_sym_surface%finalize()
1192 call this%bc_vel_res%finalize()
1193 call this%bc_du%finalize()
1194 call this%bc_dv%finalize()
1195 call this%bc_dw%finalize()
1196 call this%bc_dp%finalize()
1198 call this%bclst_vel_res%append(this%bc_vel_res)
1199 call this%bclst_du%append(this%bc_du)
1200 call this%bclst_dv%append(this%bc_dv)
1201 call this%bclst_dw%append(this%bc_dw)
1202 call this%bclst_dp%append(this%bc_dp)
1205 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1206 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1207 mpi_logical, mpi_lor, neko_comm)
1210 if (
allocated(marked_zones))
then
1211 deallocate(marked_zones)
1214 if (
allocated(zone_indices))
then
1215 deallocate(zone_indices)
1221 subroutine fluid_pnpn_write_boundary_conditions(this)
1228 class(fluid_pnpn_t),
target,
intent(inout) :: this
1229 type(dirichlet_t) :: bdry_mask
1230 type(field_t),
pointer :: bdry_field
1231 type(file_t) :: bdry_file
1232 integer :: temp_index, i
1233 class(bc_t),
pointer :: bci
1234 character(len=LOG_SIZE) :: log_buf
1236 write(log_buf,
'(A)')
'Marking using integer keys in bdry0.f00000'
1237 call neko_log%message(log_buf)
1238 write(log_buf,
'(A)')
'Condition-value pairs: '
1239 call neko_log%message(log_buf)
1240 write(log_buf,
'(A)')
' no_slip (stationary wall) = 1'
1241 call neko_log%message(log_buf)
1242 write(log_buf,
'(A)')
' velocity_value = 2'
1243 call neko_log%message(log_buf)
1244 write(log_buf,
'(A)')
' outflow, normal_outflow (+dong) = 3'
1245 call neko_log%message(log_buf)
1246 write(log_buf,
'(A)')
' symmetry = 4'
1247 call neko_log%message(log_buf)
1248 write(log_buf,
'(A)')
' periodic = 6'
1249 call neko_log%message(log_buf)
1250 write(log_buf,
'(A)')
' user_velocity = 7'
1251 call neko_log%message(log_buf)
1252 write(log_buf,
'(A)')
' user_pressure = 8'
1253 call neko_log%message(log_buf)
1254 write(log_buf,
'(A)')
' shear_stress = 9'
1255 call neko_log%message(log_buf)
1256 write(log_buf,
'(A)')
' wall_modelling = 10'
1257 call neko_log%message(log_buf)
1258 write(log_buf,
'(A)')
' blasius_profile = 11'
1259 call neko_log%message(log_buf)
1260 write(log_buf,
'(A)')
' no_slip (moving wall) = 12'
1261 call neko_log%message(log_buf)
1263 call neko_scratch_registry%request_field(bdry_field, temp_index, .true.)
1267 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1268 call bdry_mask%mark_zone(this%msh%periodic)
1269 call bdry_mask%finalize()
1270 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1271 call bdry_mask%free()
1273 do i = 1, this%bcs_prs%size()
1274 bci => this%bcs_prs%get(i)
1275 select type (
bc => bci)
1276 type is (zero_dirichlet_t)
1277 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1278 call bdry_mask%mark_facets(bci%marked_facet)
1279 call bdry_mask%finalize()
1280 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1281 call bdry_mask%free()
1283 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1284 call bdry_mask%mark_facets(bci%marked_facet)
1285 call bdry_mask%finalize()
1286 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1287 call bdry_mask%free()
1289 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1290 call bdry_mask%mark_facets(bci%marked_facet)
1291 call bdry_mask%finalize()
1292 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1293 call bdry_mask%free()
1297 do i = 1, this%bcs_vel%size()
1298 bci => this%bcs_vel%get(i)
1299 select type (
bc => bci)
1301 if (
bc%is_moving)
then
1303 call bdry_mask%init_from_components(this%c_Xh, 12.0_rp)
1306 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1308 call bdry_mask%mark_facets(bci%marked_facet)
1309 call bdry_mask%finalize()
1310 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1311 call bdry_mask%free()
1313 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1314 call bdry_mask%mark_facets(bci%marked_facet)
1315 call bdry_mask%finalize()
1316 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1317 call bdry_mask%free()
1318 type is (symmetry_t)
1319 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1320 call bdry_mask%mark_facets(bci%marked_facet)
1321 call bdry_mask%finalize()
1322 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1323 call bdry_mask%free()
1325 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1326 call bdry_mask%mark_facets(bci%marked_facet)
1327 call bdry_mask%finalize()
1328 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1329 call bdry_mask%free()
1330 type is (shear_stress_t)
1331 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1332 call bdry_mask%mark_facets(bci%marked_facet)
1333 call bdry_mask%finalize()
1334 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1335 call bdry_mask%free()
1336 type is (wall_model_bc_t)
1337 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1338 call bdry_mask%mark_facets(bci%marked_facet)
1339 call bdry_mask%finalize()
1340 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1341 call bdry_mask%free()
1343 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1344 call bdry_mask%mark_facets(bci%marked_facet)
1345 call bdry_mask%finalize()
1346 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1347 call bdry_mask%free()
1352 call bdry_file%init(
'bdry.fld')
1353 call bdry_file%write(bdry_field)
1355 call neko_scratch_registry%relinquish_field(temp_index)