diff --git a/src/LIB/SURCOUCHE/src/modd_field.f90 b/src/LIB/SURCOUCHE/src/modd_field.f90
index 277e8ad0825863688d2afe959f25982c74a21d25..8eecb7d3e17b39f41d784869d3f3837c880081d6 100644
--- a/src/LIB/SURCOUCHE/src/modd_field.f90
+++ b/src/LIB/SURCOUCHE/src/modd_field.f90
@@ -14,6 +14,7 @@
 !  P. Wautelet 10/11/2020: new data structures for netCDF dimensions
 !  P. Wautelet 24/09/2021: add Fill_tfielddata and use it as a custom constructor for tfielddata type
 !  P. Wautelet 08/10/2021: add 2 new dimensions: LW_bands (NMNHDIM_NLWB) and SW_bands (NMNHDIM_NSWB)
+!  P. Wautelet 14/10/2021: dynamically allocate tfieldlist (+ reallocate if necessary)
 !-----------------------------------------------------------------
 module modd_field
 
@@ -26,9 +27,11 @@ use NETCDF,          only: NF90_FILL_INT, NF90_FILL_REAL
 implicit none
 
 integer, parameter :: NMNHMAXDIMS = 6 ! Cannot be less than 6
-INTEGER,PARAMETER :: MAXFIELDS = 250
 INTEGER,PARAMETER :: TYPEUNDEF = -1, TYPEINT = 1, TYPELOG = 2, TYPEREAL = 3, TYPECHAR = 4, TYPEDATE = 5
 !
+INTEGER, PARAMETER :: NMAXFIELDINIT = 200 !Initial maximum number of fields in tfieldlist
+INTEGER, PARAMETER :: NMAXFIELDSTEP = 50  !Number of fields to add each time tfieldlist is too small
+
 integer, parameter :: NMNHDIM_UNKNOWN             = -2
 
 !For efficient use of memory, it is better that all values for real dimensions be contiguous
@@ -115,6 +118,8 @@ integer, dimension(0:8,3), parameter :: NMNHDIM_ARAKAWA = reshape( [ &
   NMNHDIM_NI_U,    NMNHDIM_NJ_V,    NMNHDIM_LEVEL_W] & ! fw point (=uvw point)
   , shape = [ 9, 3 ], order = [ 2, 1 ] )
 
+INTEGER, SAVE :: NMAXFIELDS !Maximum number of fields in tfieldlist (value is automatically increased if too small)
+
 TYPE TFIELDPTR_C0D
   CHARACTER(LEN=:),     POINTER :: DATA => NULL()
 END TYPE TFIELDPTR_C0D
@@ -216,7 +221,7 @@ TYPE, extends( tfield_metadata_base ) :: TFIELDDATA
   CHARACTER(LEN=4)   :: CLBTYPE   = 'NONE' !Type of the lateral boundary (LBX,LBY,LBXU,LBYV)
   LOGICAL            :: LTIMEDEP  = .FALSE. !Is the field time-dependent?
   !
-  INTEGER :: NMODELMAX = -1 !Number of models for which the field has been allocated
+  INTEGER :: NMODELMAX = -1 !Number of models for which the field has been allocated (default value must be negative)
   !
   TYPE(TFIELDPTR_C0D),DIMENSION(:),ALLOCATABLE :: TFIELD_C0D !Pointer to the character string fields (one per nested mesh)
   TYPE(TFIELDPTR_C1D),DIMENSION(:),ALLOCATABLE :: TFIELD_C1D !Pointer to the character string 1D fields (one per nested mesh)
@@ -243,7 +248,7 @@ END TYPE TFIELDDATA
 !
 integer, save :: NFIELDS_USED = 0
 LOGICAL, SAVE :: LFIELDLIST_ISINIT = .FALSE.
-TYPE(TFIELDDATA),DIMENSION(MAXFIELDS),SAVE :: TFIELDLIST
+TYPE(TFIELDDATA), ALLOCATABLE, DIMENSION(:), SAVE :: TFIELDLIST
 
 interface TFIELDDATA
   module procedure :: Fill_tfielddata
diff --git a/src/LIB/SURCOUCHE/src/mode_field.f90 b/src/LIB/SURCOUCHE/src/mode_field.f90
index 42b139e95cc2a9c32bf998c9be172fc26299c64e..326562a42b637fab2bd0889bed7b0960b192f98e 100644
--- a/src/LIB/SURCOUCHE/src/mode_field.f90
+++ b/src/LIB/SURCOUCHE/src/mode_field.f90
@@ -16,6 +16,7 @@
 !  P. Wautelet 23/01/2020: split in modd_field.f90 and mode_field.f90
 !  JL Redelsperger 03/2021: add variables for Ocean LES and auto-coupled version
 !  P. Wautelet 08/10/2021: add Goto_model_1field + Add_field2list procedures + remove Fieldlist_nmodel_resize
+!  P. Wautelet 14/10/2021: dynamically allocate tfieldlist (+ reallocate if necessary)
 !-----------------------------------------------------------------
 module mode_field
 
@@ -70,17 +71,14 @@ INTEGER,INTENT(IN),OPTIONAL :: KMODEL
 INTEGER :: IMODEL
 CHARACTER(LEN=42) :: YMSG
 !
-!F90/95: TFIELDLIST(1) = TFIELDDATA('UT','x_wind','m s-1','XY','X_Y_Z_U component of wind',2)
-!F2003:
-!TFIELDLIST(1) = TFIELDDATA(CMNHNAME='UT',CSTDNAME='x_wind',CUNITS='m s-1',CDIR='XY',&
-!                           CCOMMENT='X_Y_Z_U component of wind',NGRID=2)
-!
 CALL PRINT_MSG(NVERB_DEBUG,'GEN','INI_FIELD_LIST','called')
 IF (LFIELDLIST_ISINIT) THEN
   CALL PRINT_MSG(NVERB_ERROR,'GEN','INI_FIELD_LIST','already called')
   RETURN
 END IF
 LFIELDLIST_ISINIT = .TRUE.
+Allocate( tfieldlist(NMAXFIELDINIT) )
+NMAXFIELDS = NMAXFIELDINIT
 !
 IF (PRESENT(KMODEL)) THEN
   IMODEL = KMODEL
@@ -3098,7 +3096,7 @@ call Add_field2list( TFIELDDATA( &
   LTIMEDEP   = .FALSE.           ) )
 !
 !
-WRITE(YMSG,'("number of used fields=",I4," out of ",I4)') nfields_used-1,MAXFIELDS
+WRITE(YMSG,'("number of used fields=",I4," out of ",I4)') nfields_used-1,NMAXFIELDS
 CALL PRINT_MSG(NVERB_INFO,'GEN','INI_FIELD_LIST',TRIM(YMSG))
 !
 #if 0
@@ -3749,11 +3747,17 @@ implicit none
 
 type(tfielddata) :: tpfield
 
-CHARACTER(LEN=42) :: YMSG
-
-if ( nfields_used >= MAXFIELDS ) then
-  WRITE(YMSG,'( "nfields_used>=MAXFIELDS (",I5,")" )') MAXFIELDS
-  CALL PRINT_MSG(NVERB_FATAL,'GEN','INI_FIELD_LIST',TRIM(YMSG))
+character(len=64) :: ymsg
+type(tfielddata), allocatable, dimension(:) :: tzfieldlistnew
+
+!Check if tfieldlist big enough and enlarge it if necessary
+if ( nfields_used >= NMAXFIELDS ) then
+  Allocate( tzfieldlistnew(nmaxfields + nmaxfieldstep) )
+  tzfieldlistnew(1 : nmaxfields) = tfieldlist(1 : nmaxfields)
+  call Move_alloc( from = tzfieldlistnew, to = tfieldlist )
+  nmaxfields = nmaxfields + nmaxfieldstep
+  Write( ymsg, '( "nmaxfields increased from ", i5, " to ", i5 )') nmaxfields - nmaxfieldstep, nmaxfields
+  call Print_msg( NVERB_DEBUG, 'GEN', 'Add_field2list', Trim( ymsg ) )
 end if
 
 nfields_used = nfields_used + 1