Skip to content
6 changes: 3 additions & 3 deletions src/common/m_model.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -226,7 +226,7 @@ contains
character(LEN=*), intent(in) :: filepath
type(t_model), intent(out) :: model

integer :: i, j, k, l, iunit, iostat, nVertices
integer :: i, j, k, l, iv3, iunit, iostat, nVertices

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

Expand Down Expand Up @@ -275,10 +275,10 @@ contains
read (line(3:), *) vertices(i, :)
i = i + 1
case ("f ")
read (line(3:), *) k, l, j
read (line(3:), *) k, l, iv3
model%trs(j)%v(1, :) = vertices(k, :)
model%trs(j)%v(2, :) = vertices(l, :)
model%trs(j)%v(3, :) = vertices(j, :)
model%trs(j)%v(3, :) = vertices(iv3, :)
j = j + 1
case default
print *, "Error: unknown line type in OBJ file ", filepath
Expand Down
11 changes: 5 additions & 6 deletions src/post_process/m_data_output.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -1102,6 +1102,7 @@ contains
call MPI_BCAST(file_time, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
call MPI_BCAST(file_dt, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
call MPI_BCAST(file_num_procs, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
time_real = file_time

allocate (proc_bubble_counts(file_num_procs))

Expand Down Expand Up @@ -1271,6 +1272,7 @@ contains
call MPI_BCAST(file_time, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
call MPI_BCAST(file_dt, 1, mpi_p, 0, MPI_COMM_WORLD, ierr)
call MPI_BCAST(file_num_procs, 1, MPI_INTEGER, 0, MPI_COMM_WORLD, ierr)
time_real = file_time

allocate (proc_bubble_counts(file_num_procs))

Expand Down Expand Up @@ -1546,19 +1548,16 @@ contains
counter = counter + 1
x_d1(counter) = x_cc(j)
y_d1(counter) = y_cc(k)
euc_d = sqrt((x_cc(j) - x_d1(i))**2 + (y_cc(k) - y_d1(i))**2)
tgp = sqrt(dx(j)**2 + dy(k)**2)
else
euc_d = sqrt((x_cc(j) - x_d1(i))**2 + (y_cc(k) - y_d1(i))**2)
tgp = sqrt(dx(j)**2 + dy(k)**2)
do i = 1, counter
euc_d = sqrt((x_cc(j) - x_d1(i))**2 + (y_cc(k) - y_d1(i))**2)
if (euc_d < tgp) then
cycle
elseif (euc_d > tgp .and. i == counter) then
exit
elseif (i == counter) then
counter = counter + 1
x_d1(counter) = x_cc(j)
y_d1(counter) = y_cc(k)

end if
end do
end if
Expand Down
4 changes: 2 additions & 2 deletions src/pre_process/m_simplex_noise.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -205,8 +205,8 @@ contains
x2 = x0 - 1._wp + 2._wp*G2
y2 = y0 - 1._wp + 2._wp*G2

ii = mod(i, 255)
jj = mod(j, 255)
ii = iand(i, 255)
jj = iand(j, 255)

gi0 = mod(p_vec(ii + p_vec(jj)), 10) + 1
gi1 = mod(p_vec(ii + i1 + p_vec(jj + j1)), 10) + 1
Expand Down
4 changes: 2 additions & 2 deletions src/pre_process/m_start_up.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -355,7 +355,7 @@ contains
! the time-step directory that will contain the new grid and initial
! condition data are also generated.
if (old_ic .neqv. .true.) then
call s_delete_directory(trim(proc_rank_dir)//'/*')
call s_delete_directory(trim(proc_rank_dir))
call s_create_directory(trim(proc_rank_dir)//'/0')
end if

Expand Down Expand Up @@ -505,7 +505,7 @@ contains
! process may be cleaned out to make room for new pre-process data.
! In addition, the time-step folder that will contain the new grid
! and initial condition data are also generated.
call s_create_directory(trim(proc_rank_dir)//'/*')
call s_delete_directory(trim(proc_rank_dir))
call s_create_directory(trim(proc_rank_dir)//'/0')

end subroutine s_read_serial_ic_data_files
Expand Down
4 changes: 2 additions & 2 deletions src/simulation/m_global_parameters.fpp
Original file line number Diff line number Diff line change
Expand Up @@ -814,8 +814,8 @@ contains
integral(i)%xmax = dflt_real
integral(i)%ymin = dflt_real
integral(i)%ymax = dflt_real
integral(i)%ymin = dflt_real
integral(i)%ymax = dflt_real
integral(i)%zmin = dflt_real
integral(i)%zmax = dflt_real
end do

! GRCBC flags
Expand Down
Loading