Skip to content
This repository was archived by the owner on Aug 28, 2024. It is now read-only.

Commit

Permalink
Fixed issue that was causing a segfault when an append operation was …
Browse files Browse the repository at this point in the history
…done using an unallocated source array
  • Loading branch information
daminton committed May 10, 2023
1 parent e1bc876 commit 678e812
Show file tree
Hide file tree
Showing 3 changed files with 23 additions and 9 deletions.
20 changes: 15 additions & 5 deletions src/base/base_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -295,12 +295,14 @@ subroutine base_util_append_arr_char_string(arr, source, nold, lsource_mask)
implicit none
! Arguments
character(len=STRMAX), dimension(:), allocatable, intent(inout) :: arr !! Destination array
character(len=STRMAX), dimension(:), intent(in) :: source !! Array to append
character(len=STRMAX), dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
! Internals
integer(I4B) :: nnew, nsrc, nend_orig

if (.not.allocated(source)) return

if (present(lsource_mask)) then
nsrc = count(lsource_mask(:))
else
Expand Down Expand Up @@ -338,12 +340,14 @@ subroutine base_util_append_arr_DP(arr, source, nold, lsource_mask)
implicit none
! Arguments
real(DP), dimension(:), allocatable, intent(inout) :: arr !! Destination array
real(DP), dimension(:), intent(in) :: source !! Array to append
real(DP), dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
! Internals
integer(I4B) :: nnew, nsrc, nend_orig

if (.not.allocated(source)) return

if (present(lsource_mask)) then
nsrc = count(lsource_mask(:))
else
Expand Down Expand Up @@ -381,12 +385,14 @@ subroutine base_util_append_arr_DPvec(arr, source, nold, lsource_mask)
implicit none
! Arguments
real(DP), dimension(:,:), allocatable, intent(inout) :: arr !! Destination array
real(DP), dimension(:,:), intent(in) :: source !! Array to append
real(DP), dimension(:,:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
! Internals
integer(I4B) :: nnew, nsrc, nend_orig

if (.not.allocated(source)) return

if (present(lsource_mask)) then
nsrc = count(lsource_mask(:))
else
Expand Down Expand Up @@ -426,12 +432,14 @@ subroutine base_util_append_arr_I4B(arr, source, nold, lsource_mask)
implicit none
! Arguments
integer(I4B), dimension(:), allocatable, intent(inout) :: arr !! Destination array
integer(I4B), dimension(:), intent(in) :: source !! Array to append
integer(I4B), dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
! Internals
integer(I4B) :: nnew, nsrc, nend_orig

if (.not.allocated(source)) return

if (present(lsource_mask)) then
nsrc = count(lsource_mask(:))
else
Expand Down Expand Up @@ -469,12 +477,14 @@ subroutine base_util_append_arr_logical(arr, source, nold, lsource_mask)
implicit none
! Arguments
logical, dimension(:), allocatable, intent(inout) :: arr !! Destination array
logical, dimension(:), intent(in) :: source !! Array to append
logical, dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
! Internals
integer(I4B) :: nnew, nsrc, nend_orig

if (.not.allocated(source)) return

if (present(lsource_mask)) then
nsrc = count(lsource_mask(:))
else
Expand Down
4 changes: 2 additions & 2 deletions src/swiftest/swiftest_module.f90
Original file line number Diff line number Diff line change
Expand Up @@ -1132,15 +1132,15 @@ end subroutine swiftest_user_kick_getacch_body
module subroutine swiftest_util_append_arr_info(arr, source, nold, lsource_mask)
implicit none
type(swiftest_particle_info), dimension(:), allocatable, intent(inout) :: arr !! Destination array
type(swiftest_particle_info), dimension(:), intent(in) :: source !! Array to append
type(swiftest_particle_info), dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
end subroutine swiftest_util_append_arr_info

module subroutine swiftest_util_append_arr_kin(arr, source, nold, lsource_mask)
implicit none
type(swiftest_kinship), dimension(:), allocatable, intent(inout) :: arr !! Destination array
type(swiftest_kinship), dimension(:), intent(in) :: source !! Array to append
type(swiftest_kinship), dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
end subroutine swiftest_util_append_arr_kin
Expand Down
8 changes: 6 additions & 2 deletions src/swiftest/swiftest_util.f90
Original file line number Diff line number Diff line change
Expand Up @@ -23,13 +23,15 @@ module subroutine swiftest_util_append_arr_info(arr, source, nold, lsource_mask)
implicit none
! Arguments
type(swiftest_particle_info), dimension(:), allocatable, intent(inout) :: arr !! Destination array
type(swiftest_particle_info), dimension(:), intent(in) :: source !! Array to append
type(swiftest_particle_info), dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
! Internals
integer(I4B) :: nnew, nsrc, nend_orig, i
integer(I4B), dimension(:), allocatable :: idx

if (.not.allocated(source)) return

if (present(lsource_mask)) then
nsrc = count(lsource_mask(:))
else
Expand Down Expand Up @@ -70,12 +72,14 @@ module subroutine swiftest_util_append_arr_kin(arr, source, nold, lsource_mask)
implicit none
! Arguments
type(swiftest_kinship), dimension(:), allocatable, intent(inout) :: arr !! Destination array
type(swiftest_kinship), dimension(:), intent(in) :: source !! Array to append
type(swiftest_kinship), dimension(:), allocatable, intent(in) :: source !! Array to append
integer(I4B), intent(in), optional :: nold !! Extent of original array. If passed, the source array will begin at arr(nold+1). Otherwise, the size of arr will be used.
logical, dimension(:), intent(in), optional :: lsource_mask !! Logical mask indicating which elements to append to
! Internals
integer(I4B) :: nnew, nsrc, nend_orig

if (.not.allocated(source)) return

if (present(lsource_mask)) then
nsrc = count(lsource_mask(:))
else
Expand Down

0 comments on commit 678e812

Please sign in to comment.