From 2e37f8af8a13f2a16da9e42e1b304f92896853f1 Mon Sep 17 00:00:00 2001
From: Philippe WAUTELET <philippe.wautelet@aero.obs-mip.fr>
Date: Mon, 26 Mar 2018 17:13:54 +0200
Subject: [PATCH] Philippe 26/03/2018: IO: time-dependent scalar (integer,
 logical or real) variables are now allowed

---
 src/LIB/SURCOUCHE/src/mode_netcdf.f90 | 80 +++++++++++++++++++--------
 1 file changed, 56 insertions(+), 24 deletions(-)

diff --git a/src/LIB/SURCOUCHE/src/mode_netcdf.f90 b/src/LIB/SURCOUCHE/src/mode_netcdf.f90
index 49f7fc09d..38ac393a8 100644
--- a/src/LIB/SURCOUCHE/src/mode_netcdf.f90
+++ b/src/LIB/SURCOUCHE/src/mode_netcdf.f90
@@ -966,7 +966,7 @@ TYPE(DIMCDF), POINTER :: PTDIM
 !
 CALL PRINT_MSG(NVERB_DEBUG,'IO','FILLVDIMS','called for '//TRIM(TPFIELD%CMNHNAME))
 !
-IF (SIZE(KSHAPE) < 1) CALL PRINT_MSG(NVERB_FATAL,'IO','FILLVDIMS','empty KSHAPE')
+IF (SIZE(KSHAPE) < 1 .AND. .NOT.TPFIELD%LTIMEDEP) CALL PRINT_MSG(NVERB_FATAL,'IO','FILLVDIMS','empty KSHAPE')
 !
 IGRID  =  TPFIELD%NGRID
 YDIR   =  TPFIELD%CDIR
@@ -1189,6 +1189,7 @@ INTEGER(KIND=IDCDF_KIND) :: STATUS
 INTEGER(KIND=IDCDF_KIND) :: INCID
 CHARACTER(LEN=LEN(TPFIELD%CMNHNAME)) :: YVARNAME
 INTEGER(KIND=IDCDF_KIND) :: IVARID
+INTEGER(KIND=IDCDF_KIND), DIMENSION(:), ALLOCATABLE :: IVDIMS
 INTEGER                  :: IRESP
 LOGICAL                  :: GEXISTED !True if variable was already defined
 !
@@ -1202,16 +1203,21 @@ GEXISTED = .FALSE.
 !
 CALL CLEANMNHNAME(TPFIELD%CMNHNAME,YVARNAME)
 !
-IF (TPFIELD%LTIMEDEP) &
-  CALL PRINT_MSG(NVERB_WARNING,'IO','IO_WRITE_FIELD_NC4_X0',TRIM(TPFILE%CNAME)// &
-                 ': time dependent variable not (yet) possible for 0D variable '//TRIM(TPFIELD%CMNHNAME))
-!
 ! The variable should not already exist but who knows ?
 STATUS = NF90_INQ_VARID(INCID, YVARNAME, IVARID)
 IF (STATUS /= NF90_NOERR) THEN
-   ! Define the scalar variable 
-   STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_DOUBLE, IVARID)
-   IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_X0[NF90_DEF_VAR]')
+   IF (TPFIELD%LTIMEDEP) THEN
+     ! Get the netcdf dimensions
+     CALL FILLVDIMS(TPFILE, TPFIELD, INT(SHAPE(PFIELD),KIND=IDCDF_KIND), IVDIMS)
+     ! Define the variable
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_DOUBLE,  IVDIMS, IVARID)
+     IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_X0[NF90_DEF_VAR]')
+     DEALLOCATE(IVDIMS)
+   ELSE
+     ! Define the scalar variable
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_DOUBLE, IVARID)
+     IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_X0[NF90_DEF_VAR]')
+   END IF
 ELSE
    GEXISTED = .TRUE.
    CALL PRINT_MSG(NVERB_WARNING,'IO','IO_WRITE_FIELD_NC4_X0',TRIM(TPFILE%CNAME)//': '//TRIM(YVARNAME)//' already defined')
@@ -1669,6 +1675,7 @@ INTEGER(KIND=IDCDF_KIND) :: STATUS
 INTEGER(KIND=IDCDF_KIND) :: INCID
 CHARACTER(LEN=LEN(TPFIELD%CMNHNAME)) :: YVARNAME
 INTEGER(KIND=IDCDF_KIND) :: IVARID
+INTEGER(KIND=IDCDF_KIND), DIMENSION(:), ALLOCATABLE :: IVDIMS
 INTEGER                  :: IRESP
 LOGICAL                  :: GEXISTED !True if variable was already defined
 !
@@ -1682,20 +1689,29 @@ GEXISTED = .FALSE.
 !
 CALL CLEANMNHNAME(TPFIELD%CMNHNAME,YVARNAME)
 !
-IF (TPFIELD%LTIMEDEP) &
-  CALL PRINT_MSG(NVERB_WARNING,'IO','IO_WRITE_FIELD_NC4_N0',TRIM(TPFILE%CNAME)// &
-                 ': time dependent variable not (yet) possible for 0D variable '//TRIM(TPFIELD%CMNHNAME))
-!
 ! The variable should not already exist but who knows ?
 STATUS = NF90_INQ_VARID(INCID, YVARNAME, IVARID)
 IF (STATUS /= NF90_NOERR) THEN
-   ! Define the scalar variable 
+   IF (TPFIELD%LTIMEDEP) THEN
+     ! Get the netcdf dimensions
+     CALL FILLVDIMS(TPFILE, TPFIELD, INT(SHAPE(KFIELD),KIND=IDCDF_KIND), IVDIMS)
+     ! Define the variable
 #if ( MNH_INT == 4 )
-   STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT, IVARID)
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT,   IVDIMS, IVARID)
 #else
-   STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT64, IVARID)
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT64, IVDIMS, IVARID)
 #endif
-   IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_N0[NF90_DEF_VAR]')
+     IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_N0[NF90_DEF_VAR]')
+     DEALLOCATE(IVDIMS)
+   ELSE
+     ! Define the scalar variable
+#if ( MNH_INT == 4 )
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT,   IVARID)
+#else
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT64, IVARID)
+#endif
+     IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_N0[NF90_DEF_VAR]')
+   END IF
 ELSE
    GEXISTED = .TRUE.
    CALL PRINT_MSG(NVERB_WARNING,'IO','IO_WRITE_FIELD_NC4_N0',TRIM(TPFILE%CNAME)//': '//TRIM(YVARNAME)//' already defined')
@@ -1941,6 +1957,7 @@ INTEGER(KIND=IDCDF_KIND) :: STATUS
 INTEGER(KIND=IDCDF_KIND) :: INCID
 CHARACTER(LEN=LEN(TPFIELD%CMNHNAME)) :: YVARNAME
 INTEGER(KIND=IDCDF_KIND) :: IVARID
+INTEGER(KIND=IDCDF_KIND), DIMENSION(:), ALLOCATABLE      :: IVDIMS
 INTEGER                  :: IRESP
 LOGICAL                  :: GEXISTED !True if variable was already defined
 !
@@ -1954,17 +1971,23 @@ GEXISTED = .FALSE.
 !
 CALL CLEANMNHNAME(TPFIELD%CMNHNAME,YVARNAME)
 !
-IF (TPFIELD%LTIMEDEP) &
-  CALL PRINT_MSG(NVERB_WARNING,'IO','IO_WRITE_FIELD_NC4_L0',TRIM(TPFILE%CNAME)// &
-                 ': time dependent variable not (yet) possible for 0D variable '//TRIM(TPFIELD%CMNHNAME))
-!
 ! The variable should not already exist but who knows ?
 STATUS = NF90_INQ_VARID(INCID, YVARNAME, IVARID)
 IF (STATUS /= NF90_NOERR) THEN
-   ! Define the scalar variable 
-   ! Use of NF90_INT1 datatype (=NF90_BYTE) that is enough to store a boolean
-   STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT1, IVARID)
-   IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_L0[NF90_DEF_VAR]')
+   IF (TPFIELD%LTIMEDEP) THEN
+     ! Get the netcdf dimensions
+     CALL FILLVDIMS(TPFILE, TPFIELD, INT(SHAPE(OFIELD),KIND=IDCDF_KIND), IVDIMS)
+     ! Define the variable
+     ! Use of NF90_INT1 datatype (=NF90_BYTE) that is enough to store a boolean
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT1, IVDIMS, IVARID)
+     IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_L0[NF90_DEF_VAR]')
+     DEALLOCATE(IVDIMS)
+   ELSE
+     ! Define the scalar variable
+     ! Use of NF90_INT1 datatype (=NF90_BYTE) that is enough to store a boolean
+     STATUS = NF90_DEF_VAR(INCID, YVARNAME, NF90_INT1, IVARID)
+     IF (status /= NF90_NOERR) CALL HANDLE_ERR(status,__LINE__,'IO_WRITE_FIELD_NC4_L0[NF90_DEF_VAR]')
+   END IF
 ELSE
    GEXISTED = .TRUE.
    CALL PRINT_MSG(NVERB_WARNING,'IO','IO_WRITE_FIELD_NC4_L0',TRIM(TPFILE%CNAME)//': '//TRIM(YVARNAME)//' already defined')
@@ -2542,6 +2565,9 @@ END IF
 STATUS = NF90_INQUIRE_VARIABLE(INCID, IVARID, XTYPE=ITYPE, NDIMS=IDIMS)
 IF (STATUS /= NF90_NOERR) CALL HANDLE_ERR(STATUS,__LINE__,'IO_READ_FIELD_NC4_X0[NF90_INQUIRE_VARIABLE] '//TRIM(YVARNAME))
 
+!Neglect the time dimension (of size 1)
+IF (TPFIELD%LTIMEDEP) IDIMS=IDIMS-1
+
 IF (IDIMS == 0 .AND. ITYPE == NF90_DOUBLE) THEN
   ! Read variable
   STATUS = NF90_GET_VAR(INCID, IVARID, PFIELD)
@@ -3044,6 +3070,9 @@ END IF
 STATUS = NF90_INQUIRE_VARIABLE(INCID, IVARID, XTYPE=ITYPE, NDIMS=IDIMS)
 IF (STATUS /= NF90_NOERR) CALL HANDLE_ERR(STATUS,__LINE__,'IO_READ_FIELD_NC4_N0[NF90_INQUIRE_VARIABLE] '//TRIM(YVARNAME))
 
+!Neglect the time dimension (of size 1)
+IF (TPFIELD%LTIMEDEP) IDIMS=IDIMS-1
+
 #if ( MNH_INT == 4 )
 IF (IDIMS == 0 .AND. (ITYPE == NF90_INT) ) THEN
 #else
@@ -3260,6 +3289,9 @@ END IF
 STATUS = NF90_INQUIRE_VARIABLE(INCID, IVARID, XTYPE=ITYPE, NDIMS=IDIMS)
 IF (STATUS /= NF90_NOERR) CALL HANDLE_ERR(STATUS,__LINE__,'IO_READ_FIELD_NC4_L0[NF90_INQUIRE_VARIABLE] '//TRIM(YVARNAME))
 
+!Neglect the time dimension (of size 1)
+IF (TPFIELD%LTIMEDEP) IDIMS=IDIMS-1
+
 !NF90_INT1 is for the case a boolean was written
 !Accept also INT and INT64 (for backward compatibility)
 IF (IDIMS == 0 .AND. (ITYPE == NF90_INT1 .OR. ITYPE == NF90_INT .OR. ITYPE == NF90_INT64)  ) THEN
-- 
GitLab