diff --git a/Common/CMakeLists.txt b/Common/CMakeLists.txt
index 5a51c07d790ee782650521ba7f154f3e8fccdc71..ad46dbebd64d2d3699edf234ea8a43d7955628a0 100644
--- a/Common/CMakeLists.txt
+++ b/Common/CMakeLists.txt
@@ -12,6 +12,8 @@ target_sources(test_spher_polyline_area PRIVATE spher_polyline_area.f90)
 target_sources(inst_eddies PRIVATE write_eddy.f90 spher_polyline_area.f90
   derived_types.f90 shpc_create.f90 shpc_close.f90 shpc_open.f90
   read_field_indices.f90)
+target_sources(test_write_eddy PRIVATE derived_types.f90 shpc_open.f90
+  shpc_close.f90 read_field_indices.f90 shpc_create.f90)
 
 if(MPI_Fortran_HAVE_F08_MODULE)
   target_sources(test_overlap PRIVATE derived_types.f90 read_snapshot.f90
diff --git a/Inst_eddies/Tests/CMakeLists.txt b/Inst_eddies/Tests/CMakeLists.txt
index cdb8de46b06fe08ac57f77ef03c6b4fca46e1814..34cf7391c7a3db8d75d9784ea0d5cb3bffaa74e5 100644
--- a/Inst_eddies/Tests/CMakeLists.txt
+++ b/Inst_eddies/Tests/CMakeLists.txt
@@ -52,10 +52,15 @@ target_link_libraries(test_set_max_speed PRIVATE Geometry::geometry
   Shapelib_03::shapelib_03 Contour_531::contour_531 Jumble::jumble
   NetCDF_Fortran::netcdff gpc_f)
 
+# test_write_eddy
+add_executable(test_write_eddy test_write_eddy.f90)
+target_link_libraries(test_write_eddy PRIVATE Shapelib_03::shapelib_03
+  Contour_531::contour_531 Jumble::jumble)
+
 foreach(my_target IN ITEMS test_get_1_outerm test_set_all_outerm
     test_good_contour test_inside_4 test_mean_speed
     test_max_speed_contour_ssh test_nearby_extr test_local_extrema
-    test_set_max_speed)
+    test_set_max_speed test_write_eddy)
   set_target_properties(${my_target} PROPERTIES Fortran_MODULE_DIRECTORY
     ${CMAKE_CURRENT_BINARY_DIR}/${my_target}_modules)
   target_include_directories(${my_target} PRIVATE
diff --git a/Inst_eddies/Tests/test_write_eddy.f90 b/Inst_eddies/Tests/test_write_eddy.f90
new file mode 100644
index 0000000000000000000000000000000000000000..64130aba87c62816066d4556b6640cd041ccaaf7
--- /dev/null
+++ b/Inst_eddies/Tests/test_write_eddy.f90
@@ -0,0 +1,108 @@
+program test_write_eddy
+
+  ! This program is a performance test for output of a shapefile
+  ! collection. It writes null shapes.
+
+  ! Libraries:
+  use shapelib_03, only: shp_append_point_03, dbf_write_attribute_03, &
+       shp_append_null_03
+
+  use derived_types, only: shpc
+  use shpc_close_m, only: shpc_close
+  use shpc_create_m, only: shpc_create
+  use shpc_open_m, only: shpc_open
+
+  implicit none
+
+  TYPE(shpc) hshpc
+  integer i, ishape
+  integer:: n_eddies = 50000
+  logical exist
+  logical:: mixed = .true.
+  namelist /main_nml/ n_eddies, mixed
+
+  !--------------------------------------------------------------
+
+  inquire(file = "SHPC_cyclo/extremum.shp", exist = exist)
+  
+  if (exist) then
+     call shpc_open(hshpc, shpc_dir = "SHPC_cyclo", pszaccess = "rb+")
+  else
+     call shpc_create(hshpc, shpc_dir = "SHPC_cyclo", cyclone = .true.)
+  end if
+
+  print *, "Enter namelist main_nml."
+  read(unit = *, nml = main_nml)
+
+  if (mixed) then
+     ! Alternate between shapefiles inside the loop on eddies:
+     do i = 1, n_eddies
+        call shp_append_point_03(ishape, hshpc%extremum, [0., 0.])
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_ssh, 0.)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_date, 0)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, &
+             hshpc%extr_eddy_index, i)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_valid, 1)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_speed, &
+             0.)
+        call shp_append_null_03(ishape, hshpc%outermost)
+        call shp_append_null_03(ishape, hshpc%max_speed)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, &
+             hshpc%out_r_eq_area, - 100.)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, hshpc%out_ssh, 0.)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, hshpc%out_date, 0)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, &
+             hshpc%out_eddy_index, i)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, &
+             hshpc%out_radius4, 0)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_r_eq_area, - 100.)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_ssh, 0.)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_date, 0)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_eddy_index, i)
+     end do
+  else
+     ! Loop on eddies for each shapefile:
+     
+     do i = 1, n_eddies
+        call shp_append_point_03(ishape, hshpc%extremum, [0., 0.])
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_ssh, 0.)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_date, 0)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, &
+             hshpc%extr_eddy_index, i)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_valid, 1)
+        call dbf_write_attribute_03(hshpc%extremum, ishape, hshpc%extr_speed, &
+             0.)
+     end do
+
+     do i = 1, n_eddies
+        call shp_append_null_03(ishape, hshpc%outermost)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, &
+             hshpc%out_r_eq_area, - 100.)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, hshpc%out_ssh, 0.)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, hshpc%out_date, 0)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, &
+             hshpc%out_eddy_index, i)
+        call dbf_write_attribute_03(hshpc%outermost, ishape, &
+             hshpc%out_radius4, 0)
+     end do
+
+     do i = 1, n_eddies
+        call shp_append_null_03(ishape, hshpc%max_speed)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_r_eq_area, - 100.)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_ssh, 0.)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_date, 0)
+        call dbf_write_attribute_03(hshpc%max_speed, ishape, &
+             hshpc%max_speed_eddy_index, i)
+     end do
+  end if
+
+  CALL shpc_close(hshpc)
+
+end program test_write_eddy