diff --git a/Common/read_snapshot.f90 b/Common/read_snapshot.f90
index 9d71e3eb0317e1623768358e5a26083f986bc524..2e21b7a3dd63e801ae1913c77d5b431c98dfaa69 100644
--- a/Common/read_snapshot.f90
+++ b/Common/read_snapshot.f90
@@ -31,7 +31,7 @@ contains
 
     real, intent(in):: step(:) ! (2) longitude and latitude steps, in rad
 
-    integer, intent(in):: copy, ishape_last(:)
+    integer, intent(in):: copy, ishape_last(d_init:)
 
     ! Local:
     integer ishape, ishape_first
@@ -45,13 +45,13 @@ contains
        ishape_first = 0
     else
        ! {k > d_init}
-       ishape_first = ishape_last(k - d_init) + 1
+       ishape_first = ishape_last(k - 1) + 1
     end if
 
-    s%number_extr = ishape_last(k - d_init + 1) - ishape_first + 1
+    s%number_extr = ishape_last(k) - ishape_first + 1
     allocate(s%list(s%number_extr))
 
-    do ishape = ishape_first, ishape_last(k - d_init + 1)
+    do ishape = ishape_first, ishape_last(k)
        call read_eddy(e, date_read, eddy_i, hshp, ishape)
 
        ! Check that all the eddies have the same date index:
diff --git a/Overlap/Tests/test_get_dispatch_snap.f90 b/Overlap/Tests/test_get_dispatch_snap.f90
index 9fabee31537cd446a96d8d4b15f2838468294be6..91c25b54289ff603a64840301504cd0b2a7371ad 100644
--- a/Overlap/Tests/test_get_dispatch_snap.f90
+++ b/Overlap/Tests/test_get_dispatch_snap.f90
@@ -43,7 +43,7 @@ program test_get_dispatch_snap
 
   integer:: k = 0
 
-  integer, allocatable:: ishape_last(:)
+  integer, allocatable:: ishape_last(:) ! (d_init:d_init + n_dates - 1)
   ! shape index (0-based) in the collection of shapefiles of the last
   ! shape at a given date index
 
@@ -83,7 +83,7 @@ program test_get_dispatch_snap
      call dbf_read_attribute_03(d_init, hshp%extremum, hshp%extr_date, &
           ishape = 0)
      call read_column(ishape_last, &
-          file = trim(shpc_dir) // "/ishape_last.txt")
+          file = trim(shpc_dir) // "/ishape_last.txt", my_lbound = d_init)
      n_dates = size(ishape_last)
      call new_unit(unit_isolated)
      open(unit_isolated, file = "isolated_nodes.txt", status = "replace", &
@@ -97,7 +97,7 @@ program test_get_dispatch_snap
   call ezmpi_bcast(k, root = 0)
   call ezmpi_bcast(copy, root = 0)
   call ezmpi_bcast(n_dates, root = 0)
-  if (rank /= 0) allocate(ishape_last(n_dates))
+  if (rank /= 0) allocate(ishape_last(d_init:d_init + n_dates - 1))
   call ezmpi_bcast(ishape_last, root = 0)
   call ezmpi_bcast(d_init, root = 0)
   k_begin = d_init + (rank * n_dates) / n_proc
diff --git a/Overlap/Tests/test_overlap.f90 b/Overlap/Tests/test_overlap.f90
index 27d3f2db0425aad472435ded98eec983a926287c..71998889adcbb98ac0f5d49cd14538308e6a3cfd 100644
--- a/Overlap/Tests/test_overlap.f90
+++ b/Overlap/Tests/test_overlap.f90
@@ -68,11 +68,12 @@ program test_overlap
   allocate(flow(max_delta + 1))
   call shpc_open(hshp, trim(shpc_dir), rank = 0, pszaccess = "rb")
   call dbf_read_attribute_03(d_init, hshp%extremum, hshp%extr_date, ishape = 0)
-  call read_column(ishape_last, file = trim(shpc_dir) // "/ishape_last.txt")
+  call read_column(ishape_last, file = trim(shpc_dir) // "/ishape_last.txt", &
+       my_lbound = d_init)
   n_dates = size(ishape_last)
   call assert(d_init <= [k_test_1, k_test_2] .and. [k_test_1, k_test_2] &
        < d_init + n_dates, "test_overlap k_test_1, k_test_2")
-  e_overestim = maxval([ishape_last(1) + 1, ediff1d(ishape_last)])
+  e_overestim = maxval([ishape_last(d_init) + 1, ediff1d(ishape_last)])
   open(unit, file = "node_id_param.json", status = "replace", &
        action = "write")
   write(unit, fmt = *) '{"e_overestim": ', e_overestim, '}'
diff --git a/Overlap/Tests/test_send_recv.f90 b/Overlap/Tests/test_send_recv.f90
index fe160df0f3417b810853f617047baa91107727d4..d644bbab9bf875787b74c409784c8c0fd548bd2f 100644
--- a/Overlap/Tests/test_send_recv.f90
+++ b/Overlap/Tests/test_send_recv.f90
@@ -43,7 +43,7 @@ program test_send_recv
   logical periodic ! grid is periodic in longitude
   TYPE(shpc) hshp
 
-  integer, allocatable:: ishape_last(:)
+  integer, allocatable:: ishape_last(:) ! (d_init:d_init + n_dates - 1)
   ! shape index (0-based) in the collection of shapefiles of the last
   ! shape at a given date index
 
@@ -92,7 +92,7 @@ program test_send_recv
      call dbf_read_attribute_03(d_init, hshp%extremum, hshp%extr_date, &
           ishape = 0)
      call read_column(ishape_last, &
-          file = trim(shpc_dir) // "/ishape_last.txt")
+          file = trim(shpc_dir) // "/ishape_last.txt", my_lbound = d_init)
      n_dates = size(ishape_last)
   end if
 
@@ -102,7 +102,7 @@ program test_send_recv
   call ezmpi_bcast(nlat, root = 0)
   call ezmpi_bcast(copy, root = 0)
   call ezmpi_bcast(n_dates, root = 0)
-  if (rank == 1) allocate(ishape_last(n_dates))
+  if (rank == 1) allocate(ishape_last(d_init:d_init + n_dates - 1))
   call ezmpi_bcast(ishape_last, root = 0)
 
   call ezmpi_bcast(d_init, root = 0)
diff --git a/Overlap/eddy_graph.f90 b/Overlap/eddy_graph.f90
index 11d536ce7c41f56964d524c23f905799b1aa975e..55b0c00e5faffff324e551a3bb170ea22dbaa683 100644
--- a/Overlap/eddy_graph.f90
+++ b/Overlap/eddy_graph.f90
@@ -50,7 +50,7 @@ program eddy_graph
 
   logical periodic ! grid is periodic in longitude
 
-  integer, allocatable:: ishape_last(:)
+  integer, allocatable:: ishape_last(:) ! (d_init:d_init + n_dates - 1)
   ! shape index (0-based) in the collection of shapefiles of the last
   ! shape at a given date index
 
@@ -99,13 +99,14 @@ program eddy_graph
      call dbf_read_attribute_03(d_init, hshp%extremum, hshp%extr_date, &
           ishape = 0)
      call read_column(ishape_last, &
-          file = trim(shpc_dir) // "/ishape_last.txt", last = n_dates)
+          file = trim(shpc_dir) // "/ishape_last.txt", last = n_dates, &
+          my_lbound = d_init)
      if (n_dates == huge(0)) n_dates = size(ishape_last)
      call assert(n_dates >= max_delta + 1, &
           "eddy_graph: n_dates should be >= max_delta + 1")
      call assert(n_proc <= n_dates / (max_delta + 1), &
           "eddy_graph: n_proc should be <= n_dates / (max_delta + 1)")
-     e_overestim = maxval([ishape_last(1) + 1, ediff1d(ishape_last)])
+     e_overestim = maxval([ishape_last(d_init) + 1, ediff1d(ishape_last)])
      corner = corner_deg * deg_to_rad
      step = step_deg * deg_to_rad
   end if
@@ -119,7 +120,7 @@ program eddy_graph
   call ezmpi_bcast(periodic, root = 0)
   call ezmpi_bcast(copy, root = 0)
   call ezmpi_bcast(n_dates, root = 0)
-  if (rank /= 0) allocate(ishape_last(n_dates))
+  if (rank /= 0) allocate(ishape_last(d_init:d_init + n_dates - 1))
   call ezmpi_bcast(ishape_last, root = 0)
   call ezmpi_bcast(e_overestim, root = 0)
   call ezmpi_bcast(d_init, root = 0)
diff --git a/Overlap/get_snapshot.f90 b/Overlap/get_snapshot.f90
index b2fa3cf3182de0516007776272d791993eabece3..e1ed2c7938d09e7c4cd280fbea72e961e0de6ab0 100644
--- a/Overlap/get_snapshot.f90
+++ b/Overlap/get_snapshot.f90
@@ -13,8 +13,9 @@ contains
 
     type(snapshot), intent(out):: s
     integer, intent(in):: nlon, nlat
+    integer, intent(in):: d_init ! first date in the collection of shapefiles
     
-    integer, intent(in):: ishape_last(:)
+    integer, intent(in):: ishape_last(d_init:)
     ! shape index (0-based) in the collection of shapefiles of the last
     ! shape at a given date index
 
@@ -24,7 +25,6 @@ contains
     real, intent(in):: step(:) ! (2) longitude and latitude steps, in rad
     integer, intent(in):: copy
     TYPE(shpc), intent(in):: hshp
-    integer, intent(in):: d_init ! first date in the collection of shapefiles
     integer, intent(in):: k ! date index
     integer, intent(in):: k_end ! last date index analyzed by this MPI process
     integer, intent(in):: rank ! of MPI process