Skip to content

Commit 07a4ad6

Browse files
Merge pull request #12 from jacobwilliams/stack-overflow-fix
some changes to avoid stack overflow on windows for large arrays
2 parents 09abfdc + b474a13 commit 07a4ad6

1 file changed

Lines changed: 48 additions & 3 deletions

File tree

src/numerical_differentiation_module.f90

Lines changed: 48 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -222,6 +222,7 @@ module numerical_differentiation_module
222222

223223
procedure,public :: failed !! to check if an exception was raised.
224224
procedure,public :: get_error_status !! the status of error condition
225+
procedure,public :: set_dpert !! change the `dpert` values
225226

226227
! internal routines:
227228
procedure :: destroy_sparsity_pattern !! destroy the sparsity pattern
@@ -1260,6 +1261,26 @@ subroutine set_numdiff_sparsity_bounds(me,xlow,xhigh)
12601261
end subroutine set_numdiff_sparsity_bounds
12611262
!*******************************************************************************
12621263

1264+
!*******************************************************************************
1265+
!>
1266+
! Change the `dpert` vector. Can be used after the class has been initialized
1267+
! to change the perturbation step sizes (e.g., after an iteration).
1268+
1269+
subroutine set_dpert(me,dpert)
1270+
1271+
class(numdiff_type),intent(inout) :: me
1272+
real(wp),dimension(:),intent(in) :: dpert !! perturbation vector for `x`
1273+
1274+
if (size(dpert)/=me%n) then
1275+
call me%raise_exception(29,'set_dpert',&
1276+
'incorrect size of dpert array')
1277+
else
1278+
me%dpert = abs(dpert) ! update
1279+
end if
1280+
1281+
end subroutine set_dpert
1282+
!*******************************************************************************
1283+
12631284
!*******************************************************************************
12641285
!>
12651286
! Initialize a [[numdiff_type]] class. This must be called first.
@@ -1606,7 +1627,17 @@ subroutine columns_in_partition_group(me,igroup,n_cols,cols,nonzero_rows,indices
16061627
num_nonzero_elements_in_col
16071628
if (allocated(col_indices)) deallocate(col_indices)
16081629
allocate(col_indices(num_nonzero_elements_in_col))
1609-
col_indices = pack(me%indices,mask=me%icol==cols(i))
1630+
!col_indices = pack(me%indices,mask=me%icol==cols(i))
1631+
block
1632+
integer :: j,n
1633+
n = 0
1634+
do j = 1, size(me%icol)
1635+
if (me%icol(j)==cols(i)) then
1636+
n = n + 1
1637+
col_indices(n) = j
1638+
end if
1639+
end do
1640+
end block
16101641
if (allocated(nonzero_rows)) then
16111642
nonzero_rows = [nonzero_rows,me%irow(col_indices)]
16121643
indices = [indices,col_indices]
@@ -1652,7 +1683,10 @@ subroutine compute_indices(me)
16521683
integer :: i !! counter
16531684

16541685
allocate(me%indices(me%num_nonzero_elements))
1655-
me%indices = [(i,i=1,me%num_nonzero_elements)]
1686+
!me%indices = [(i,i=1,me%num_nonzero_elements)]
1687+
do i = 1, me%num_nonzero_elements
1688+
me%indices(i) = i
1689+
end do
16561690

16571691
end subroutine compute_indices
16581692
!*******************************************************************************
@@ -2780,7 +2814,18 @@ subroutine compute_jacobian_partitioned(me,x,dx,jac)
27802814
num_nonzero_elements_in_col = count(me%sparsity%icol==cols(i))
27812815
if (allocated(col_indices)) deallocate(col_indices)
27822816
allocate(col_indices(num_nonzero_elements_in_col))
2783-
col_indices = pack(me%sparsity%indices,mask=me%sparsity%icol==cols(i))
2817+
! col_indices = pack(me%sparsity%indices,mask=me%sparsity%icol==cols(i))
2818+
block
2819+
integer :: j,n
2820+
n = 0
2821+
do j = 1, size(me%sparsity%icol)
2822+
if (me%sparsity%icol(j)==cols(i)) then
2823+
n = n + 1
2824+
col_indices(n) = j
2825+
end if
2826+
end do
2827+
end block
2828+
27842829
df(me%sparsity%irow(col_indices)) = df(me%sparsity%irow(col_indices)) / &
27852830
(me%meth(1)%df_den_factor*dx(cols(i)))
27862831
end do

0 commit comments

Comments
 (0)