247 subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
248 class(fluid_pnpn_t),
target,
intent(inout) :: this
249 type(
mesh_t),
target,
intent(inout) :: msh
250 integer,
intent(in) :: lx
251 type(json_file),
target,
intent(inout) :: params
252 type(
user_t),
target,
intent(in) :: user
253 type(
chkp_t),
target,
intent(inout) :: chkp
254 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
256 class(
bc_t),
pointer :: bc_i, vel_bc
257 real(kind=
rp) :: abs_tol
258 character(len=LOG_SIZE) :: log_buf
259 integer :: ierr, integer_val, solver_maxiter
260 character(len=:),
allocatable :: solver_type, precon_type
261 logical :: monitor, found
263 type(json_file) :: numerics_params, precon_params
268 call this%init_base(msh, lx, params, scheme,
user, .true.)
280 allocate(this%ext_bdf)
281 call this%ext_bdf%init(integer_val)
284 this%full_stress_formulation, .false.)
288 call this%c_Xh%check_cyclic()
290 if (this%full_stress_formulation)
then
292 call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
295 call pnpn_prs_res_stress_factory(this%prs_res)
298 call pnpn_vel_res_stress_factory(this%vel_res)
301 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
304 call pnpn_prs_res_factory(this%prs_res)
307 call pnpn_vel_res_factory(this%vel_res)
312 if (params%valid_path(
'case.fluid.nut_field'))
then
313 if (.not. this%full_stress_formulation)
then
314 call neko_error(
"You need to set full_stress_formulation to " // &
315 "true for the fluid to have a spatially varying " // &
318 call json_get(params,
'case.fluid.nut_field', this%nut_field_name)
320 this%nut_field_name =
""
324 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
328 call rhs_maker_sumab_fctry(this%sumab)
331 call rhs_maker_ext_fctry(this%makeabf)
334 call rhs_maker_bdf_fctry(this%makebdf)
337 call rhs_maker_oifs_fctry(this%makeoifs)
340 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
341 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
343 call this%p_res%init(dm_xh,
"p_res")
344 call this%u_res%init(dm_xh,
"u_res")
345 call this%v_res%init(dm_xh,
"v_res")
346 call this%w_res%init(dm_xh,
"w_res")
347 call this%abx1%init(dm_xh,
"abx1")
348 call this%aby1%init(dm_xh,
"aby1")
349 call this%abz1%init(dm_xh,
"abz1")
350 call this%abx2%init(dm_xh,
"abx2")
351 call this%aby2%init(dm_xh,
"aby2")
352 call this%abz2%init(dm_xh,
"abz2")
353 call this%advx%init(dm_xh,
"advx")
354 call this%advy%init(dm_xh,
"advy")
355 call this%advz%init(dm_xh,
"advz")
358 call this%du%init(this%dm_Xh,
'du')
359 call this%dv%init(this%dm_Xh,
'dv')
360 call this%dw%init(this%dm_Xh,
'dw')
361 call this%dp%init(this%dm_Xh,
'dp')
363 call neko_log%section(
"Fluid boundary conditions")
365 call this%setup_bcs(
user, params)
369 if (found)
call this%write_boundary_conditions()
372 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
373 this%pr_projection_activ_step)
375 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
376 this%vel_projection_activ_step)
382 if (params%valid_path(
'case.fluid.flow_rate_force'))
then
383 call this%vol_flow%init(this%dm_Xh, params)
387 call neko_log%section(
"Pressure solver")
390 'case.fluid.pressure_solver.max_iterations', &
392 call json_get(params,
'case.fluid.pressure_solver.type', solver_type)
393 call json_get(params,
'case.fluid.pressure_solver.preconditioner.type', &
396 'case.fluid.pressure_solver.preconditioner', precon_params)
398 'case.fluid.pressure_solver.absolute_tolerance', &
402 call neko_log%message(
'Type : ('// trim(solver_type) // &
403 ', ' // trim(precon_type) //
')')
404 write(log_buf,
'(A,ES13.6)')
'Abs tol :', abs_tol
407 call this%solver_factory(this%ksp_prs, this%dm_Xh%size(), &
408 solver_type, solver_maxiter, abs_tol, monitor)
409 call this%precon_factory_(this%pc_prs, this%ksp_prs, &
410 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, &
411 precon_type, precon_params)
416 call json_get(params,
'case.numerics', numerics_params)
417 call advection_factory(this%adv, numerics_params, this%c_Xh, &
418 this%ulag, this%vlag, this%wlag, &
419 chkp%dtlag, chkp%tlag, this%ext_bdf, &
424 call this%chkp%add_fluid(this%u, this%v, this%w, this%p)
426 this%chkp%abx1 => this%abx1
427 this%chkp%abx2 => this%abx2
428 this%chkp%aby1 => this%aby1
429 this%chkp%aby2 => this%aby2
430 this%chkp%abz1 => this%abz1
431 this%chkp%abz2 => this%abz2
432 call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
440 subroutine fluid_pnpn_restart(this, chkp)
441 class(fluid_pnpn_t),
target,
intent(inout) :: this
442 type(chkp_t),
intent(inout) :: chkp
443 real(kind=rp) :: dtlag(10), tlag(10)
444 type(field_t) :: u_temp, v_temp, w_temp
450 n = this%u%dof%size()
451 if (
allocated(chkp%previous_mesh%elements) .or. &
452 chkp%previous_Xh%lx .ne. this%Xh%lx)
then
453 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
454 c_xh => this%c_Xh, ulag => this%ulag, vlag => this%vlag, &
456 do concurrent(j = 1:n)
457 u%x(j,1,1,1) = u%x(j,1,1,1) * c_xh%mult(j,1,1,1)
458 v%x(j,1,1,1) = v%x(j,1,1,1) * c_xh%mult(j,1,1,1)
459 w%x(j,1,1,1) = w%x(j,1,1,1) * c_xh%mult(j,1,1,1)
460 p%x(j,1,1,1) = p%x(j,1,1,1) * c_xh%mult(j,1,1,1)
462 do i = 1, this%ulag%size()
463 do concurrent(j = 1:n)
464 ulag%lf(i)%x(j,1,1,1) = ulag%lf(i)%x(j,1,1,1) &
466 vlag%lf(i)%x(j,1,1,1) = vlag%lf(i)%x(j,1,1,1) &
468 wlag%lf(i)%x(j,1,1,1) = wlag%lf(i)%x(j,1,1,1) &
475 if (neko_bcknd_device .eq. 1)
then
476 associate(u => this%u, v => this%v, w => this%w, &
477 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
479 call device_memcpy(u%x, u%x_d, u%dof%size(), &
480 host_to_device, sync = .false.)
481 call device_memcpy(v%x, v%x_d, v%dof%size(), &
482 host_to_device, sync = .false.)
483 call device_memcpy(w%x, w%x_d, w%dof%size(), &
484 host_to_device, sync = .false.)
485 call device_memcpy(p%x, p%x_d, p%dof%size(), &
486 host_to_device, sync = .false.)
487 call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
488 u%dof%size(), host_to_device, sync = .false.)
489 call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
490 u%dof%size(), host_to_device, sync = .false.)
492 call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
493 v%dof%size(), host_to_device, sync = .false.)
494 call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
495 v%dof%size(), host_to_device, sync = .false.)
497 call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
498 w%dof%size(), host_to_device, sync = .false.)
499 call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
500 w%dof%size(), host_to_device, sync = .false.)
501 call device_memcpy(this%abx1%x, this%abx1%x_d, &
502 w%dof%size(), host_to_device, sync = .false.)
503 call device_memcpy(this%abx2%x, this%abx2%x_d, &
504 w%dof%size(), host_to_device, sync = .false.)
505 call device_memcpy(this%aby1%x, this%aby1%x_d, &
506 w%dof%size(), host_to_device, sync = .false.)
507 call device_memcpy(this%aby2%x, this%aby2%x_d, &
508 w%dof%size(), host_to_device, sync = .false.)
509 call device_memcpy(this%abz1%x, this%abz1%x_d, &
510 w%dof%size(), host_to_device, sync = .false.)
511 call device_memcpy(this%abz2%x, this%abz2%x_d, &
512 w%dof%size(), host_to_device, sync = .false.)
513 call device_memcpy(this%advx%x, this%advx%x_d, &
514 w%dof%size(), host_to_device, sync = .false.)
515 call device_memcpy(this%advy%x, this%advy%x_d, &
516 w%dof%size(), host_to_device, sync = .false.)
517 call device_memcpy(this%advz%x, this%advz%x_d, &
518 w%dof%size(), host_to_device, sync = .false.)
525 if (
allocated(chkp%previous_mesh%elements) &
526 .or. chkp%previous_Xh%lx .ne. this%Xh%lx)
then
528 call rotate_cyc(this%u%x, this%v%x, this%w%x, 1, this%c_Xh)
529 call this%gs_Xh%op(this%u, gs_op_add)
530 call this%gs_Xh%op(this%v, gs_op_add)
531 call this%gs_Xh%op(this%w, gs_op_add)
532 call this%gs_Xh%op(this%p, gs_op_add)
533 call rotate_cyc(this%u%x, this%v%x, this%w%x, 0, this%c_Xh)
535 do i = 1, this%ulag%size()
536 call rotate_cyc(this%ulag%lf(i)%x, this%vlag%lf(i)%x, &
537 this%wlag%lf(i)%x, 1, this%c_Xh)
538 call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
539 call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
540 call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
541 call rotate_cyc(this%ulag%lf(i)%x, this%vlag%lf(i)%x, &
542 this%wlag%lf(i)%x, 0, this%c_Xh)
648 subroutine fluid_pnpn_step(this, time, dt_controller)
649 class(fluid_pnpn_t),
target,
intent(inout) :: this
650 type(time_state_t),
intent(in) :: time
651 type(time_step_controller_t),
intent(in) :: dt_controller
655 type(ksp_monitor_t) :: ksp_results(4)
657 type(file_t) :: dump_file
658 class(bc_t),
pointer :: bc_i
659 type(non_normal_t),
pointer :: bc_j
661 if (this%freeze)
return
663 n = this%dm_Xh%size()
665 call profiler_start_region(
'Fluid', 1)
666 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
667 u_e => this%u_e, v_e => this%v_e, w_e => this%w_e, &
668 du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
669 u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
670 p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
672 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
673 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
674 msh => this%msh, prs_res => this%prs_res, &
675 source_term => this%source_term, vel_res => this%vel_res, &
676 sumab => this%sumab, makeoifs => this%makeoifs, &
677 makeabf => this%makeabf, makebdf => this%makebdf, &
678 vel_projection_dim => this%vel_projection_dim, &
679 pr_projection_dim => this%pr_projection_dim, &
681 rho => this%rho, mu_tot => this%mu_tot, &
682 f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
683 t => time%t, tstep => time%tstep, dt => time%dt, &
684 ext_bdf => this%ext_bdf, event => glb_cmd_event)
687 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
688 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
691 call this%source_term%compute(time)
694 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
695 this%dm_Xh%size(), time, strong = .false.)
699 call this%adv%compute(u, v, w, &
700 this%advx, this%advy, this%advz, &
701 xh, this%c_Xh, dm_xh%size(), dt)
708 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
709 this%abx2, this%aby2, this%abz2, &
710 f_x%x, f_y%x, f_z%x, &
711 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
714 call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
715 f_x%x, f_y%x, f_z%x, &
716 rho%x(1,1,1,1), dt, n)
719 call this%adv%compute(u, v, w, &
721 xh, this%c_Xh, dm_xh%size())
727 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
728 this%abx2, this%aby2, this%abz2, &
729 f_x%x, f_y%x, f_z%x, &
730 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
733 call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
734 u, v, w, c_xh%B, rho%x(1,1,1,1), dt, &
735 ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
742 call this%bc_apply_vel(time, strong = .true.)
743 call this%bc_apply_prs(time)
746 call this%update_material_properties(time)
749 call profiler_start_region(
'Pressure_residual', 18)
751 call prs_res%compute(p, p_res,&
756 this%bc_prs_surface, this%bc_sym_surface,&
757 ax_prs, ext_bdf%diffusion_coeffs(1), dt, &
761 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
762 call device_ortho(p_res%x_d, this%glb_n_points, n)
763 else if (.not. this%prs_dirichlet)
then
764 call ortho(p_res%x, this%glb_n_points, n)
767 call gs_xh%op(p_res, gs_op_add, event)
768 call device_event_sync(event)
771 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
774 call profiler_end_region(
'Pressure_residual', 18)
777 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
780 call this%pc_prs%update()
782 call profiler_start_region(
'Pressure_solve', 3)
786 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
787 this%bclst_dp, gs_xh)
788 ksp_results(1)%name =
'Pressure'
791 call profiler_end_region(
'Pressure_solve', 3)
793 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
794 this%bclst_dp, gs_xh, n, tstep, dt_controller)
797 call field_add2(p, dp, n)
798 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1)
then
799 call device_ortho(p%x_d, this%glb_n_points, n)
800 else if (.not. this%prs_dirichlet)
then
801 call ortho(p%x, this%glb_n_points, n)
805 call profiler_start_region(
'Velocity_residual', 19)
806 call vel_res%compute(ax_vel, u, v, w, &
807 u_res, v_res, w_res, &
811 mu_tot, rho, ext_bdf%diffusion_coeffs(1), &
814 call rotate_cyc(u_res%x, v_res%x, w_res%x, 1, c_xh)
815 call gs_xh%op(u_res, gs_op_add, event)
816 call device_event_sync(event)
817 call gs_xh%op(v_res, gs_op_add, event)
818 call device_event_sync(event)
819 call gs_xh%op(w_res, gs_op_add, event)
820 call device_event_sync(event)
821 call rotate_cyc(u_res%x, v_res%x, w_res%x, 0, c_xh)
824 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
827 call profiler_end_region(
'Velocity_residual', 19)
829 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
830 tstep, c_xh, n, dt_controller,
'Velocity')
832 call this%pc_vel%update()
834 call profiler_start_region(
"Velocity_solve", 4)
835 ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
836 u_res%x, v_res%x, w_res%x, n, c_xh, &
837 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
838 this%ksp_vel%max_iter)
839 call profiler_end_region(
"Velocity_solve", 4)
840 if (this%full_stress_formulation)
then
841 ksp_results(2)%name =
'Momentum'
843 ksp_results(2)%name =
'X-Velocity'
844 ksp_results(3)%name =
'Y-Velocity'
845 ksp_results(4)%name =
'Z-Velocity'
848 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
849 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
852 if (neko_bcknd_device .eq. 1)
then
853 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
854 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
856 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
859 if (this%forced_flow_rate)
then
861 call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
862 c_xh, gs_xh, ext_bdf, rho%x(1,1,1,1), mu_tot, &
863 dt, this%bclst_dp, this%bclst_du, this%bclst_dv, &
864 this%bclst_dw, this%bclst_vel_res, ax_vel, ax_prs, this%ksp_prs, &
865 this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
866 this%ksp_vel%max_iter)
869 call fluid_step_info(time, ksp_results, &
870 this%full_stress_formulation, this%strict_convergence, &
871 this%allow_stabilization)
874 call profiler_end_region(
'Fluid', 1)
879 subroutine fluid_pnpn_setup_bcs(this, user, params)
880 class(fluid_pnpn_t),
target,
intent(inout) :: this
881 type(user_t),
target,
intent(in) :: user
882 type(json_file),
intent(inout) :: params
883 integer :: i, n_bcs, zone_index, j, zone_size, global_zone_size, ierr
884 class(bc_t),
pointer :: bc_i
885 type(json_core) :: core
886 type(json_value),
pointer :: bc_object
887 type(json_file) :: bc_subdict
890 logical,
allocatable :: marked_zones(:)
891 integer,
allocatable :: zone_indices(:)
894 call this%bclst_vel_res%init()
895 call this%bclst_du%init()
896 call this%bclst_dv%init()
897 call this%bclst_dw%init()
898 call this%bclst_dp%init()
900 call this%bc_vel_res%init_from_components(this%c_Xh)
901 call this%bc_du%init_from_components(this%c_Xh)
902 call this%bc_dv%init_from_components(this%c_Xh)
903 call this%bc_dw%init_from_components(this%c_Xh)
904 call this%bc_dp%init_from_components(this%c_Xh)
907 call this%bc_prs_surface%init_from_components(this%c_Xh)
908 call this%bc_sym_surface%init_from_components(this%c_Xh)
911 if (params%valid_path(
'case.fluid.boundary_conditions'))
then
912 call params%info(
'case.fluid.boundary_conditions', n_children = n_bcs)
913 call params%get_core(core)
914 call params%get(
'case.fluid.boundary_conditions', bc_object, found)
919 call this%bcs_vel%init(n_bcs)
921 allocate(marked_zones(
size(this%msh%labeled_zones)))
922 marked_zones = .false.
926 call json_extract_item(core, bc_object, i, bc_subdict)
928 call json_get_or_lookup(bc_subdict,
"zone_indices", zone_indices)
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)
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 ", &
942 "is not in the mesh. Check fluid boundary condition ", &
947 if (marked_zones(zone_indices(j)))
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 "fluid and make sure that each zone index appears only ", &
953 "in a single boundary condition."
956 marked_zones(zone_indices(j)) = .true.
961 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
965 if (
associated(bc_i))
then
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)
984 call this%bcs_vel%append(bc_i)
986 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
987 type is (non_normal_t)
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)
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)
1001 call this%bcs_vel%append(bc_i)
1002 type is (wall_model_bc_t)
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)
1009 call this%bcs_vel%append(bc_i)
1016 if (bc_i%strong)
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)
1022 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
1025 call this%bcs_vel%append(bc_i)
1031 do i = 1,
size(this%msh%labeled_zones)
1032 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1033 (.not. marked_zones(i)))
then
1034 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
1035 "No fluid boundary condition assigned to zone ", i
1043 call this%bcs_prs%init(n_bcs)
1047 call json_extract_item(core, bc_object, i, bc_subdict)
1049 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
1053 if (
associated(bc_i))
then
1054 call this%bcs_prs%append(bc_i)
1057 if (bc_i%strong)
then
1058 call this%bc_dp%mark_facets(bc_i%marked_facet)
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!")
1074 call this%bcs_vel%init()
1075 call this%bcs_prs%init()
1079 call this%bc_prs_surface%finalize()
1080 call this%bc_sym_surface%finalize()
1082 call this%bc_vel_res%finalize()
1083 call this%bc_du%finalize()
1084 call this%bc_dv%finalize()
1085 call this%bc_dw%finalize()
1086 call this%bc_dp%finalize()
1088 call this%bclst_vel_res%append(this%bc_vel_res)
1089 call this%bclst_du%append(this%bc_du)
1090 call this%bclst_dv%append(this%bc_dv)
1091 call this%bclst_dw%append(this%bc_dw)
1092 call this%bclst_dp%append(this%bc_dp)
1095 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1096 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1097 mpi_logical, mpi_lor, neko_comm)
1100 if (
allocated(marked_zones))
then
1101 deallocate(marked_zones)
1104 if (
allocated(zone_indices))
then
1105 deallocate(zone_indices)
1111 subroutine fluid_pnpn_write_boundary_conditions(this)
1117 class(fluid_pnpn_t),
target,
intent(inout) :: this
1118 type(dirichlet_t) :: bdry_mask
1119 type(field_t),
pointer :: bdry_field
1120 type(file_t) :: bdry_file
1121 integer :: temp_index, i
1122 class(bc_t),
pointer :: bci
1123 character(len=LOG_SIZE) :: log_buf
1125 write(log_buf,
'(A)')
'Marking using integer keys in bdry0.f00000'
1126 call neko_log%message(log_buf)
1127 write(log_buf,
'(A)')
'Condition-value pairs: '
1128 call neko_log%message(log_buf)
1129 write(log_buf,
'(A)')
' no_slip = 1'
1130 call neko_log%message(log_buf)
1131 write(log_buf,
'(A)')
' velocity_value = 2'
1132 call neko_log%message(log_buf)
1133 write(log_buf,
'(A)')
' outflow, normal_outflow (+dong) = 3'
1134 call neko_log%message(log_buf)
1135 write(log_buf,
'(A)')
' symmetry = 4'
1136 call neko_log%message(log_buf)
1137 write(log_buf,
'(A)')
' periodic = 6'
1138 call neko_log%message(log_buf)
1139 write(log_buf,
'(A)')
' user_velocity = 7'
1140 call neko_log%message(log_buf)
1141 write(log_buf,
'(A)')
' user_pressure = 8'
1142 call neko_log%message(log_buf)
1143 write(log_buf,
'(A)')
' shear_stress = 9'
1144 call neko_log%message(log_buf)
1145 write(log_buf,
'(A)')
' wall_modelling = 10'
1146 call neko_log%message(log_buf)
1147 write(log_buf,
'(A)')
' blasius_profile = 11'
1148 call neko_log%message(log_buf)
1150 call neko_scratch_registry%request_field(bdry_field, temp_index, .true.)
1154 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1155 call bdry_mask%mark_zone(this%msh%periodic)
1156 call bdry_mask%finalize()
1157 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1158 call bdry_mask%free()
1160 do i = 1, this%bcs_prs%size()
1161 bci => this%bcs_prs%get(i)
1162 select type (
bc => bci)
1163 type is (zero_dirichlet_t)
1164 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1165 call bdry_mask%mark_facets(bci%marked_facet)
1166 call bdry_mask%finalize()
1167 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1168 call bdry_mask%free()
1170 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1171 call bdry_mask%mark_facets(bci%marked_facet)
1172 call bdry_mask%finalize()
1173 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1174 call bdry_mask%free()
1176 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1177 call bdry_mask%mark_facets(bci%marked_facet)
1178 call bdry_mask%finalize()
1179 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1180 call bdry_mask%free()
1184 do i = 1, this%bcs_vel%size()
1185 bci => this%bcs_vel%get(i)
1186 select type (
bc => bci)
1187 type is (zero_dirichlet_t)
1188 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1189 call bdry_mask%mark_facets(bci%marked_facet)
1190 call bdry_mask%finalize()
1191 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1192 call bdry_mask%free()
1194 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1195 call bdry_mask%mark_facets(bci%marked_facet)
1196 call bdry_mask%finalize()
1197 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1198 call bdry_mask%free()
1199 type is (symmetry_t)
1200 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1201 call bdry_mask%mark_facets(bci%marked_facet)
1202 call bdry_mask%finalize()
1203 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1204 call bdry_mask%free()
1206 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1207 call bdry_mask%mark_facets(bci%marked_facet)
1208 call bdry_mask%finalize()
1209 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1210 call bdry_mask%free()
1211 type is (shear_stress_t)
1212 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1213 call bdry_mask%mark_facets(bci%marked_facet)
1214 call bdry_mask%finalize()
1215 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1216 call bdry_mask%free()
1217 type is (wall_model_bc_t)
1218 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1219 call bdry_mask%mark_facets(bci%marked_facet)
1220 call bdry_mask%finalize()
1221 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1222 call bdry_mask%free()
1224 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1225 call bdry_mask%mark_facets(bci%marked_facet)
1226 call bdry_mask%finalize()
1227 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1228 call bdry_mask%free()
1233 call bdry_file%init(
'bdry.fld')
1234 call bdry_file%write(bdry_field)
1236 call neko_scratch_registry%relinquish_field(temp_index)