Skip to content

Commit d7c0801

Browse files
committed
correctness fixes for projection
1 parent 96d5f2f commit d7c0801

11 files changed

Lines changed: 466 additions & 201 deletions

src/common/m_boundary_common.fpp

Lines changed: 114 additions & 18 deletions
Original file line numberDiff line numberDiff line change
@@ -1212,6 +1212,8 @@ contains
12121212
select case (bc_x%beg)
12131213
case (BC_PERIODIC)
12141214
call s_beta_periodic(q_beta, 1, -1, k, l, nvar)
1215+
case (BC_REFLECTIVE)
1216+
call s_beta_reflective(q_beta, 1, -1, k, l, nvar)
12151217
case default
12161218
end select
12171219
end do
@@ -1228,6 +1230,8 @@ contains
12281230
select case (bc_x%end)
12291231
case (BC_PERIODIC)
12301232
call s_beta_periodic(q_beta, 1, 1, k, l, nvar)
1233+
case (BC_REFLECTIVE)
1234+
call s_beta_reflective(q_beta, 1, 1, k, l, nvar)
12311235
case default
12321236
end select
12331237
end do
@@ -1245,13 +1249,35 @@ contains
12451249
select case (bc_y%beg)
12461250
case (BC_PERIODIC)
12471251
call s_beta_periodic(q_beta, 2, -1, k, l, nvar)
1252+
case (BC_REFLECTIVE)
1253+
call s_beta_reflective(q_beta, 2, -1, k, l, nvar)
12481254
case default
12491255
end select
12501256
end do
12511257
end do
12521258
$:END_GPU_PARALLEL_LOOP()
12531259
end if
12541260

1261+
if (bc_y%end >= 0) then
1262+
call s_mpi_reduce_beta_variables_buffers(q_beta, 2, 1, nvar)
1263+
else
1264+
$:GPU_PARALLEL_LOOP(private='[l,k]', collapse=2)
1265+
do l = beta_bc_bounds(3)%beg, beta_bc_bounds(3)%end
1266+
do k = beta_bc_bounds(1)%beg, beta_bc_bounds(1)%end
1267+
select case (bc_y%end)
1268+
case (BC_PERIODIC)
1269+
call s_beta_periodic(q_beta, 2, 1, k, l, nvar)
1270+
case (BC_REFLECTIVE)
1271+
call s_beta_reflective(q_beta, 2, 1, k, l, nvar)
1272+
case default
1273+
end select
1274+
end do
1275+
end do
1276+
$:END_GPU_PARALLEL_LOOP()
1277+
end if
1278+
1279+
if (num_dims == 2) return
1280+
12551281
#:if not MFC_CASE_OPTIMIZATION or num_dims > 2
12561282
!< z-direction
12571283
if (bc_z%beg >= 0) then
@@ -1263,31 +1289,15 @@ contains
12631289
select case (bc_type(3, 1)%sf(k, l, 0))
12641290
case (BC_PERIODIC)
12651291
call s_beta_periodic(q_beta, 3, -1, k, l, nvar)
1292+
case (BC_REFLECTIVE)
1293+
call s_beta_reflective(q_beta, 3, -1, k, l, nvar)
12661294
case default
12671295
end select
12681296
end do
12691297
end do
12701298
$:END_GPU_PARALLEL_LOOP()
12711299
end if
1272-
#:endif
12731300

1274-
if (bc_y%end >= 0) then
1275-
call s_mpi_reduce_beta_variables_buffers(q_beta, 2, 1, nvar)
1276-
else
1277-
$:GPU_PARALLEL_LOOP(private='[l,k]', collapse=2)
1278-
do l = beta_bc_bounds(3)%beg, beta_bc_bounds(3)%end
1279-
do k = beta_bc_bounds(1)%beg, beta_bc_bounds(1)%end
1280-
select case (bc_y%end)
1281-
case (BC_PERIODIC)
1282-
call s_beta_periodic(q_beta, 2, 1, k, l, nvar)
1283-
case default
1284-
end select
1285-
end do
1286-
end do
1287-
$:END_GPU_PARALLEL_LOOP()
1288-
end if
1289-
1290-
#:if not MFC_CASE_OPTIMIZATION or num_dims > 2
12911301
if (bc_z%end >= 0) then
12921302
call s_mpi_reduce_beta_variables_buffers(q_beta, 3, 1, nvar)
12931303
else
@@ -1297,6 +1307,8 @@ contains
12971307
select case (bc_type(3, 2)%sf(k, l, 0))
12981308
case (BC_PERIODIC)
12991309
call s_beta_periodic(q_beta, 3, 1, k, l, nvar)
1310+
case (BC_REFLECTIVE)
1311+
call s_beta_reflective(q_beta, 3, 1, k, l, nvar)
13001312
case default
13011313
end select
13021314
end do
@@ -1424,6 +1436,90 @@ contains
14241436

14251437
end subroutine s_beta_extrapolation
14261438

1439+
subroutine s_beta_reflective(q_beta, bc_dir, bc_loc, k, l, nvar)
1440+
$:GPU_ROUTINE(function_name='s_beta_reflective', &
1441+
& parallelism='[seq]', cray_inline=True)
1442+
type(scalar_field), dimension(num_dims + 1), intent(inout) :: q_beta
1443+
integer, intent(in) :: bc_dir, bc_loc
1444+
integer, intent(in) :: k, l
1445+
integer, intent(in) :: nvar
1446+
1447+
integer :: j, i
1448+
1449+
! Reflective BC for void fraction:
1450+
! 1) Fold ghost-cell contributions back onto their mirror interior cells
1451+
! 2) Set ghost cells = mirror of (now-folded) interior values
1452+
1453+
if (bc_dir == 1) then !< x-direction
1454+
if (bc_loc == -1) then !< bc_x%beg
1455+
do i = 1, nvar
1456+
do j = 1, mapCells + 1
1457+
q_beta(beta_vars(i))%sf(j - 1, k, l) = q_beta(beta_vars(i))%sf(j - 1, k, l) + &
1458+
q_beta(beta_vars(i))%sf(-j, k, l)
1459+
end do
1460+
do j = 1, mapCells + 1
1461+
q_beta(beta_vars(i))%sf(-j, k, l) = q_beta(beta_vars(i))%sf(j - 1, k, l)
1462+
end do
1463+
end do
1464+
else !< bc_x%end
1465+
do i = 1, nvar
1466+
do j = 1, mapCells + 1
1467+
q_beta(beta_vars(i))%sf(m - (j - 1), k, l) = q_beta(beta_vars(i))%sf(m - (j - 1), k, l) + &
1468+
q_beta(beta_vars(i))%sf(m + j, k, l)
1469+
end do
1470+
do j = 1, mapCells + 1
1471+
q_beta(beta_vars(i))%sf(m + j, k, l) = q_beta(beta_vars(i))%sf(m - (j - 1), k, l)
1472+
end do
1473+
end do
1474+
end if
1475+
elseif (bc_dir == 2) then !< y-direction
1476+
if (bc_loc == -1) then !< bc_y%beg
1477+
do i = 1, nvar
1478+
do j = 1, mapCells + 1
1479+
q_beta(beta_vars(i))%sf(k, j - 1, l) = q_beta(beta_vars(i))%sf(k, j - 1, l) + &
1480+
q_beta(beta_vars(i))%sf(k, -j, l)
1481+
end do
1482+
do j = 1, mapCells + 1
1483+
q_beta(beta_vars(i))%sf(k, -j, l) = q_beta(beta_vars(i))%sf(k, j - 1, l)
1484+
end do
1485+
end do
1486+
else !< bc_y%end
1487+
do i = 1, nvar
1488+
do j = 1, mapCells + 1
1489+
q_beta(beta_vars(i))%sf(k, n - (j - 1), l) = q_beta(beta_vars(i))%sf(k, n - (j - 1), l) + &
1490+
q_beta(beta_vars(i))%sf(k, n + j, l)
1491+
end do
1492+
do j = 1, mapCells + 1
1493+
q_beta(beta_vars(i))%sf(k, n + j, l) = q_beta(beta_vars(i))%sf(k, n - (j - 1), l)
1494+
end do
1495+
end do
1496+
end if
1497+
elseif (bc_dir == 3) then !< z-direction
1498+
if (bc_loc == -1) then !< bc_z%beg
1499+
do i = 1, nvar
1500+
do j = 1, mapCells + 1
1501+
q_beta(beta_vars(i))%sf(k, l, j - 1) = q_beta(beta_vars(i))%sf(k, l, j - 1) + &
1502+
q_beta(beta_vars(i))%sf(k, l, -j)
1503+
end do
1504+
do j = 1, mapCells + 1
1505+
q_beta(beta_vars(i))%sf(k, l, -j) = q_beta(beta_vars(i))%sf(k, l, j - 1)
1506+
end do
1507+
end do
1508+
else !< bc_z%end
1509+
do i = 1, nvar
1510+
do j = 1, mapCells + 1
1511+
q_beta(beta_vars(i))%sf(k, l, p - (j - 1)) = q_beta(beta_vars(i))%sf(k, l, p - (j - 1)) + &
1512+
q_beta(beta_vars(i))%sf(k, l, p + j)
1513+
end do
1514+
do j = 1, mapCells + 1
1515+
q_beta(beta_vars(i))%sf(k, l, p + j) = q_beta(beta_vars(i))%sf(k, l, p - (j - 1))
1516+
end do
1517+
end do
1518+
end if
1519+
end if
1520+
1521+
end subroutine s_beta_reflective
1522+
14271523
impure subroutine s_populate_capillary_buffers(c_divs, bc_type)
14281524

14291525
type(scalar_field), dimension(num_dims + 1), intent(inout) :: c_divs

src/common/m_derived_types.fpp

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -492,6 +492,7 @@ module m_derived_types
492492
character(LEN=pathlen_max) :: input_path !< Path to lag_bubbles.dat
493493
real(wp) :: epsilonb !< Standard deviation scaling for the gaussian function
494494
real(wp) :: charwidth !< Domain virtual depth (z direction, for 2D simulations)
495+
integer :: charNz !< Number of grid cells in characteristic depth
495496
real(wp) :: valmaxvoid !< Maximum void fraction permitted
496497
497498
end type bubbles_lagrange_parameters

src/common/m_mpi_common.fpp

Lines changed: 10 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -1169,12 +1169,16 @@ contains
11691169
call nvtxStartRange("BETA-COMM-PACKBUF")
11701170
11711171
! Set bounds for each dimension
1172-
comm_coords(1)%beg = merge(-mapcells - 1, 0, bc_x%beg >= 0)
1173-
comm_coords(1)%end = merge(m + mapcells + 1, m, bc_x%end >= 0)
1174-
comm_coords(2)%beg = merge(-mapcells - 1, 0, bc_y%beg >= 0)
1175-
comm_coords(2)%end = merge(n + mapcells + 1, n, bc_y%end >= 0)
1176-
comm_coords(3)%beg = merge(-mapcells - 1, 0, (bc_z%beg >= 0 .and. p > 0))
1177-
comm_coords(3)%end = merge(p + mapcells + 1, p, (bc_z%end >= 0 .and. p > 0))
1172+
! Always include the full buffer range for each existing dimension.
1173+
! The Gaussian smearing kernel writes to buffer cells even at physical
1174+
! boundaries, and these contributions must be communicated to neighbors
1175+
! in other directions via ADD operations.
1176+
comm_coords(1)%beg = -mapcells - 1
1177+
comm_coords(1)%end = m + mapcells + 1
1178+
comm_coords(2)%beg = merge(-mapcells - 1, 0, n > 0)
1179+
comm_coords(2)%end = merge(n + mapcells + 1, n, n > 0)
1180+
comm_coords(3)%beg = merge(-mapcells - 1, 0, p > 0)
1181+
comm_coords(3)%end = merge(p + mapcells + 1, p, p > 0)
11781182
11791183
! Compute sizes
11801184
comm_size(1) = comm_coords(1)%end - comm_coords(1)%beg + 1

src/pre_process/m_global_parameters.fpp

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -461,6 +461,7 @@ contains
461461
lag_params%drag_model = dflt_int
462462
lag_params%epsilonb = 1._wp
463463
lag_params%charwidth = dflt_real
464+
lag_params%charNz = dflt_int
464465
lag_params%valmaxvoid = dflt_real
465466
466467
do i = 1, num_patches_max

0 commit comments

Comments
 (0)