Skip to content

Commit 0545695

Browse files
sbryngelsonclaude
andcommitted
Fix bugs detected by new Fortran/Fypp static analysis linter
- m_mpi_proxy.fpp: fix duplicate 'bc_x%ve2' → 'bc_x%ve3' in MPI broadcast Fypp list; bc_x%ve3 was never broadcast to all ranks - m_assign_variables.fpp: remove duplicate R3bar accumulation line in QBMM path; bubble radius cubed was summed twice per iteration - m_rhs.fpp: fix cylindrical viscous boundary call passing dq_prim_dy_vf twice; 4th arg should be dq_prim_dz_vf (z-gradient) - m_data_input.f90, m_data_output.fpp, m_start_up.fpp (all targets): replace hardcoded WP_MOK = int(8._wp,...) with storage_size(0._stp)/8 for correct mixed-precision MPI-IO offset computation Co-Authored-By: Claude Sonnet 4.6 <noreply@anthropic.com>
1 parent 10e0984 commit 0545695

8 files changed

Lines changed: 12 additions & 13 deletions

File tree

src/post_process/m_data_input.f90

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -133,7 +133,7 @@ impure subroutine s_setup_mpi_io_params(data_size, m_MOK, n_MOK, p_MOK, WP_MOK,
133133
m_MOK = int(m_glb + 1, MPI_OFFSET_KIND)
134134
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
135135
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
136-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
136+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
137137
MOK = int(1._wp, MPI_OFFSET_KIND)
138138
str_MOK = int(name_len, MPI_OFFSET_KIND)
139139
NVARS_MOK = int(sys_size, MPI_OFFSET_KIND)
@@ -177,7 +177,7 @@ impure subroutine s_read_ib_data_files(file_loc_base, t_step)
177177
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
178178
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
179179
MOK = int(1._wp, MPI_OFFSET_KIND)
180-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
180+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
181181
save_index = t_step/t_step_save ! get the number of saves done to this point
182182

183183
data_size = (m + 1)*(n + 1)*(p + 1)
@@ -517,7 +517,7 @@ impure subroutine s_read_parallel_conservative_data(t_step, m_MOK, n_MOK, p_MOK,
517517
m_MOK = int(m_glb + 1, MPI_OFFSET_KIND)
518518
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
519519
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
520-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
520+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
521521
MOK = int(1._wp, MPI_OFFSET_KIND)
522522
str_MOK = int(name_len, MPI_OFFSET_KIND)
523523
NVARS_MOK = int(sys_size, MPI_OFFSET_KIND)

src/pre_process/m_assign_variables.fpp

Lines changed: 0 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -233,7 +233,6 @@ contains
233233
if (qbmm) then
234234
do i = 1, nb
235235
R3bar = R3bar + weight(i)*0.5_wp*(q_prim_vf(bubxb + 1 + (i - 1)*nmom)%sf(j, k, l))**3._wp
236-
R3bar = R3bar + weight(i)*0.5_wp*(q_prim_vf(bubxb + 1 + (i - 1)*nmom)%sf(j, k, l))**3._wp
237236
end do
238237
else
239238
do i = 1, nb

src/pre_process/m_data_output.fpp

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -549,7 +549,7 @@ contains
549549
m_MOK = int(m_glb_save, MPI_OFFSET_KIND)
550550
n_MOK = int(n_glb_save, MPI_OFFSET_KIND)
551551
p_MOK = int(p_glb_save, MPI_OFFSET_KIND)
552-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
552+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
553553
MOK = int(1._wp, MPI_OFFSET_KIND)
554554
str_MOK = int(name_len, MPI_OFFSET_KIND)
555555
NVARS_MOK = int(sys_size, MPI_OFFSET_KIND)
@@ -615,7 +615,7 @@ contains
615615
m_MOK = int(m_glb + 1, MPI_OFFSET_KIND)
616616
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
617617
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
618-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
618+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
619619
MOK = int(1._wp, MPI_OFFSET_KIND)
620620
str_MOK = int(name_len, MPI_OFFSET_KIND)
621621
NVARS_MOK = int(sys_size, MPI_OFFSET_KIND)

src/pre_process/m_start_up.fpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -660,7 +660,7 @@ contains
660660
m_MOK = int(m_glb + 1, MPI_OFFSET_KIND)
661661
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
662662
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
663-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
663+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
664664
MOK = int(1._wp, MPI_OFFSET_KIND)
665665
str_MOK = int(name_len, MPI_OFFSET_KIND)
666666
NVARS_MOK = int(sys_size, MPI_OFFSET_KIND)

src/simulation/m_data_output.fpp

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -895,7 +895,7 @@ contains
895895
m_MOK = int(m_glb_save + 1, MPI_OFFSET_KIND)
896896
n_MOK = int(n_glb_save + 1, MPI_OFFSET_KIND)
897897
p_MOK = int(p_glb_save + 1, MPI_OFFSET_KIND)
898-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
898+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
899899
MOK = int(1._wp, MPI_OFFSET_KIND)
900900
str_MOK = int(name_len, MPI_OFFSET_KIND)
901901
NVARS_MOK = int(sys_size, MPI_OFFSET_KIND)
@@ -963,7 +963,7 @@ contains
963963
m_MOK = int(m_glb + 1, MPI_OFFSET_KIND)
964964
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
965965
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
966-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
966+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
967967
MOK = int(1._wp, MPI_OFFSET_KIND)
968968
str_MOK = int(name_len, MPI_OFFSET_KIND)
969969
NVARS_MOK = int(alt_sys, MPI_OFFSET_KIND)
@@ -1087,7 +1087,7 @@ contains
10871087
m_MOK = int(m_glb + 1, MPI_OFFSET_KIND)
10881088
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
10891089
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
1090-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
1090+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
10911091
MOK = int(1._wp, MPI_OFFSET_KIND)
10921092

10931093
write (file_loc, '(A)') 'ib.dat'

src/simulation/m_mpi_proxy.fpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -145,7 +145,7 @@ contains
145145

146146
#:for VAR in [ 'dt','weno_eps','teno_CT','pref','rhoref','R0ref','Web','Ca', 'sigma', &
147147
& 'Re_inv', 'poly_sigma', 'palpha_eps', 'ptgalpha_eps', 'pi_fac', &
148-
& 'bc_x%vb1','bc_x%vb2','bc_x%vb3','bc_x%ve1','bc_x%ve2','bc_x%ve2', &
148+
& 'bc_x%vb1','bc_x%vb2','bc_x%vb3','bc_x%ve1','bc_x%ve2','bc_x%ve3', &
149149
& 'bc_y%vb1','bc_y%vb2','bc_y%vb3','bc_y%ve1','bc_y%ve2','bc_y%ve3', &
150150
& 'bc_z%vb1','bc_z%vb2','bc_z%vb3','bc_z%ve1','bc_z%ve2','bc_z%ve3', &
151151
& 'bc_x%pres_in','bc_x%pres_out','bc_y%pres_in','bc_y%pres_out', 'bc_z%pres_in','bc_z%pres_out', &

src/simulation/m_rhs.fpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -1653,7 +1653,7 @@ contains
16531653
call s_compute_viscous_stress_cylindrical_boundary(q_prim_vf, &
16541654
dq_prim_dx_vf(mom_idx%beg:mom_idx%end), &
16551655
dq_prim_dy_vf(mom_idx%beg:mom_idx%end), &
1656-
dq_prim_dy_vf(mom_idx%beg:mom_idx%end), &
1656+
dq_prim_dz_vf(mom_idx%beg:mom_idx%end), &
16571657
tau_Re_vf, &
16581658
idwbuff(1), idwbuff(2), idwbuff(3))
16591659
end if

src/simulation/m_start_up.fpp

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -652,7 +652,7 @@ contains
652652
m_MOK = int(m_glb + 1, MPI_OFFSET_KIND)
653653
n_MOK = int(n_glb + 1, MPI_OFFSET_KIND)
654654
p_MOK = int(p_glb + 1, MPI_OFFSET_KIND)
655-
WP_MOK = int(8._wp, MPI_OFFSET_KIND)
655+
WP_MOK = int(storage_size(0._stp)/8, MPI_OFFSET_KIND)
656656
MOK = int(1._wp, MPI_OFFSET_KIND)
657657
str_MOK = int(name_len, MPI_OFFSET_KIND)
658658
NVARS_MOK = int(sys_size, MPI_OFFSET_KIND)

0 commit comments

Comments
 (0)