245 subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
246 class(fluid_pnpn_t),
target,
intent(inout) :: this
247 type(
mesh_t),
target,
intent(inout) :: msh
248 integer,
intent(in) :: lx
249 type(json_file),
target,
intent(inout) :: params
250 type(
user_t),
target,
intent(in) :: user
251 type(
chkp_t),
target,
intent(inout) :: chkp
252 character(len=15),
parameter :: scheme =
'Modular (Pn/Pn)'
254 class(
bc_t),
pointer :: bc_i, vel_bc
255 real(kind=
rp) :: abs_tol
256 character(len=LOG_SIZE) :: log_buf
257 integer :: ierr, integer_val, solver_maxiter
258 character(len=:),
allocatable :: solver_type, precon_type
259 logical :: monitor, found
261 type(json_file) :: numerics_params, precon_params
266 call this%init_base(msh, lx, params, scheme,
user, .true.)
277 call json_get(params,
'case.numerics.time_order', integer_val)
278 allocate(this%ext_bdf)
279 call this%ext_bdf%init(integer_val)
282 this%full_stress_formulation, .false.)
286 call this%c_Xh%generate_cyclic_bc()
288 if (this%full_stress_formulation .eqv. .true.)
then
290 call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
293 call pnpn_prs_res_stress_factory(this%prs_res)
296 call pnpn_vel_res_stress_factory(this%vel_res)
299 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
302 call pnpn_prs_res_factory(this%prs_res)
305 call pnpn_vel_res_factory(this%vel_res)
310 if (params%valid_path(
'case.fluid.nut_field'))
then
311 if (this%full_stress_formulation .eqv. .false.)
then
312 call neko_error(
"You need to set full_stress_formulation to " // &
313 "true for the fluid to have a spatially varying " // &
316 call json_get(params,
'case.fluid.nut_field', this%nut_field_name)
318 this%nut_field_name =
""
322 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
326 call rhs_maker_sumab_fctry(this%sumab)
329 call rhs_maker_ext_fctry(this%makeabf)
332 call rhs_maker_bdf_fctry(this%makebdf)
335 call rhs_maker_oifs_fctry(this%makeoifs)
338 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
339 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
341 call this%p_res%init(dm_xh,
"p_res")
342 call this%u_res%init(dm_xh,
"u_res")
343 call this%v_res%init(dm_xh,
"v_res")
344 call this%w_res%init(dm_xh,
"w_res")
345 call this%abx1%init(dm_xh,
"abx1")
346 call this%aby1%init(dm_xh,
"aby1")
347 call this%abz1%init(dm_xh,
"abz1")
348 call this%abx2%init(dm_xh,
"abx2")
349 call this%aby2%init(dm_xh,
"aby2")
350 call this%abz2%init(dm_xh,
"abz2")
351 call this%advx%init(dm_xh,
"advx")
352 call this%advy%init(dm_xh,
"advy")
353 call this%advz%init(dm_xh,
"advz")
356 call this%du%init(this%dm_Xh,
'du')
357 call this%dv%init(this%dm_Xh,
'dv')
358 call this%dw%init(this%dm_Xh,
'dw')
359 call this%dp%init(this%dm_Xh,
'dp')
362 call this%setup_bcs(
user, params)
366 if (found)
call this%write_boundary_conditions()
368 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
369 this%pr_projection_activ_step)
371 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
372 this%vel_projection_activ_step)
378 if (params%valid_path(
'case.fluid.flow_rate_force'))
then
379 call this%vol_flow%init(this%dm_Xh, params)
383 call neko_log%section(
"Pressure solver")
386 'case.fluid.pressure_solver.max_iterations', &
388 call json_get(params,
'case.fluid.pressure_solver.type', solver_type)
389 call json_get(params,
'case.fluid.pressure_solver.preconditioner.type', &
392 'case.fluid.pressure_solver.preconditioner', precon_params)
393 call json_get(params,
'case.fluid.pressure_solver.absolute_tolerance', &
397 call neko_log%message(
'Type : ('// trim(solver_type) // &
398 ', ' // trim(precon_type) //
')')
399 write(log_buf,
'(A,ES13.6)')
'Abs tol :', abs_tol
402 call this%solver_factory(this%ksp_prs, this%dm_Xh%size(), &
403 solver_type, solver_maxiter, abs_tol, monitor)
404 call this%precon_factory_(this%pc_prs, this%ksp_prs, &
405 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, &
406 precon_type, precon_params)
411 call json_get(params,
'case.numerics', numerics_params)
412 call advection_factory(this%adv, numerics_params, this%c_Xh, &
413 this%ulag, this%vlag, this%wlag, &
414 chkp%dtlag, chkp%tlag, this%ext_bdf, &
420 call this%chkp%init(this%u, this%v, this%w, this%p)
422 this%chkp%abx1 => this%abx1
423 this%chkp%abx2 => this%abx2
424 this%chkp%aby1 => this%aby1
425 this%chkp%aby2 => this%aby2
426 this%chkp%abz1 => this%abz1
427 this%chkp%abz2 => this%abz2
428 call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
436 subroutine fluid_pnpn_restart(this, chkp)
437 class(fluid_pnpn_t),
target,
intent(inout) :: this
438 type(chkp_t),
intent(inout) :: chkp
439 real(kind=rp) :: dtlag(10), tlag(10)
440 type(field_t) :: u_temp, v_temp, w_temp
446 n = this%u%dof%size()
447 if (
allocated(chkp%previous_mesh%elements) .or. &
448 chkp%previous_Xh%lx .ne. this%Xh%lx)
then
449 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
450 c_xh => this%c_Xh, ulag => this%ulag, vlag => this%vlag, &
452 do concurrent(j = 1:n)
453 u%x(j,1,1,1) = u%x(j,1,1,1) * c_xh%mult(j,1,1,1)
454 v%x(j,1,1,1) = v%x(j,1,1,1) * c_xh%mult(j,1,1,1)
455 w%x(j,1,1,1) = w%x(j,1,1,1) * c_xh%mult(j,1,1,1)
456 p%x(j,1,1,1) = p%x(j,1,1,1) * c_xh%mult(j,1,1,1)
458 do i = 1, this%ulag%size()
459 do concurrent(j = 1:n)
460 ulag%lf(i)%x(j,1,1,1) = ulag%lf(i)%x(j,1,1,1) &
462 vlag%lf(i)%x(j,1,1,1) = vlag%lf(i)%x(j,1,1,1) &
464 wlag%lf(i)%x(j,1,1,1) = wlag%lf(i)%x(j,1,1,1) &
471 if (neko_bcknd_device .eq. 1)
then
472 associate(u => this%u, v => this%v, w => this%w, &
473 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
475 call device_memcpy(u%x, u%x_d, u%dof%size(), &
476 host_to_device, sync = .false.)
477 call device_memcpy(v%x, v%x_d, v%dof%size(), &
478 host_to_device, sync = .false.)
479 call device_memcpy(w%x, w%x_d, w%dof%size(), &
480 host_to_device, sync = .false.)
481 call device_memcpy(p%x, p%x_d, p%dof%size(), &
482 host_to_device, sync = .false.)
483 call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
484 u%dof%size(), host_to_device, sync = .false.)
485 call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
486 u%dof%size(), host_to_device, sync = .false.)
488 call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
489 v%dof%size(), host_to_device, sync = .false.)
490 call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
491 v%dof%size(), host_to_device, sync = .false.)
493 call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
494 w%dof%size(), host_to_device, sync = .false.)
495 call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
496 w%dof%size(), host_to_device, sync = .false.)
497 call device_memcpy(this%abx1%x, this%abx1%x_d, &
498 w%dof%size(), host_to_device, sync = .false.)
499 call device_memcpy(this%abx2%x, this%abx2%x_d, &
500 w%dof%size(), host_to_device, sync = .false.)
501 call device_memcpy(this%aby1%x, this%aby1%x_d, &
502 w%dof%size(), host_to_device, sync = .false.)
503 call device_memcpy(this%aby2%x, this%aby2%x_d, &
504 w%dof%size(), host_to_device, sync = .false.)
505 call device_memcpy(this%abz1%x, this%abz1%x_d, &
506 w%dof%size(), host_to_device, sync = .false.)
507 call device_memcpy(this%abz2%x, this%abz2%x_d, &
508 w%dof%size(), host_to_device, sync = .false.)
509 call device_memcpy(this%advx%x, this%advx%x_d, &
510 w%dof%size(), host_to_device, sync = .false.)
511 call device_memcpy(this%advy%x, this%advy%x_d, &
512 w%dof%size(), host_to_device, sync = .false.)
513 call device_memcpy(this%advz%x, this%advz%x_d, &
514 w%dof%size(), host_to_device, sync = .false.)
521 if (
allocated(chkp%previous_mesh%elements) &
522 .or. chkp%previous_Xh%lx .ne. this%Xh%lx)
then
524 call rotate_cyc(this%u%x, this%v%x, this%w%x, 1, this%c_Xh)
525 call this%gs_Xh%op(this%u, gs_op_add)
526 call this%gs_Xh%op(this%v, gs_op_add)
527 call this%gs_Xh%op(this%w, gs_op_add)
528 call this%gs_Xh%op(this%p, gs_op_add)
529 call rotate_cyc(this%u%x, this%v%x, this%w%x, 0, this%c_Xh)
531 do i = 1, this%ulag%size()
532 call rotate_cyc(this%ulag%lf(i)%x, this%vlag%lf(i)%x, this%wlag%lf(i)%x, 1, this%c_Xh)
533 call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
534 call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
535 call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
536 call rotate_cyc(this%ulag%lf(i)%x, this%vlag%lf(i)%x, this%wlag%lf(i)%x, 0, this%c_Xh)
639 subroutine fluid_pnpn_step(this, time, dt_controller)
640 class(fluid_pnpn_t),
target,
intent(inout) :: this
641 type(time_state_t),
intent(in) :: time
642 type(time_step_controller_t),
intent(in) :: dt_controller
646 type(ksp_monitor_t) :: ksp_results(4)
648 type(file_t) :: dump_file
649 class(bc_t),
pointer :: bc_i
650 type(non_normal_t),
pointer :: bc_j
652 if (this%freeze)
return
654 n = this%dm_Xh%size()
656 call profiler_start_region(
'Fluid', 1)
657 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
658 u_e => this%u_e, v_e => this%v_e, w_e => this%w_e, &
659 du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
660 u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
661 p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
663 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
664 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
665 msh => this%msh, prs_res => this%prs_res, &
666 source_term => this%source_term, vel_res => this%vel_res, &
667 sumab => this%sumab, makeoifs => this%makeoifs, &
668 makeabf => this%makeabf, makebdf => this%makebdf, &
669 vel_projection_dim => this%vel_projection_dim, &
670 pr_projection_dim => this%pr_projection_dim, &
672 rho => this%rho, mu_tot => this%mu_tot, &
673 f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
674 t => time%t, tstep => time%tstep, dt => time%dt, &
675 ext_bdf => this%ext_bdf, event => glb_cmd_event)
678 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
679 ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
682 call this%source_term%compute(time)
685 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
686 this%dm_Xh%size(), time, strong = .false.)
690 call this%adv%compute(u, v, w, &
691 this%advx, this%advy, this%advz, &
692 xh, this%c_Xh, dm_xh%size(), dt)
699 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
700 this%abx2, this%aby2, this%abz2, &
701 f_x%x, f_y%x, f_z%x, &
702 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
705 call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
706 f_x%x, f_y%x, f_z%x, &
707 rho%x(1,1,1,1), dt, n)
710 call this%adv%compute(u, v, w, &
712 xh, this%c_Xh, dm_xh%size())
718 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
719 this%abx2, this%aby2, this%abz2, &
720 f_x%x, f_y%x, f_z%x, &
721 rho%x(1,1,1,1), ext_bdf%advection_coeffs, n)
724 call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
725 u, v, w, c_xh%B, rho%x(1,1,1,1), dt, &
726 ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
733 call this%bc_apply_vel(time, strong = .true.)
734 call this%bc_apply_prs(time)
737 call this%update_material_properties(time)
740 call profiler_start_region(
'Pressure_residual', 18)
742 call prs_res%compute(p, p_res,&
747 this%bc_prs_surface, this%bc_sym_surface,&
748 ax_prs, ext_bdf%diffusion_coeffs(1), dt, &
752 if (.not. this%prs_dirichlet)
call ortho(p_res%x, this%glb_n_points, n)
754 call gs_xh%op(p_res, gs_op_add, event)
755 call device_event_sync(event)
758 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
761 call profiler_end_region(
'Pressure_residual', 18)
764 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
767 call this%pc_prs%update()
769 call profiler_start_region(
'Pressure_solve', 3)
773 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
774 this%bclst_dp, gs_xh)
775 ksp_results(1)%name =
'Pressure'
778 call profiler_end_region(
'Pressure_solve', 3)
780 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
781 this%bclst_dp, gs_xh, n, tstep, dt_controller)
784 call field_add2(p, dp, n)
785 if (.not. this%prs_dirichlet)
call ortho(p%x, this%glb_n_points, n)
788 call profiler_start_region(
'Velocity_residual', 19)
789 call vel_res%compute(ax_vel, u, v, w, &
790 u_res, v_res, w_res, &
794 mu_tot, rho, ext_bdf%diffusion_coeffs(1), &
797 call rotate_cyc(u_res%x, v_res%x, w_res%x, 1, c_xh)
798 call gs_xh%op(u_res, gs_op_add, event)
799 call device_event_sync(event)
800 call gs_xh%op(v_res, gs_op_add, event)
801 call device_event_sync(event)
802 call gs_xh%op(w_res, gs_op_add, event)
803 call device_event_sync(event)
804 call rotate_cyc(u_res%x, v_res%x, w_res%x, 0, c_xh)
807 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
810 call profiler_end_region(
'Velocity_residual', 19)
812 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
813 tstep, c_xh, n, dt_controller,
'Velocity')
815 call this%pc_vel%update()
817 call profiler_start_region(
"Velocity_solve", 4)
818 ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
819 u_res%x, v_res%x, w_res%x, n, c_xh, &
820 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
821 this%ksp_vel%max_iter)
822 call profiler_end_region(
"Velocity_solve", 4)
823 if (this%full_stress_formulation)
then
824 ksp_results(2)%name =
'Momentum'
826 ksp_results(2)%name =
'X-Velocity'
827 ksp_results(3)%name =
'Y-Velocity'
828 ksp_results(4)%name =
'Z-Velocity'
831 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
832 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
835 if (neko_bcknd_device .eq. 1)
then
836 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
837 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
839 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
842 if (this%forced_flow_rate)
then
844 call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
845 c_xh, gs_xh, ext_bdf, rho%x(1,1,1,1), mu_tot, &
846 dt, this%bclst_dp, this%bclst_du, this%bclst_dv, &
847 this%bclst_dw, this%bclst_vel_res, ax_vel, ax_prs, this%ksp_prs, &
848 this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
849 this%ksp_vel%max_iter)
852 call fluid_step_info(time, ksp_results, &
853 this%full_stress_formulation, this%strict_convergence, &
854 this%allow_stabilization)
857 call profiler_end_region(
'Fluid', 1)
862 subroutine fluid_pnpn_setup_bcs(this, user, params)
863 class(fluid_pnpn_t),
target,
intent(inout) :: this
864 type(user_t),
target,
intent(in) :: user
865 type(json_file),
intent(inout) :: params
866 integer :: i, n_bcs, zone_index, j, zone_size, global_zone_size, ierr
867 class(bc_t),
pointer :: bc_i
868 type(json_core) :: core
869 type(json_value),
pointer :: bc_object
870 type(json_file) :: bc_subdict
873 logical,
allocatable :: marked_zones(:)
874 integer,
allocatable :: zone_indices(:)
877 call this%bclst_vel_res%init()
878 call this%bclst_du%init()
879 call this%bclst_dv%init()
880 call this%bclst_dw%init()
881 call this%bclst_dp%init()
883 call this%bc_vel_res%init_from_components(this%c_Xh)
884 call this%bc_du%init_from_components(this%c_Xh)
885 call this%bc_dv%init_from_components(this%c_Xh)
886 call this%bc_dw%init_from_components(this%c_Xh)
887 call this%bc_dp%init_from_components(this%c_Xh)
890 call this%bc_prs_surface%init_from_components(this%c_Xh)
891 call this%bc_sym_surface%init_from_components(this%c_Xh)
894 if (params%valid_path(
'case.fluid.boundary_conditions'))
then
895 call params%info(
'case.fluid.boundary_conditions', n_children = n_bcs)
896 call params%get_core(core)
897 call params%get(
'case.fluid.boundary_conditions', bc_object, found)
902 call this%bcs_vel%init(n_bcs)
904 allocate(marked_zones(
size(this%msh%labeled_zones)))
905 marked_zones = .false.
909 call json_extract_item(core, bc_object, i, bc_subdict)
911 call json_get(bc_subdict,
"zone_indices", zone_indices)
916 do j = 1,
size(zone_indices)
917 zone_size = this%msh%labeled_zones(zone_indices(j))%size
918 call mpi_allreduce(zone_size, global_zone_size, 1, &
919 mpi_integer, mpi_max, neko_comm, ierr)
921 if (global_zone_size .eq. 0)
then
922 write(error_unit,
'(A, A, I0, A, A, I0, A)')
"*** ERROR ***: ",&
923 "Zone index ", zone_indices(j), &
924 " is invalid as this zone has 0 size, meaning it ", &
925 "is not in the mesh. Check fluid boundary condition ", &
930 if (marked_zones(zone_indices(j)) .eqv. .true.)
then
931 write(error_unit,
'(A, A, I0, A, A, A, A)')
"*** ERROR ***: ", &
932 "Zone with index ", zone_indices(j), &
933 " has already been assigned a boundary condition. ", &
934 "Please check your boundary_conditions entry for the ", &
935 "fluid and make sure that each zone index appears only ", &
936 "in a single boundary condition."
939 marked_zones(zone_indices(j)) = .true.
944 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
948 if (
associated(bc_i))
then
962 call this%bclst_vel_res%append(bc_i)
963 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
964 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
965 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
967 call this%bcs_vel%append(bc_i)
969 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
970 type is (non_normal_t)
973 call this%bclst_vel_res%append(bc_i)
974 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
975 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
976 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
977 type is (shear_stress_t)
979 call this%bclst_vel_res%append(bc_i%symmetry)
980 call this%bclst_du%append(bc_i%symmetry%bc_x)
981 call this%bclst_dv%append(bc_i%symmetry%bc_y)
982 call this%bclst_dw%append(bc_i%symmetry%bc_z)
984 call this%bcs_vel%append(bc_i)
985 type is (wall_model_bc_t)
987 call this%bclst_vel_res%append(bc_i%symmetry)
988 call this%bclst_du%append(bc_i%symmetry%bc_x)
989 call this%bclst_dv%append(bc_i%symmetry%bc_y)
990 call this%bclst_dw%append(bc_i%symmetry%bc_z)
992 call this%bcs_vel%append(bc_i)
999 if (bc_i%strong .eqv. .true.)
then
1000 call this%bc_vel_res%mark_facets(bc_i%marked_facet)
1001 call this%bc_du%mark_facets(bc_i%marked_facet)
1002 call this%bc_dv%mark_facets(bc_i%marked_facet)
1003 call this%bc_dw%mark_facets(bc_i%marked_facet)
1005 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
1008 call this%bcs_vel%append(bc_i)
1014 do i = 1,
size(this%msh%labeled_zones)
1015 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1016 (marked_zones(i) .eqv. .false.))
then
1017 write(error_unit,
'(A, A, I0)')
"*** ERROR ***: ", &
1018 "No fluid boundary condition assigned to zone ", i
1026 call this%bcs_prs%init(n_bcs)
1030 call json_extract_item(core, bc_object, i, bc_subdict)
1032 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh,
user)
1036 if (
associated(bc_i))
then
1037 call this%bcs_prs%append(bc_i)
1040 if (bc_i%strong .eqv. .true.)
then
1041 call this%bc_dp%mark_facets(bc_i%marked_facet)
1049 do i = 1,
size(this%msh%labeled_zones)
1050 if (this%msh%labeled_zones(i)%size .gt. 0)
then
1051 call neko_error(
"No boundary_conditions entry in the case file!")
1057 call this%bcs_vel%init()
1058 call this%bcs_prs%init()
1062 call this%bc_prs_surface%finalize()
1063 call this%bc_sym_surface%finalize()
1065 call this%bc_vel_res%finalize()
1066 call this%bc_du%finalize()
1067 call this%bc_dv%finalize()
1068 call this%bc_dw%finalize()
1069 call this%bc_dp%finalize()
1071 call this%bclst_vel_res%append(this%bc_vel_res)
1072 call this%bclst_du%append(this%bc_du)
1073 call this%bclst_dv%append(this%bc_dv)
1074 call this%bclst_dw%append(this%bc_dw)
1075 call this%bclst_dp%append(this%bc_dp)
1078 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1079 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1080 mpi_logical, mpi_lor, neko_comm)
1083 if (
allocated(marked_zones))
then
1084 deallocate(marked_zones)
1087 if (
allocated(zone_indices))
then
1088 deallocate(zone_indices)
1094 subroutine fluid_pnpn_write_boundary_conditions(this)
1100 class(fluid_pnpn_t),
target,
intent(inout) :: this
1101 type(dirichlet_t) :: bdry_mask
1102 type(field_t),
pointer :: bdry_field
1103 type(file_t) :: bdry_file
1104 integer :: temp_index, i
1105 class(bc_t),
pointer :: bci
1106 character(len=LOG_SIZE) :: log_buf
1108 call neko_log%section(
"Fluid boundary conditions")
1109 write(log_buf,
'(A)')
'Marking using integer keys in bdry0.f00000'
1110 call neko_log%message(log_buf)
1111 write(log_buf,
'(A)')
'Condition-value pairs: '
1112 call neko_log%message(log_buf)
1113 write(log_buf,
'(A)')
' no_slip = 1'
1114 call neko_log%message(log_buf)
1115 write(log_buf,
'(A)')
' velocity_value = 2'
1116 call neko_log%message(log_buf)
1117 write(log_buf,
'(A)')
' outflow, normal_outflow (+dong) = 3'
1118 call neko_log%message(log_buf)
1119 write(log_buf,
'(A)')
' symmetry = 4'
1120 call neko_log%message(log_buf)
1121 write(log_buf,
'(A)')
' periodic = 6'
1122 call neko_log%message(log_buf)
1123 write(log_buf,
'(A)')
' user_velocity = 7'
1124 call neko_log%message(log_buf)
1125 write(log_buf,
'(A)')
' user_pressure = 8'
1126 call neko_log%message(log_buf)
1127 write(log_buf,
'(A)')
' shear_stress = 9'
1128 call neko_log%message(log_buf)
1129 write(log_buf,
'(A)')
' wall_modelling = 10'
1130 call neko_log%message(log_buf)
1131 write(log_buf,
'(A)')
' blasius_profile = 11'
1132 call neko_log%message(log_buf)
1133 call neko_log%end_section()
1135 call neko_scratch_registry%request_field(bdry_field, temp_index, .true.)
1139 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1140 call bdry_mask%mark_zone(this%msh%periodic)
1141 call bdry_mask%finalize()
1142 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1143 call bdry_mask%free()
1145 do i = 1, this%bcs_prs%size()
1146 bci => this%bcs_prs%get(i)
1147 select type (
bc => bci)
1148 type is (zero_dirichlet_t)
1149 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1150 call bdry_mask%mark_facets(bci%marked_facet)
1151 call bdry_mask%finalize()
1152 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1153 call bdry_mask%free()
1155 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1156 call bdry_mask%mark_facets(bci%marked_facet)
1157 call bdry_mask%finalize()
1158 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1159 call bdry_mask%free()
1161 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1162 call bdry_mask%mark_facets(bci%marked_facet)
1163 call bdry_mask%finalize()
1164 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1165 call bdry_mask%free()
1169 do i = 1, this%bcs_vel%size()
1170 bci => this%bcs_vel%get(i)
1171 select type (
bc => bci)
1172 type is (zero_dirichlet_t)
1173 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1174 call bdry_mask%mark_facets(bci%marked_facet)
1175 call bdry_mask%finalize()
1176 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1177 call bdry_mask%free()
1179 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1180 call bdry_mask%mark_facets(bci%marked_facet)
1181 call bdry_mask%finalize()
1182 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1183 call bdry_mask%free()
1184 type is (symmetry_t)
1185 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1186 call bdry_mask%mark_facets(bci%marked_facet)
1187 call bdry_mask%finalize()
1188 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1189 call bdry_mask%free()
1191 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1192 call bdry_mask%mark_facets(bci%marked_facet)
1193 call bdry_mask%finalize()
1194 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1195 call bdry_mask%free()
1196 type is (shear_stress_t)
1197 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1198 call bdry_mask%mark_facets(bci%marked_facet)
1199 call bdry_mask%finalize()
1200 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1201 call bdry_mask%free()
1202 type is (wall_model_bc_t)
1203 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1204 call bdry_mask%mark_facets(bci%marked_facet)
1205 call bdry_mask%finalize()
1206 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1207 call bdry_mask%free()
1209 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1210 call bdry_mask%mark_facets(bci%marked_facet)
1211 call bdry_mask%finalize()
1212 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1213 call bdry_mask%free()
1218 call bdry_file%init(
'bdry.fld')
1219 call bdry_file%write(bdry_field)
1221 call neko_scratch_registry%relinquish_field(temp_index)