Skip to content

Commit

Permalink
boundary conditions refactor and MPI fix for elliptic smoothing
Browse files Browse the repository at this point in the history
  • Loading branch information
wilfonba committed Mar 8, 2025
1 parent a0ed080 commit 2e42cc7
Show file tree
Hide file tree
Showing 12 changed files with 89 additions and 303 deletions.
Original file line number Diff line number Diff line change
Expand Up @@ -16,9 +16,11 @@ module m_boundary_conditions

implicit none

private;
public :: s_populate_variables_buffers, &
s_populate_capillary_buffers
#ifdef MFC_SIMULATION
private; public :: s_populate_variables_buffers, s_populate_capillary_buffers
#else
private; public :: s_populate_variables_buffers
#endif

contains

Expand All @@ -28,7 +30,8 @@ contains
subroutine s_populate_variables_buffers(q_prim_vf, pb, mv)

type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv

real(wp), optional, dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv

integer :: bc_loc, bc_dir

Expand Down Expand Up @@ -66,6 +69,7 @@ contains
q_prim_vf, pb, mv, 1, 1)
end select

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
select case (bc_x%beg)
case (-13:-3) ! Ghost-cell extrap. BC at beginning
Expand All @@ -85,8 +89,7 @@ contains
call s_qbmm_extrapolation(pb, mv, 1, 1)
end select
end if

! END: Population of Buffers in x-direction
#endif

! Population of Buffers in y-direction

Expand Down Expand Up @@ -126,8 +129,8 @@ contains
q_prim_vf, pb, mv, 2, 1)
end select

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then

select case (bc_y%beg)
case (-13:-3) ! Ghost-cell extrap. BC at beginning
call s_qbmm_extrapolation(pb, mv, 2, -1)
Expand All @@ -145,10 +148,8 @@ contains
case (-16) ! No-slip wall BC at end
call s_qbmm_extrapolation(pb, mv, 2, 1)
end select

end if

! END: Population of Buffers in y-direction
#endif

! Population of Buffers in z-direction

Expand Down Expand Up @@ -186,8 +187,8 @@ contains
q_prim_vf, pb, mv, 3, 1)
end select

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then

select case (bc_z%beg)
case (-13:-3) ! Ghost-cell extrap. BC at beginning
call s_qbmm_extrapolation(pb, mv, 3, -1)
Expand All @@ -205,17 +206,16 @@ contains
case (-16) ! No-slip wall BC at end
call s_qbmm_extrapolation(pb, mv, 3, 1)
end select

end if

#endif
! END: Population of Buffers in z-direction

end subroutine s_populate_variables_buffers

subroutine s_ghost_cell_extrapolation(q_prim_vf, pb, mv, bc_dir, bc_loc)

type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv
real(wp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv
integer, intent(in) :: bc_dir, bc_loc
integer :: j, k, l, q, i

Expand Down Expand Up @@ -325,7 +325,7 @@ contains
subroutine s_symmetry(q_prim_vf, pb, mv, bc_dir, bc_loc)

type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv
real(wp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv
integer, intent(in) :: bc_dir, bc_loc

integer :: j, k, l, q, i
Expand Down Expand Up @@ -362,7 +362,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -380,6 +380,7 @@ contains
end do
end do
end if
#endif

else !< bc_x%end

Expand Down Expand Up @@ -411,7 +412,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -429,7 +430,7 @@ contains
end do
end do
end if

#endif
end if

!< y-direction
Expand Down Expand Up @@ -463,7 +464,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -481,7 +482,7 @@ contains
end do
end do
end if

#endif
else !< bc_y%end

!$acc parallel loop collapse(3) gang vector default(present)
Expand Down Expand Up @@ -510,7 +511,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -528,7 +529,7 @@ contains
end do
end do
end if

#endif
end if

!< z-direction
Expand Down Expand Up @@ -562,7 +563,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -580,7 +581,7 @@ contains
end do
end do
end if

#endif
else !< bc_z%end

!$acc parallel loop collapse(3) gang vector default(present)
Expand Down Expand Up @@ -609,7 +610,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -627,7 +628,7 @@ contains
end do
end do
end if

#endif
end if

end if
Expand All @@ -637,7 +638,7 @@ contains
subroutine s_periodic(q_prim_vf, pb, mv, bc_dir, bc_loc)

type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv
real(wp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv
integer, intent(in) :: bc_dir, bc_loc

integer :: j, k, l, q, i
Expand All @@ -658,7 +659,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -676,7 +677,7 @@ contains
end do
end do
end if

#endif
else !< bc_x%end

!$acc parallel loop collapse(4) gang vector default(present)
Expand All @@ -690,7 +691,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -708,7 +709,7 @@ contains
end do
end do
end if

#endif
end if

!< y-direction
Expand All @@ -727,7 +728,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(4) gang vector default(present)
do i = 1, nb
Expand All @@ -745,7 +746,7 @@ contains
end do
end do
end if

#endif
else !< bc_y%end

!$acc parallel loop collapse(4) gang vector default(present)
Expand All @@ -759,7 +760,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -777,7 +778,7 @@ contains
end do
end do
end if

#endif
end if

!< z-direction
Expand All @@ -796,7 +797,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -814,7 +815,7 @@ contains
end do
end do
end if

#endif
else !< bc_z%end

!$acc parallel loop collapse(4) gang vector default(present)
Expand All @@ -828,7 +829,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -846,7 +847,7 @@ contains
end do
end do
end if

#endif
end if

end if
Expand All @@ -856,7 +857,7 @@ contains
subroutine s_axis(q_prim_vf, pb, mv, bc_dir, bc_loc)

type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv
real(wp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv
integer, intent(in) :: bc_dir, bc_loc

integer :: j, k, l, q, i
Expand Down Expand Up @@ -905,7 +906,7 @@ contains
end do
end do
end do

#ifdef MFC_SIMULATION
if (qbmm .and. .not. polytropic) then
!$acc parallel loop collapse(5) gang vector default(present)
do i = 1, nb
Expand All @@ -923,13 +924,13 @@ contains
end do
end do
end if

#endif
end subroutine s_axis

subroutine s_slip_wall(q_prim_vf, pb, mv, bc_dir, bc_loc)

type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv
real(wp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv
integer, intent(in) :: bc_dir, bc_loc

integer :: j, k, l, q, i
Expand Down Expand Up @@ -1070,7 +1071,7 @@ contains
subroutine s_no_slip_wall(q_prim_vf, pb, mv, bc_dir, bc_loc)

type(scalar_field), dimension(sys_size), intent(inout) :: q_prim_vf
real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv
real(wp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv
integer, intent(in) :: bc_dir, bc_loc

integer :: j, k, l, q, i
Expand Down Expand Up @@ -1246,7 +1247,7 @@ contains

subroutine s_qbmm_extrapolation(pb, mv, bc_dir, bc_loc)

real(wp), dimension(startx:, starty:, startz:, 1:, 1:), intent(inout) :: pb, mv
real(wp), dimension(idwbuff(1)%beg:, idwbuff(2)%beg:, idwbuff(3)%beg:, 1:, 1:), intent(inout) :: pb, mv
integer, intent(in) :: bc_dir, bc_loc

integer :: j, k, l, q, i
Expand Down Expand Up @@ -1378,6 +1379,7 @@ contains

end subroutine s_qbmm_extrapolation

#ifdef MFC_SIMULATION
subroutine s_populate_capillary_buffers(c_divs)

type(scalar_field), dimension(num_dims + 1), intent(inout) :: c_divs
Expand Down Expand Up @@ -1659,5 +1661,6 @@ contains
end if

end subroutine s_populate_capillary_buffers
#endif

end module m_boundary_conditions
Loading

0 comments on commit 2e42cc7

Please sign in to comment.