Skip to content

Commit 25a074e

Browse files
sbryngelsonclaude
andauthored
Fix 6 safe pre/post-process bugs (batch) (#1240)
Co-authored-by: Claude Opus 4.6 <noreply@anthropic.com>
1 parent e331547 commit 25a074e

File tree

5 files changed

+15
-16
lines changed

5 files changed

+15
-16
lines changed

src/common/m_model.fpp

Lines changed: 4 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -241,14 +241,14 @@ contains
241241
end subroutine s_read_stl
242242

243243
!> This procedure reads an OBJ file.
244-
!! @param filepath Path to the odj file.
244+
!! @param filepath Path to the obj file.
245245
!! @param model The obj file.
246246
impure subroutine s_read_obj(filepath, model)
247247

248248
character(LEN=*), intent(in) :: filepath
249249
type(t_model), intent(out) :: model
250250

251-
integer :: i, j, k, l, iunit, iostat, nVertices
251+
integer :: i, j, k, l, iv3, iunit, iostat, nVertices
252252

253253
real(wp), dimension(1:3), allocatable :: vertices(:, :)
254254

@@ -297,10 +297,10 @@ contains
297297
read (line(3:), *) vertices(i, :)
298298
i = i + 1
299299
case ("f ")
300-
read (line(3:), *) k, l, j
300+
read (line(3:), *) k, l, iv3
301301
model%trs(j)%v(1, :) = vertices(k, :)
302302
model%trs(j)%v(2, :) = vertices(l, :)
303-
model%trs(j)%v(3, :) = vertices(j, :)
303+
model%trs(j)%v(3, :) = vertices(iv3, :)
304304
j = j + 1
305305
case default
306306
print *, "Error: unknown line type in OBJ file ", filepath

src/post_process/m_data_output.fpp

Lines changed: 5 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -1102,6 +1102,7 @@ contains
11021102
call MPI_BCAST(file_time, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
11031103
call MPI_BCAST(file_dt, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
11041104
call MPI_BCAST(file_num_procs, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
1105+
time_real = file_time
11051106
11061107
allocate (proc_bubble_counts(file_num_procs))
11071108
@@ -1271,6 +1272,7 @@ contains
12711272
call MPI_BCAST(file_time, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
12721273
call MPI_BCAST(file_dt, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
12731274
call MPI_BCAST(file_num_procs, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
1275+
time_real = file_time
12741276
12751277
allocate (proc_bubble_counts(file_num_procs))
12761278
@@ -1546,19 +1548,16 @@ contains
15461548
counter = counter + 1
15471549
x_d1(counter) = x_cc(j)
15481550
y_d1(counter) = y_cc(k)
1549-
euc_d = sqrt((x_cc(j) - x_d1(i))**2 + (y_cc(k) - y_d1(i))**2)
1550-
tgp = sqrt(dx(j)**2 + dy(k)**2)
15511551
else
1552-
euc_d = sqrt((x_cc(j) - x_d1(i))**2 + (y_cc(k) - y_d1(i))**2)
15531552
tgp = sqrt(dx(j)**2 + dy(k)**2)
15541553
do i = 1, counter
1554+
euc_d = sqrt((x_cc(j) - x_d1(i))**2 + (y_cc(k) - y_d1(i))**2)
15551555
if (euc_d < tgp) then
1556-
cycle
1557-
elseif (euc_d > tgp .and. i == counter) then
1556+
exit
1557+
elseif (i == counter) then
15581558
counter = counter + 1
15591559
x_d1(counter) = x_cc(j)
15601560
y_d1(counter) = y_cc(k)
1561-
15621561
end if
15631562
end do
15641563
end if

src/pre_process/m_simplex_noise.fpp

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -205,8 +205,8 @@ contains
205205
x2 = x0 - 1._wp + 2._wp*G2
206206
y2 = y0 - 1._wp + 2._wp*G2
207207

208-
ii = mod(i, 255)
209-
jj = mod(j, 255)
208+
ii = iand(i, 255)
209+
jj = iand(j, 255)
210210

211211
gi0 = mod(p_vec(ii + p_vec(jj)), 10) + 1
212212
gi1 = mod(p_vec(ii + i1 + p_vec(jj + j1)), 10) + 1

src/pre_process/m_start_up.fpp

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -357,7 +357,7 @@ contains
357357
! the time-step directory that will contain the new grid and initial
358358
! condition data are also generated.
359359
if (old_ic .neqv. .true.) then
360-
call s_delete_directory(trim(proc_rank_dir)//'/*')
360+
call s_delete_directory(trim(proc_rank_dir))
361361
call s_create_directory(trim(proc_rank_dir)//'/0')
362362
end if
363363
@@ -507,7 +507,7 @@ contains
507507
! process may be cleaned out to make room for new pre-process data.
508508
! In addition, the time-step folder that will contain the new grid
509509
! and initial condition data are also generated.
510-
call s_create_directory(trim(proc_rank_dir)//'/*')
510+
call s_delete_directory(trim(proc_rank_dir))
511511
call s_create_directory(trim(proc_rank_dir)//'/0')
512512
513513
end subroutine s_read_serial_ic_data_files

src/simulation/m_global_parameters.fpp

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -815,8 +815,8 @@ contains
815815
integral(i)%xmax = dflt_real
816816
integral(i)%ymin = dflt_real
817817
integral(i)%ymax = dflt_real
818-
integral(i)%ymin = dflt_real
819-
integral(i)%ymax = dflt_real
818+
integral(i)%zmin = dflt_real
819+
integral(i)%zmax = dflt_real
820820
end do
821821
822822
! GRCBC flags

0 commit comments

Comments
 (0)