Skip to content

Commit 7f2cd02

Browse files
committed
refactor: move dx/dy/dz min spacing into grid_axis%min_spacing
Remove module-level dx, dy, dz scalars from pre_process m_global_parameters. Add min_spacing field to the grid_axis derived type so each axis carries its own minimum cell width. Update all call sites in m_grid, m_start_up, m_icpp_patches, m_mpi_common, and 2dHardcodedIC.
1 parent 15f3a64 commit 7f2cd02

7 files changed

Lines changed: 65 additions & 53 deletions

File tree

src/common/include/2dHardcodedIC.fpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -157,7 +157,7 @@
157157
if (x%cc(i) <= 0.7_wp*lam) then
158158
d = x%cc(i) - lam*(0.4_wp - 0.1_wp*sin(2.0_wp*pi*(y%cc(j)/lam + 0.25_wp)))
159159
#ifdef MFC_PRE_PROCESS
160-
fsm = 0.5_wp*(1.0_wp + erf(d/(ei*sqrt(dx*dy))))
160+
fsm = 0.5_wp*(1.0_wp + erf(d/(ei*sqrt(x%min_spacing*y%min_spacing))))
161161
#else
162162
fsm = 0.5_wp*(1.0_wp + erf(d/(ei*sqrt(x%spacing(i)*y%spacing(j)))))
163163
#endif

src/common/m_derived_types.fpp

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -13,9 +13,10 @@ module m_derived_types
1313

1414
implicit none
1515

16-
!> Derived type for a single spatial grid axis: cell-boundary, cell-center, and per-cell spacing arrays
16+
!> Derived type for a single spatial grid axis: cell-boundary, cell-center, per-cell spacing arrays, and minimum spacing scalar
1717
type grid_axis
1818
real(wp), allocatable, dimension(:) :: cb, cc, spacing
19+
real(wp) :: min_spacing = 0._wp
1920
end type grid_axis
2021

2122
!> Derived type adding the field position (fp) as an attribute

src/common/m_mpi_common.fpp

Lines changed: 9 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -1022,12 +1022,15 @@ contains
10221022
integer :: num_procs_x, num_procs_y, num_procs_z !< Optimal number of processors in the x-, y- and z-directions
10231023
!> Non-optimal number of processors in the x-, y- and z-directions
10241024
real(wp) :: tmp_num_procs_x, tmp_num_procs_y, tmp_num_procs_z
1025-
real(wp) :: fct_min !< Processor factorization (fct) minimization parameter
1026-
integer :: MPI_COMM_CART !< Cartesian processor topology communicator
1027-
integer :: rem_cells !< Remaining cells after distribution among processors
1028-
integer :: recon_order !< WENO or MUSCL reconstruction order
1029-
integer :: i, j !< Generic loop iterators
1030-
integer :: ierr !< Generic flag used to identify and report MPI errors
1025+
real(wp) :: fct_min !< Processor factorization (fct) minimization parameter
1026+
#ifdef MFC_PRE_PROCESS
1027+
real(wp) :: dx, dy, dz !< uniform spacing temporaries for domain decomposition
1028+
#endif
1029+
integer :: MPI_COMM_CART !< Cartesian processor topology communicator
1030+
integer :: rem_cells !< Remaining cells after distribution among processors
1031+
integer :: recon_order !< WENO or MUSCL reconstruction order
1032+
integer :: i, j !< Generic loop iterators
1033+
integer :: ierr !< Generic flag used to identify and report MPI errors
10311034

10321035
if (recon_type == WENO_TYPE) then
10331036
recon_order = weno_order

src/pre_process/m_global_parameters.fpp

Lines changed: 1 addition & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -41,9 +41,8 @@ module m_global_parameters
4141
integer :: num_vels !< Number of velocity components (different from num_dims for mhd)
4242
logical :: cyl_coord
4343
integer :: grid_geometry !< Cylindrical coordinates (either axisymmetric or full 3D)
44-
!> Cell-boundary (cb) and cell-center (cc) arrays per direction
44+
!> Cell-boundary (cb), cell-center (cc), per-cell spacing, and minimum spacing per direction
4545
type(grid_axis) :: x, y, z
46-
real(wp) :: dx, dy, dz !< Minimum cell-widths in the x-, y- and z-coordinate directions
4746
type(bounds_info) :: x_domain, y_domain, z_domain !< Locations of the domain bounds in the x-, y- and z-coordinate directions
4847
logical :: stretch_x, stretch_y, stretch_z !< Grid stretching flags for the x-, y- and z-coordinate directions
4948
! Grid stretching: a_x/a_y/a_z = rate, x_stretch%beg/end = locations

src/pre_process/m_grid.f90

Lines changed: 8 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -34,8 +34,9 @@ end subroutine s_generate_abstract_grid
3434
impure subroutine s_generate_serial_grid
3535

3636
! Generic loop iterator
37-
integer :: i, j !< generic loop operators
38-
real(wp) :: length !< domain lengths
37+
integer :: i, j !< generic loop operators
38+
real(wp) :: length !< domain lengths
39+
real(wp) :: dx, dy, dz !< local uniform spacing temporaries
3940
! Uniform grid: dx = (x_end - x_beg) / (m + 1)
4041

4142
dx = (x_domain%end - x_domain%beg)/real(m + 1, wp)
@@ -68,6 +69,7 @@ impure subroutine s_generate_serial_grid
6869
print *, 'Stretched grid: min/max x grid: ', minval(x%cc(:)), maxval(x%cc(:))
6970
if (num_procs > 1) call s_mpi_reduce_min(dx)
7071
end if
72+
x%min_spacing = dx
7173

7274
! Grid Generation in the y-direction
7375
if (n == 0) return
@@ -115,6 +117,7 @@ impure subroutine s_generate_serial_grid
115117

116118
if (num_procs > 1) call s_mpi_reduce_min(dy)
117119
end if
120+
y%min_spacing = dy
118121

119122
! Grid Generation in the z-direction
120123
if (p == 0) return
@@ -149,14 +152,16 @@ impure subroutine s_generate_serial_grid
149152

150153
if (num_procs > 1) call s_mpi_reduce_min(dz)
151154
end if
155+
z%min_spacing = dz
152156

153157
end subroutine s_generate_serial_grid
154158

155159
!> Generate a uniform or stretched rectilinear grid in parallel from user parameters.
156160
impure subroutine s_generate_parallel_grid
157161

158162
#ifdef MFC_MPI
159-
real(wp) :: length !< domain lengths
163+
real(wp) :: length !< domain lengths
164+
real(wp) :: dx, dy, dz !< local uniform spacing temporaries
160165
! Locations of cell boundaries
161166
real(wp), allocatable, dimension(:) :: x_cb_glb, y_cb_glb, z_cb_glb !< Locations of cell boundaries
162167
character(LEN=path_len + name_len) :: file_loc !< Generic string used to store the address of a file

src/pre_process/m_icpp_patches.fpp

Lines changed: 32 additions & 28 deletions
Original file line numberDiff line numberDiff line change
@@ -312,8 +312,9 @@ contains
312312
do i = 0, m
313313
if (patch_icpp(patch_id)%smoothen) then
314314
! Smooth Heaviside via hyperbolic tangent; smooth_coeff controls interface sharpness
315-
eta = tanh(smooth_coeff/min(dx, &
316-
& dy)*(sqrt((x%cc(i) - x_centroid)**2 + (y%cc(j) - y_centroid)**2) - radius))*(-0.5_wp) + 0.5_wp
315+
eta = tanh(smooth_coeff/min(x%min_spacing, &
316+
& y%min_spacing)*(sqrt((x%cc(i) - x_centroid)**2 + (y%cc(j) - y_centroid)**2) - radius))*(-0.5_wp) &
317+
& + 0.5_wp
317318
end if
318319
319320
if (((x%cc(i) - x_centroid)**2 + (y%cc(j) - y_centroid)**2 <= radius**2 &
@@ -485,9 +486,9 @@ contains
485486
do j = 0, n
486487
do i = 0, m
487488
if (patch_icpp(patch_id)%smoothen) then
488-
eta = tanh(smooth_coeff/min(dx, &
489-
& dy)*(sqrt(((x%cc(i) - x_centroid)/a)**2 + ((y%cc(j) - y_centroid)/b)**2) - 1._wp))*(-0.5_wp) &
490-
& + 0.5_wp
489+
eta = tanh(smooth_coeff/min(x%min_spacing, &
490+
& y%min_spacing)*(sqrt(((x%cc(i) - x_centroid)/a)**2 + ((y%cc(j) - y_centroid)/b)**2) - 1._wp)) &
491+
& *(-0.5_wp) + 0.5_wp
491492
end if
492493

493494
if ((((x%cc(i) - x_centroid)/a)**2 + ((y%cc(j) - y_centroid)/b)**2 <= 1._wp &
@@ -555,8 +556,8 @@ contains
555556
end if
556557
557558
if (patch_icpp(patch_id)%smoothen) then
558-
eta = tanh(smooth_coeff/min(dx, dy, &
559-
& dz)*(sqrt(((x%cc(i) - x_centroid)/a)**2 + ((cart_y - y_centroid)/b)**2 + ((cart_z &
559+
eta = tanh(smooth_coeff/min(x%min_spacing, y%min_spacing, &
560+
& z%min_spacing)*(sqrt(((x%cc(i) - x_centroid)/a)**2 + ((cart_y - y_centroid)/b)**2 + ((cart_z &
560561
& - z_centroid)/c)**2) - 1._wp))*(-0.5_wp) + 0.5_wp
561562
end if
562563
@@ -688,7 +689,8 @@ contains
688689
do j = 0, n
689690
do i = 0, m
690691
if (patch_icpp(patch_id)%smoothen) then
691-
eta = 5.e-1_wp + 5.e-1_wp*tanh(smooth_coeff/min(dx, dy)*(a*x%cc(i) + b*y%cc(j) + c)/sqrt(a**2 + b**2))
692+
eta = 5.e-1_wp + 5.e-1_wp*tanh(smooth_coeff/min(x%min_spacing, &
693+
& y%min_spacing)*(a*x%cc(i) + b*y%cc(j) + c)/sqrt(a**2 + b**2))
692694
end if
693695

694696
if ((a*x%cc(i) + b*y%cc(j) + c >= 0._wp .and. patch_icpp(patch_id)%alter_patch(patch_id_fp(i, j, &
@@ -878,7 +880,7 @@ contains
878880
end if
879881
end if
880882
if (patch_icpp(patch_id)%smoothen) then
881-
eta = 0.5_wp + 0.5_wp*tanh(smooth_coeff/min(dx, dy)*(R_boundary - r))
883+
eta = 0.5_wp + 0.5_wp*tanh(smooth_coeff/min(x%min_spacing, y%min_spacing)*(R_boundary - r))
882884
end if
883885
if ((r <= R_boundary .and. patch_icpp(patch_id)%alter_patch(patch_id_fp(i, j, 0))) .or. patch_id_fp(i, j, &
884886
& 0) == smooth_patch_id) then
@@ -940,7 +942,8 @@ contains
940942
end do
941943
end do
942944
if (patch_icpp(patch_id)%smoothen) then
943-
eta_local = 0.5_wp + 0.5_wp*tanh(smooth_coeff/min(dx, dy, dz)*(R_surface - r))
945+
eta_local = 0.5_wp + 0.5_wp*tanh(smooth_coeff/min(x%min_spacing, y%min_spacing, &
946+
& z%min_spacing)*(R_surface - r))
944947
end if
945948
if ((r <= R_surface .and. patch_icpp(patch_id)%alter_patch(patch_id_fp(i, j, k))) .or. patch_id_fp(i, j, &
946949
& k) == smooth_patch_id) then
@@ -998,9 +1001,9 @@ contains
9981001
end if
9991002
10001003
if (patch_icpp(patch_id)%smoothen) then
1001-
eta = tanh(smooth_coeff/min(dx, dy, &
1002-
& dz)*(sqrt((x%cc(i) - x_centroid)**2 + (cart_y - y_centroid)**2 + (cart_z - z_centroid)**2) &
1003-
& - radius))*(-0.5_wp) + 0.5_wp
1004+
eta = tanh(smooth_coeff/min(x%min_spacing, y%min_spacing, &
1005+
& z%min_spacing)*(sqrt((x%cc(i) - x_centroid)**2 + (cart_y - y_centroid)**2 + (cart_z &
1006+
& - z_centroid)**2) - radius))*(-0.5_wp) + 0.5_wp
10041007
end if
10051008
10061009
if ((((x%cc(i) - x_centroid)**2 + (cart_y - y_centroid)**2 + (cart_z - z_centroid)**2 <= radius**2) &
@@ -1146,17 +1149,17 @@ contains
11461149

11471150
if (patch_icpp(patch_id)%smoothen) then
11481151
if (.not. f_is_default(length_x)) then
1149-
eta = tanh(smooth_coeff/min(dy, &
1150-
& dz)*(sqrt((cart_y - y_centroid)**2 + (cart_z - z_centroid)**2) - radius))*(-0.5_wp) &
1151-
& + 0.5_wp
1152+
eta = tanh(smooth_coeff/min(y%min_spacing, &
1153+
& z%min_spacing)*(sqrt((cart_y - y_centroid)**2 + (cart_z - z_centroid)**2) - radius)) &
1154+
& *(-0.5_wp) + 0.5_wp
11521155
else if (.not. f_is_default(length_y)) then
1153-
eta = tanh(smooth_coeff/min(dx, &
1154-
& dz)*(sqrt((x%cc(i) - x_centroid)**2 + (cart_z - z_centroid)**2) - radius))*(-0.5_wp) &
1155-
& + 0.5_wp
1156+
eta = tanh(smooth_coeff/min(x%min_spacing, &
1157+
& z%min_spacing)*(sqrt((x%cc(i) - x_centroid)**2 + (cart_z - z_centroid)**2) - radius)) &
1158+
& *(-0.5_wp) + 0.5_wp
11561159
else
1157-
eta = tanh(smooth_coeff/min(dx, &
1158-
& dy)*(sqrt((x%cc(i) - x_centroid)**2 + (cart_y - y_centroid)**2) - radius))*(-0.5_wp) &
1159-
& + 0.5_wp
1160+
eta = tanh(smooth_coeff/min(x%min_spacing, &
1161+
& y%min_spacing)*(sqrt((x%cc(i) - x_centroid)**2 + (cart_y - y_centroid)**2) - radius)) &
1162+
& *(-0.5_wp) + 0.5_wp
11601163
end if
11611164
end if
11621165

@@ -1232,8 +1235,9 @@ contains
12321235
end if
12331236

12341237
if (patch_icpp(patch_id)%smoothen) then
1235-
eta = 5.e-1_wp + 5.e-1_wp*tanh(smooth_coeff/min(dx, dy, &
1236-
& dz)*(a*x%cc(i) + b*cart_y + c*cart_z + d)/sqrt(a**2 + b**2 + c**2))
1238+
eta = 5.e-1_wp + 5.e-1_wp*tanh(smooth_coeff/min(x%min_spacing, y%min_spacing, &
1239+
& z%min_spacing)*(a*x%cc(i) + b*cart_y + c*cart_z + d)/sqrt(a**2 + b**2 &
1240+
& + c**2))
12371241
end if
12381242

12391243
if ((a*x%cc(i) + b*cart_y + c*cart_z + d >= 0._wp .and. patch_icpp(patch_id)%alter_patch(patch_id_fp(i, j, &
@@ -1326,11 +1330,11 @@ contains
13261330
write (*, "(A, 3(2X, F20.10))") " > Cen:", (bbox%min(1:3) + bbox%max(1:3))/2._wp
13271331
write (*, "(A, 3(2X, F20.10))") " > Max:", bbox%max(1:3)
13281332

1329-
grid_mm(1,:) = (/minval(x%cc) - 0.e5_wp*dx, maxval(x%cc) + 0.e5_wp*dx/)
1330-
grid_mm(2,:) = (/minval(y%cc) - 0.e5_wp*dy, maxval(y%cc) + 0.e5_wp*dy/)
1333+
grid_mm(1,:) = (/minval(x%cc) - 0.e5_wp*x%min_spacing, maxval(x%cc) + 0.e5_wp*x%min_spacing/)
1334+
grid_mm(2,:) = (/minval(y%cc) - 0.e5_wp*y%min_spacing, maxval(y%cc) + 0.e5_wp*y%min_spacing/)
13311335

13321336
if (p > 0) then
1333-
grid_mm(3,:) = (/minval(z%cc) - 0.e5_wp*dz, maxval(z%cc) + 0.e5_wp*dz/)
1337+
grid_mm(3,:) = (/minval(z%cc) - 0.e5_wp*z%min_spacing, maxval(z%cc) + 0.e5_wp*z%min_spacing/)
13341338
else
13351339
grid_mm(3,:) = (/0._wp, 0._wp/)
13361340
end if
@@ -1356,7 +1360,7 @@ contains
13561360
point = f_convert_cyl_to_cart(point)
13571361
end if
13581362

1359-
eta = f_model_is_inside(model, point, (/dx, dy, dz/), patch_icpp(patch_id)%model_spc)
1363+
eta = f_model_is_inside(model, point, (/x%min_spacing, y%min_spacing, z%min_spacing/), patch_icpp(patch_id)%model_spc)
13601364

13611365
if (eta > patch_icpp(patch_id)%model_threshold) then
13621366
eta = 1._wp

src/pre_process/m_start_up.fpp

Lines changed: 12 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -181,8 +181,8 @@ contains
181181
182182
x%cc(0:m) = (x%cb(0:m) + x%cb(-1:(m - 1)))/2._wp
183183
184-
dx = minval(x%cb(0:m) - x%cb(-1:m - 1))
185-
if (num_procs > 1) call s_mpi_reduce_min(dx)
184+
x%min_spacing = minval(x%cb(0:m) - x%cb(-1:m - 1))
185+
if (num_procs > 1) call s_mpi_reduce_min(x%min_spacing)
186186
187187
x_domain%beg = x%cb(-1)
188188
x_domain%end = x%cb(m)
@@ -201,8 +201,8 @@ contains
201201
202202
y%cc(0:n) = (y%cb(0:n) + y%cb(-1:(n - 1)))/2._wp
203203
204-
dy = minval(y%cb(0:n) - y%cb(-1:n - 1))
205-
if (num_procs > 1) call s_mpi_reduce_min(dy)
204+
y%min_spacing = minval(y%cb(0:n) - y%cb(-1:n - 1))
205+
if (num_procs > 1) call s_mpi_reduce_min(y%min_spacing)
206206
207207
y_domain%beg = y%cb(-1)
208208
y_domain%end = y%cb(n)
@@ -221,8 +221,8 @@ contains
221221
222222
z%cc(0:p) = (z%cb(0:p) + z%cb(-1:(p - 1)))/2._wp
223223
224-
dz = minval(z%cb(0:p) - z%cb(-1:p - 1))
225-
if (num_procs > 1) call s_mpi_reduce_min(dz)
224+
z%min_spacing = minval(z%cb(0:p) - z%cb(-1:p - 1))
225+
if (num_procs > 1) call s_mpi_reduce_min(z%min_spacing)
226226
227227
z_domain%beg = z%cb(-1)
228228
z_domain%end = z%cb(p)
@@ -354,8 +354,8 @@ contains
354354
355355
x%cb(-1:m) = x_cb_glb((start_idx(1) - 1):(start_idx(1) + m))
356356
x%cc(0:m) = (x%cb(0:m) + x%cb(-1:(m - 1)))/2._wp
357-
dx = minval(x%cb(0:m) - x%cb(-1:(m - 1)))
358-
if (num_procs > 1) call s_mpi_reduce_min(dx)
357+
x%min_spacing = minval(x%cb(0:m) - x%cb(-1:(m - 1)))
358+
if (num_procs > 1) call s_mpi_reduce_min(x%min_spacing)
359359
x_domain%beg = x%cb(-1)
360360
x_domain%end = x%cb(m)
361361
@@ -374,8 +374,8 @@ contains
374374
375375
y%cb(-1:n) = y_cb_glb((start_idx(2) - 1):(start_idx(2) + n))
376376
y%cc(0:n) = (y%cb(0:n) + y%cb(-1:(n - 1)))/2._wp
377-
dy = minval(y%cb(0:n) - y%cb(-1:(n - 1)))
378-
if (num_procs > 1) call s_mpi_reduce_min(dy)
377+
y%min_spacing = minval(y%cb(0:n) - y%cb(-1:(n - 1)))
378+
if (num_procs > 1) call s_mpi_reduce_min(y%min_spacing)
379379
y_domain%beg = y%cb(-1)
380380
y_domain%end = y%cb(n)
381381
@@ -394,8 +394,8 @@ contains
394394
395395
z%cb(-1:p) = z_cb_glb((start_idx(3) - 1):(start_idx(3) + p))
396396
z%cc(0:p) = (z%cb(0:p) + z%cb(-1:(p - 1)))/2._wp
397-
dz = minval(z%cb(0:p) - z%cb(-1:(p - 1)))
398-
if (num_procs > 1) call s_mpi_reduce_min(dz)
397+
z%min_spacing = minval(z%cb(0:p) - z%cb(-1:(p - 1)))
398+
if (num_procs > 1) call s_mpi_reduce_min(z%min_spacing)
399399
z_domain%beg = z%cb(-1)
400400
z_domain%end = z%cb(p)
401401
end if

0 commit comments

Comments
 (0)