Skip to content
Snippets Groups Projects
iom.F90 143 KiB
Newer Older
Guillaume Samson's avatar
Guillaume Samson committed
MODULE iom
   !!======================================================================
   !!                    ***  MODULE  iom ***
   !! Input/Output manager :  Library to read input files
   !!======================================================================
   !! History :  2.0  ! 2005-12  (J. Belier) Original code
   !!            2.0  ! 2006-02  (S. Masson) Adaptation to NEMO
   !!            3.0  ! 2007-07  (D. Storkey) Changes to iom_gettime
   !!            3.4  ! 2012-12  (R. Bourdalle-Badie and G. Reffray)  add C1D case
   !!            3.6  ! 2014-15  DIMG format removed
   !!            3.6  ! 2015-15  (J. Harle) Added procedure to read REAL attributes
   !!            4.0  ! 2017-11  (M. Andrejczuk) Extend IOM interface to write any 3D fields
   !!----------------------------------------------------------------------

   !!----------------------------------------------------------------------
   !!   iom_open       : open a file read only
   !!   iom_close      : close a file or all files opened by iom
   !!   iom_get        : read a field (interfaced to several routines)
   !!   iom_varid      : get the id of a variable in a file
   !!   iom_rstput     : write a field in a restart file (interfaced to several routines)
   !!----------------------------------------------------------------------
   USE dom_oce         ! ocean space and time domain
   USE domutl          !
   USE flo_oce         ! floats module declarations
   USE lbclnk          ! lateal boundary condition / mpp exchanges
   USE iom_def         ! iom variables definitions
   USE iom_nf90        ! NetCDF format with native NetCDF library
   USE in_out_manager  ! I/O manager
   USE lib_mpp           ! MPP library
   USE sbc_oce  , ONLY :   nn_fsbc, ght_abl, ghw_abl, e3t_abl, e3w_abl, jpka, jpkam1
   USE icb_oce  , ONLY :   nclasses, class_num       !  !: iceberg classes
#if defined key_si3
   USE ice      , ONLY :   jpl
#endif
   USE phycst          ! physical constants
   USE dianam          ! build name of file
#if defined key_xios
   USE xios
# endif
   USE ioipsl, ONLY :  ju2ymds    ! for calendar
   USE crs             ! Grid coarsening
#if defined key_top
   USE trc, ONLY    :  profsed
#endif
   USE lib_fortran
   USE diu_bulk, ONLY : ln_diurnal_only, ln_diurnal
   USE iom_nf90
   USE netcdf

   IMPLICIT NONE
   PUBLIC   !   must be public to be able to access iom_def through iom

#if defined key_xios
   LOGICAL, PUBLIC, PARAMETER ::   lk_iomput = .TRUE.        !: iom_put flag
#else
   LOGICAL, PUBLIC, PARAMETER ::   lk_iomput = .FALSE.       !: iom_put flag
#endif
   LOGICAL, PUBLIC            ::   l_iom = .TRUE.            !: RK3 iom flag prevent writing at stage 1&2
Guillaume Samson's avatar
Guillaume Samson committed
   PUBLIC iom_init, iom_init_closedef, iom_swap, iom_open, iom_close, iom_setkt, iom_varid, iom_get, iom_get_var
   PUBLIC iom_chkatt, iom_getatt, iom_putatt, iom_getszuld, iom_rstput, iom_delay_rst, iom_put
   PUBLIC iom_use, iom_context_finalize, iom_update_file_name, iom_miss_val
   PUBLIC iom_xios_setid

   PRIVATE iom_rp0d_sp, iom_rp1d_sp, iom_rp2d_sp, iom_rp3d_sp
   PRIVATE iom_rp0d_dp, iom_rp1d_dp, iom_rp2d_dp, iom_rp3d_dp
   PRIVATE iom_get_123d
   PRIVATE iom_g0d_sp, iom_g1d_sp, iom_g2d_sp, iom_g3d_sp
   PRIVATE iom_g0d_dp, iom_g1d_dp, iom_g2d_dp, iom_g3d_dp
   PRIVATE iom_p1d_sp, iom_p2d_sp, iom_p3d_sp, iom_p4d_sp
   PRIVATE iom_p1d_dp, iom_p2d_dp, iom_p3d_dp, iom_p4d_dp
#if defined key_xios
   PRIVATE iom_set_domain_attr, iom_set_axis_attr, iom_set_field_attr, iom_set_file_attr, iom_get_file_attr, iom_set_grid_attr
   PRIVATE set_grid, set_grid_bounds, set_scalar, set_xmlatt, set_mooring, iom_sdate
   PRIVATE iom_set_rst_context, iom_set_vars_active
# endif
   PRIVATE set_xios_context
   PRIVATE iom_set_rstw_active

   INTERFACE iom_get
      MODULE PROCEDURE iom_g0d_sp, iom_g1d_sp, iom_g2d_sp, iom_g3d_sp
      MODULE PROCEDURE iom_g0d_dp, iom_g1d_dp, iom_g2d_dp, iom_g3d_dp
   END INTERFACE
   INTERFACE iom_getatt
      MODULE PROCEDURE iom_g0d_iatt, iom_g1d_iatt, iom_g0d_ratt, iom_g1d_ratt, iom_g0d_catt
   END INTERFACE
   INTERFACE iom_putatt
      MODULE PROCEDURE iom_p0d_iatt, iom_p1d_iatt, iom_p0d_ratt, iom_p1d_ratt, iom_p0d_catt
   END INTERFACE
   INTERFACE iom_rstput
      MODULE PROCEDURE iom_rp0d_sp, iom_rp1d_sp, iom_rp2d_sp, iom_rp3d_sp
      MODULE PROCEDURE iom_rp0d_dp, iom_rp1d_dp, iom_rp2d_dp, iom_rp3d_dp
   END INTERFACE
   INTERFACE iom_put
      MODULE PROCEDURE iom_p0d_sp, iom_p1d_sp, iom_p2d_sp, iom_p3d_sp, iom_p4d_sp
      MODULE PROCEDURE iom_p0d_dp, iom_p1d_dp, iom_p2d_dp, iom_p3d_dp, iom_p4d_dp
   END INTERFACE iom_put

   !! * Substitutions
#  include "do_loop_substitute.h90"
   !!----------------------------------------------------------------------
   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
   !! $Id: iom.F90 15512 2021-11-15 17:22:03Z techene $
Guillaume Samson's avatar
Guillaume Samson committed
   !! Software governed by the CeCILL license (see ./LICENSE)
   !!----------------------------------------------------------------------
CONTAINS

   SUBROUTINE iom_init( cdname, kdid, ld_closedef )
      !!----------------------------------------------------------------------
      !!                     ***  ROUTINE   ***
      !!
      !! ** Purpose :
      !!
      !!----------------------------------------------------------------------
      CHARACTER(len=*),           INTENT(in)  :: cdname
      INTEGER         , OPTIONAL, INTENT(in)  :: kdid
      LOGICAL         , OPTIONAL, INTENT(in)  :: ld_closedef
#if defined key_xios
      !
      TYPE(xios_duration) :: dtime    = xios_duration(0, 0, 0, 0, 0, 0)
      TYPE(xios_date)     :: start_date
Guillaume Samson's avatar
Guillaume Samson committed
      INTEGER             :: irefyear, irefmonth, irefday
      INTEGER           :: ji
      LOGICAL           :: llrst_context              ! is context related to restart
      LOGICAL           :: llrstr, llrstw
      INTEGER           :: inum
Guillaume Samson's avatar
Guillaume Samson committed
      !
      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zt_bnds, zw_bnds
      REAL(wp), DIMENSION(2,jpkam1)         :: za_bnds   ! ABL vertical boundaries
      LOGICAL ::   ll_closedef
      LOGICAL ::   ll_exist
      !!----------------------------------------------------------------------
      !
      ll_closedef = .TRUE.
      IF ( PRESENT(ld_closedef) ) ll_closedef = ld_closedef
      !
      ALLOCATE( zt_bnds(2,jpk), zw_bnds(2,jpk) )
      !
      clname = TRIM(cdname)
      IF ( .NOT. Agrif_Root() ) THEN
         iln    = INDEX(clname,'/', BACK=.TRUE.)
         cltmpn = clname(1:iln)
         clname = clname(iln+1:LEN_TRIM(clname))
         clname = TRIM(cltmpn)//TRIM(Agrif_CFixed())//'_'//TRIM(clname)
      ENDIF

Guillaume Samson's avatar
Guillaume Samson committed
      CALL xios_context_initialize(TRIM(clname), mpi_comm_oce)
      CALL iom_swap( cdname )

      llrstr = (cdname == cr_ocerst_cxt) .OR. (cdname == cr_icerst_cxt)
      llrstr = llrstr .OR. (cdname == cr_ablrst_cxt)
      llrstr = llrstr .OR. (cdname == cr_toprst_cxt)
      llrstr = llrstr .OR. (cdname == cr_sedrst_cxt)

      llrstw = (cdname == cw_ocerst_cxt) .OR. (cdname == cw_icerst_cxt)
      llrstw = llrstw .OR. (cdname == cw_ablrst_cxt)
      llrstw = llrstw .OR. (cdname == cw_toprst_cxt)
      llrstw = llrstw .OR. (cdname == cw_sedrst_cxt)

      llrst_context = llrstr .OR. llrstw

      ! Calendar type is now defined in xml file
      IF (.NOT.(xios_getvar('ref_year' ,irefyear ))) irefyear  = 1900
      IF (.NOT.(xios_getvar('ref_month',irefmonth))) irefmonth = 01
      IF (.NOT.(xios_getvar('ref_day'  ,irefday  ))) irefday   = 01

      SELECT CASE ( nleapy )        ! Choose calendar for IOIPSL
      CASE ( 1)   ;   CALL xios_define_calendar( TYPE = "Gregorian", time_origin = xios_date(irefyear,irefmonth,irefday,0,0,0),   &
          &                                                          start_date  = xios_date(   nyear,   nmonth,   nday,0,0,0) )
      CASE ( 0)   ;   CALL xios_define_calendar( TYPE = "NoLeap"   , time_origin = xios_date(irefyear,irefmonth,irefday,0,0,0),   &
          &                                                          start_date  = xios_date(   nyear,   nmonth,   nday,0,0,0) )
      CASE (30)   ;   CALL xios_define_calendar( TYPE = "D360"     , time_origin = xios_date(irefyear,irefmonth,irefday,0,0,0),   &
          &                                                          start_date  = xios_date(   nyear,   nmonth,   nday,0,0,0) )
      END SELECT

      ! horizontal grid definition
      IF(.NOT.llrst_context) CALL set_scalar
      !
      IF( cdname == cxios_context ) THEN
         CALL set_grid( "T", glamt, gphit, .FALSE., .FALSE. )
         CALL set_grid( "U", glamu, gphiu, .FALSE., .FALSE. )
         CALL set_grid( "V", glamv, gphiv, .FALSE., .FALSE. )
         CALL set_grid( "W", glamt, gphit, .FALSE., .FALSE. )
         CALL set_grid( "F", glamf, gphif, .FALSE., .FALSE. )
         CALL set_grid_znl( gphit )
         !
         IF( ln_cfmeta ) THEN   ! Add additional grid metadata
            CALL iom_set_domain_attr("grid_T", area = real( e1e2t(Nis0:Nie0, Njs0:Nje0), dp))
            CALL iom_set_domain_attr("grid_U", area = real( e1e2u(Nis0:Nie0, Njs0:Nje0), dp))
            CALL iom_set_domain_attr("grid_V", area = real( e1e2v(Nis0:Nie0, Njs0:Nje0), dp))
            CALL iom_set_domain_attr("grid_W", area = REAL( e1e2t(Nis0:Nie0, Njs0:Nje0), dp))
            CALL iom_set_domain_attr("grid_F", area = real( e1e2f(Nis0:Nie0, Njs0:Nje0), dp))
            CALL set_grid_bounds( "T", glamf, gphif, glamt, gphit )
            CALL set_grid_bounds( "U", glamv, gphiv, glamu, gphiu )
            CALL set_grid_bounds( "V", glamu, gphiu, glamv, gphiv )
            CALL set_grid_bounds( "W", glamf, gphif, glamt, gphit )
            CALL set_grid_bounds( "F", glamt, gphit, glamf, gphif )
         ENDIF
      ENDIF
      !
      IF( TRIM(cdname) == TRIM(cxios_context)//"_crs" ) THEN
         CALL dom_grid_crs   ! Save the parent grid information  & Switch to coarse grid domain
         !
         CALL set_grid( "T", glamt_crs, gphit_crs, .FALSE., .FALSE. )
         CALL set_grid( "U", glamu_crs, gphiu_crs, .FALSE., .FALSE. )
         CALL set_grid( "V", glamv_crs, gphiv_crs, .FALSE., .FALSE. )
         CALL set_grid( "W", glamt_crs, gphit_crs, .FALSE., .FALSE. )
         CALL set_grid_znl( gphit_crs )
          !
         CALL dom_grid_glo   ! Return to parent grid domain
         !
         IF( ln_cfmeta .AND. .NOT. llrst_context) THEN   ! Add additional grid metadata
            CALL iom_set_domain_attr("grid_T", area = real(e1e2t_crs(Nis0:Nie0, Njs0:Nje0), dp))
            CALL iom_set_domain_attr("grid_U", area = real(e1u_crs(Nis0:Nie0, Njs0:Nje0) * e2u_crs(Nis0:Nie0, Njs0:Nje0), dp))
            CALL iom_set_domain_attr("grid_V", area = real(e1v_crs(Nis0:Nie0, Njs0:Nje0) * e2v_crs(Nis0:Nie0, Njs0:Nje0), dp))
            CALL iom_set_domain_attr("grid_W", area = real(e1e2t_crs(Nis0:Nie0, Njs0:Nje0), dp))
            CALL set_grid_bounds( "T", glamf_crs, gphif_crs, glamt_crs, gphit_crs )
            CALL set_grid_bounds( "U", glamv_crs, gphiv_crs, glamu_crs, gphiu_crs )
            CALL set_grid_bounds( "V", glamu_crs, gphiu_crs, glamv_crs, gphiv_crs )
            CALL set_grid_bounds( "W", glamf_crs, gphif_crs, glamt_crs, gphit_crs )
         ENDIF
      ENDIF
      !
      ! vertical grid definition
      IF(.NOT.llrst_context) THEN
         CALL iom_set_axis_attr(  "deptht", paxis = gdept_1d )
         CALL iom_set_axis_attr(  "depthu", paxis = gdept_1d )
         CALL iom_set_axis_attr(  "depthv", paxis = gdept_1d )
         CALL iom_set_axis_attr(  "depthw", paxis = gdepw_1d )
          CALL iom_set_axis_attr(  "depthf", paxis = gdept_1d )

          ! ABL
         IF( .NOT. ALLOCATED(ght_abl) ) THEN   ! force definition for xml files (xios)
            ALLOCATE( ght_abl(jpka), ghw_abl(jpka), e3t_abl(jpka), e3w_abl(jpka) )   ! default allocation needed by iom
            ght_abl(:) = -1._wp   ;   ghw_abl(:) = -1._wp
            e3t_abl(:) = -1._wp   ;   e3w_abl(:) = -1._wp
         ENDIF
         CALL iom_set_axis_attr( "ght_abl", ght_abl(2:jpka) )
         CALL iom_set_axis_attr( "ghw_abl", ghw_abl(2:jpka) )

         ! Add vertical grid bounds
         zt_bnds(2,:      ) = gdept_1d(:)
         zt_bnds(1,2:jpk  ) = gdept_1d(1:jpkm1)
         zt_bnds(1,1      ) = gdept_1d(1) - e3w_1d(1)
         zw_bnds(1,:      ) = gdepw_1d(:)
         zw_bnds(2,1:jpkm1) = gdepw_1d(2:jpk)
         zw_bnds(2,jpk:   ) = gdepw_1d(jpk) + e3t_1d(jpk)
         CALL iom_set_axis_attr(  "deptht", bounds=zw_bnds )
         CALL iom_set_axis_attr(  "depthu", bounds=zw_bnds )
         CALL iom_set_axis_attr(  "depthv", bounds=zw_bnds )
         CALL iom_set_axis_attr(  "depthw", bounds=zt_bnds )
          CALL iom_set_axis_attr(  "depthf", bounds=zw_bnds )

         ! ABL
         za_bnds(1,:) = ghw_abl(1:jpkam1)
         za_bnds(2,:) = ghw_abl(2:jpka  )
         CALL iom_set_axis_attr( "ght_abl", bounds=za_bnds )
         za_bnds(1,:) = ght_abl(2:jpka  )
         za_bnds(2,:) = ght_abl(2:jpka  ) + e3w_abl(2:jpka)
         CALL iom_set_axis_attr( "ghw_abl", bounds=za_bnds )

         CALL iom_set_axis_attr(  "nfloat", (/ (REAL(ji,wp), ji=1,jpnfl) /) )
# if defined key_si3
         CALL iom_set_axis_attr( "ncatice", (/ (REAL(ji,wp), ji=1,jpl) /) )
         ! SIMIP diagnostics (4 main arctic straits)
         CALL iom_set_axis_attr( "nstrait", (/ (REAL(ji,wp), ji=1,4) /) )
# endif
#if defined key_top
         IF( ALLOCATED(profsed) ) CALL iom_set_axis_attr( "profsed", paxis = profsed )
#endif
         CALL iom_set_axis_attr( "icbcla", class_num )
         CALL iom_set_axis_attr( "iax_20C", (/ REAL(20,wp) /) )   ! strange syntaxe and idea...
         CALL iom_set_axis_attr( "iax_26C", (/ REAL(26,wp) /) )   ! strange syntaxe and idea...
         CALL iom_set_axis_attr( "iax_28C", (/ REAL(28,wp) /) )   ! strange syntaxe and idea...
         ! for diaprt, we need to define an axis which size can be 1 (default) or 5 (if the file subbasins.nc exists)
         INQUIRE( FILE = 'subbasins.nc', EXIST = ll_exist )
         nbasin = 1 + 4 * COUNT( (/ll_exist/) )
         CALL iom_set_axis_attr( "basin"  , (/ (REAL(ji,wp), ji=1,nbasin) /) )
      ENDIF
      !
      ! automatic definitions of some of the xml attributs
      IF(llrstr) THEN
         IF(PRESENT(kdid)) THEN
            CALL iom_set_rst_context(.TRUE.)
!set which fields will be read from restart file
            CALL iom_set_vars_active(kdid)
         ELSE
            CALL ctl_stop( 'iom_init:', 'restart read with XIOS: missing pointer to NETCDF file' )
         ENDIF
      ELSE IF(llrstw) THEN
         CALL iom_set_rstw_file(iom_file(kdid)%name)
      ELSE
         CALL set_xmlatt
      ENDIF
      !
      ! set time step length
      dtime%second = rn_Dt
      CALL xios_set_timestep( dtime )
      !
      ! conditional closure of context definition
      IF ( ll_closedef ) CALL iom_init_closedef
      !
      DEALLOCATE( zt_bnds, zw_bnds )
      !
#endif
      !
   END SUBROUTINE iom_init

   SUBROUTINE iom_init_closedef(cdname)
      !!----------------------------------------------------------------------
      !!            ***  SUBROUTINE iom_init_closedef  ***
      !!----------------------------------------------------------------------
      !!
      !! ** Purpose : Closure of context definition
      !!
      !!----------------------------------------------------------------------
      CHARACTER(len=*), OPTIONAL, INTENT(IN) :: cdname
#if defined key_xios
      LOGICAL :: llrstw

      llrstw = .FALSE.
      IF(PRESENT(cdname)) THEN
         llrstw = (cdname == cw_ocerst_cxt)
         llrstw = llrstw .OR. (cdname == cw_icerst_cxt)
         llrstw = llrstw .OR. (cdname == cw_ablrst_cxt)
         llrstw = llrstw .OR. (cdname == cw_toprst_cxt)
         llrstw = llrstw .OR. (cdname == cw_sedrst_cxt)
      ENDIF

      IF( llrstw ) THEN
!set names of the fields in restart file IF using XIOS to write data
         CALL iom_set_rst_context(.FALSE.)
         CALL xios_close_context_definition()
      ELSE
         CALL xios_close_context_definition()
         CALL xios_update_calendar( 0 )
      ENDIF
#else
      IF( .FALSE. )   WRITE(numout,*) 'iom_init_closedef: should not see this'   ! useless statement to avoid compilation warnings
#endif

   END SUBROUTINE iom_init_closedef

   SUBROUTINE iom_set_vars_active(idnum)
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE  iom_set_vars_active  ***
      !!
      !! ** Purpose :  define filename in XIOS context for reading file,
      !!               enable variables present in a file for reading with XIOS
      !!               id of the file is assumed to be rrestart.
      !!---------------------------------------------------------------------
      INTEGER, INTENT(IN) :: idnum

#if defined key_xios
      INTEGER                                    :: ndims, nvars, natts, unlimitedDimId, dimlen, xtype,mdims
      TYPE(xios_field)                           :: field_hdl
      TYPE(xios_file)                            :: file_hdl
      TYPE(xios_filegroup)                       :: filegroup_hdl
      INTEGER                                    :: dimids(4), jv,i, idim
      CHARACTER(LEN=256)                         :: clinfo               ! info character
      INTEGER, ALLOCATABLE                       :: indimlens(:)
      CHARACTER(LEN=nf90_max_name), ALLOCATABLE  :: indimnames(:)
      CHARACTER(LEN=nf90_max_name)               :: dimname, varname
      INTEGER                                    :: iln
      CHARACTER(LEN=lc)                          :: fname
      LOGICAL                                    :: lmeta
!metadata in restart file for restart read with XIOS
      INTEGER, PARAMETER                         :: NMETA = 11
      CHARACTER(LEN=lc)                          :: meta(NMETA)


      meta(1) = "nav_lat"
      meta(2) = "nav_lon"
      meta(3) = "nav_lev"
      meta(4) = "time_instant"
      meta(5) = "time_instant_bounds"
      meta(6) = "time_counter"
      meta(7) = "time_counter_bounds"
      meta(8) = "x"
      meta(9) = "y"
      meta(10) = "numcat"
      meta(11) = "nav_hgt"

      clinfo = '          iom_set_vars_active, file: '//TRIM(iom_file(idnum)%name)

      iln = INDEX( iom_file(idnum)%name, '.nc' )
!XIOS doee not need .nc
      IF(iln > 0) THEN
        fname =  iom_file(idnum)%name(1:iln-1)
      ELSE
        fname =  iom_file(idnum)%name
      ENDIF

!set name of the restart file and enable available fields
      CALL xios_get_handle("file_definition", filegroup_hdl )
      CALL xios_add_child(filegroup_hdl, file_hdl, 'rrestart')
      CALL xios_set_file_attr( "rrestart", name=fname, type="one_file",      &
           par_access="collective", enabled=.TRUE., mode="read",              &
                                                    output_freq=xios_timestep )

      CALL iom_nf90_check( nf90_inquire(iom_file(idnum)%nfid, ndims, nvars, natts ), clinfo )
      ALLOCATE(indimlens(ndims), indimnames(ndims))
      CALL iom_nf90_check( nf90_inquire(iom_file(idnum)%nfid, unlimitedDimId = unlimitedDimId ), clinfo )

      DO idim = 1, ndims
         CALL iom_nf90_check( nf90_inquire_dimension(iom_file(idnum)%nfid, idim, dimname, dimlen ), clinfo )
         indimlens(idim) = dimlen
         indimnames(idim) = dimname
      ENDDO

      DO jv =1, nvars
         lmeta = .FALSE.
         CALL iom_nf90_check( nf90_inquire_variable(iom_file(idnum)%nfid, jv, varname, xtype, ndims, dimids, natts ), clinfo )
         DO i = 1, NMETA
           IF(varname == meta(i)) THEN
             lmeta = .TRUE.
           ENDIF
         ENDDO
         IF(.NOT.lmeta) THEN
            CALL xios_add_child(file_hdl, field_hdl, varname)
            mdims = ndims

            IF(ANY(dimids(1:ndims) == unlimitedDimId)) THEN
               mdims = mdims - 1
            ENDIF

            IF(mdims == 3) THEN
               CALL xios_set_attr (field_hdl, enabled = .TRUE., name = varname,   &
                                   domain_ref="grid_N",                           &
                                   axis_ref=iom_axis(indimlens(dimids(mdims))),   &
                                   prec = 8, operation = "instant"                )
            ELSEIF(mdims == 2) THEN
               CALL xios_set_attr (field_hdl, enabled = .TRUE., name = varname,  &
                                   domain_ref="grid_N", prec = 8,                &
                                   operation = "instant"                         )
            ELSEIF(mdims == 1) THEN
               CALL xios_set_attr (field_hdl, enabled = .TRUE., name = varname, &
                                   axis_ref=iom_axis(indimlens(dimids(mdims))), &
                                   prec = 8, operation = "instant"              )
            ELSEIF(mdims == 0) THEN
               CALL xios_set_attr (field_hdl, enabled = .TRUE., name = varname, &
                                   scalar_ref = "grid_scalar", prec = 8,        &
                                   operation = "instant"                        )
            ELSE
               WRITE(ctmp1,*) 'iom_set_vars_active: variable ', TRIM(varname) ,' incorrect number of dimensions'
               CALL ctl_stop( 'iom_set_vars_active:', ctmp1 )
            ENDIF
         ENDIF
      ENDDO
      DEALLOCATE(indimlens, indimnames)
#endif
   END SUBROUTINE iom_set_vars_active

   SUBROUTINE iom_set_rstw_file(cdrst_file)
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE iom_set_rstw_file   ***
      !!
      !! ** Purpose :  define file name in XIOS context for writing restart
      !!---------------------------------------------------------------------
      CHARACTER(len=*) :: cdrst_file
#if defined key_xios
      TYPE(xios_file) :: file_hdl
      TYPE(xios_filegroup) :: filegroup_hdl

!set name of the restart file and enable available fields
      IF(lwp) WRITE(numout,*) 'Setting restart filename (for XIOS write) to: ', TRIM(cdrst_file)
      CALL xios_get_handle("file_definition", filegroup_hdl )
      CALL xios_add_child(filegroup_hdl, file_hdl, 'wrestart')
      IF(nxioso.eq.1) THEN
         CALL xios_set_file_attr( "wrestart", type="one_file", enabled=.TRUE.,&
                                       mode="write", output_freq=xios_timestep)
         IF(lwp) write(numout,*) 'OPEN ', TRIM(cdrst_file), ' in one_file mode'
Guillaume Samson's avatar
Guillaume Samson committed
      ELSE
         CALL xios_set_file_attr( "wrestart", type="multiple_file", enabled=.TRUE.,&
                                            mode="write", output_freq=xios_timestep)
         IF(lwp) write(numout,*) 'OPEN ', TRIM(cdrst_file), ' in multiple_file mode'
Guillaume Samson's avatar
Guillaume Samson committed
      ENDIF
      CALL xios_set_file_attr( "wrestart", name=TRIM(cdrst_file))
Guillaume Samson's avatar
Guillaume Samson committed
#endif
   END SUBROUTINE iom_set_rstw_file


   SUBROUTINE iom_set_rstw_active(sdfield, rd0, rs0, rd1, rs1, rd2, rs2, rd3, rs3)
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE iom_set_rstw_active   ***
      !!
      !! ** Purpose :  define file name in XIOS context for writing restart
      !!               enable variables present in restart file for writing
      !!---------------------------------------------------------------------
!sets enabled = .TRUE. for each field in restart file
      CHARACTER(len = *), INTENT(IN)                     :: sdfield
      REAL(dp), OPTIONAL, INTENT(IN)                     :: rd0
      REAL(sp), OPTIONAL, INTENT(IN)                     :: rs0
      REAL(dp), OPTIONAL, INTENT(IN), DIMENSION(:)       :: rd1
      REAL(sp), OPTIONAL, INTENT(IN), DIMENSION(:)       :: rs1
      REAL(dp), OPTIONAL, INTENT(IN), DIMENSION(:, :)    :: rd2
      REAL(sp), OPTIONAL, INTENT(IN), DIMENSION(:, :)    :: rs2
      REAL(dp), OPTIONAL, INTENT(IN), DIMENSION(:, :, :) :: rd3
      REAL(sp), OPTIONAL, INTENT(IN), DIMENSION(:, :, :) :: rs3
#if defined key_xios
      TYPE(xios_field) :: field_hdl
      TYPE(xios_file) :: file_hdl

      CALL xios_get_handle("wrestart", file_hdl)
!define fields for restart context
      CALL xios_add_child(file_hdl, field_hdl, sdfield)

      IF(PRESENT(rd3)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             domain_ref = "grid_N",                       &
                             axis_ref = iom_axis(size(rd3, 3)),           &
                             prec = 8, operation = "instant"              )
      ELSEIF(PRESENT(rs3)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             domain_ref = "grid_N",                       &
                             axis_ref = iom_axis(size(rd3, 3)),           &
                             prec = 4, operation = "instant"              )
      ELSEIF(PRESENT(rd2)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             domain_ref = "grid_N", prec = 8,             &
                             operation = "instant"                        )
      ELSEIF(PRESENT(rs2)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             domain_ref = "grid_N", prec = 4,             &
                             operation = "instant"                        )
      ELSEIF(PRESENT(rd1)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             axis_ref = iom_axis(size(rd1, 1)),           &
                             prec = 8, operation = "instant"              )
      ELSEIF(PRESENT(rs1)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             axis_ref = iom_axis(size(rd1, 1)),           &
                             prec = 4, operation = "instant"              )
      ELSEIF(PRESENT(rd0)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             scalar_ref = "grid_scalar", prec = 8,        &
                             operation = "instant"                        )
      ELSEIF(PRESENT(rs0)) THEN
         CALL xios_set_attr (field_hdl, enabled = .TRUE., name = sdfield, &
                             scalar_ref = "grid_scalar", prec = 4,        &
                             operation = "instant"                        )
      ENDIF
#endif
   END SUBROUTINE iom_set_rstw_active

   FUNCTION iom_axis(idlev) result(axis_ref)
      !!---------------------------------------------------------------------
      !!                   ***  FUNCTION  iom_axis  ***
      !!
      !! ** Purpose : Used for grid definition when XIOS is used to read/write
      !!              restart. Returns axis corresponding to the number of levels
      !!              given as an input variable. Axes are defined in routine
      !!              iom_set_rst_context
      !!---------------------------------------------------------------------
      INTEGER, INTENT(IN) :: idlev
      CHARACTER(len=lc)   :: axis_ref
      CHARACTER(len=12)   :: str
      IF(idlev == jpk) THEN
         axis_ref="nav_lev"
      ELSEIF(idlev == jpka) THEN
         axis_ref="nav_hgt"
#if defined key_si3
      ELSEIF(idlev == jpl) THEN
         axis_ref="numcat"
#endif
      ELSE
         write(str, *) idlev
         CALL ctl_stop( 'iom_axis', 'Definition for axis with '//TRIM(ADJUSTL(str))//' levels missing')
      ENDIF
   END FUNCTION iom_axis

   FUNCTION iom_xios_setid(cdname) result(kid)
     !!---------------------------------------------------------------------
      !!                   ***  FUNCTION    ***
      !!
      !! ** Purpose : this function returns first available id to keep information about file
      !!              sets filename in iom_file structure and sets name
      !!              of XIOS context depending on cdcomp
      !!              corresponds to iom_nf90_open
      !!---------------------------------------------------------------------
      CHARACTER(len=*), INTENT(in   ) :: cdname      ! File name
      INTEGER                         :: kid      ! identifier of the opened file
      INTEGER                         :: jl

      kid = 0
      DO jl = jpmax_files, 1, -1
         IF( iom_file(jl)%nfid == 0 )   kid = jl
      ENDDO

      iom_file(kid)%name   = TRIM(cdname)
      iom_file(kid)%nfid   = 1
      iom_file(kid)%nvars  = 0
      iom_file(kid)%irec   = -1

   END FUNCTION iom_xios_setid

   SUBROUTINE iom_set_rst_context(ld_rstr)
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE  iom_set_rst_context  ***
      !!
      !! ** Purpose : Define domain, axis and grid for restart (read/write)
      !!              context
      !!
      !!---------------------------------------------------------------------
      LOGICAL, INTENT(IN)               :: ld_rstr
      INTEGER :: ji
#if defined key_xios
      TYPE(xios_domaingroup)            :: domaingroup_hdl
      TYPE(xios_domain)                 :: domain_hdl
      TYPE(xios_axisgroup)              :: axisgroup_hdl
      TYPE(xios_axis)                   :: axis_hdl
      TYPE(xios_scalar)                 :: scalar_hdl
      TYPE(xios_scalargroup)            :: scalargroup_hdl

      CALL xios_get_handle("domain_definition",domaingroup_hdl)
      CALL xios_add_child(domaingroup_hdl, domain_hdl, "grid_N")
      CALL set_grid("N", glamt, gphit, .TRUE., ld_rstr)

      CALL xios_get_handle("axis_definition",axisgroup_hdl)
      CALL xios_add_child(axisgroup_hdl, axis_hdl, "nav_lev")
!AGRIF fails to compile when unit= is in call to xios_set_axis_attr
!     CALL xios_set_axis_attr( "nav_lev", long_name="Vertical levels",  unit="m", positive="down")
      CALL xios_set_axis_attr( "nav_lev", long_name = "Vertical levels in meters", positive = "down")
      CALL iom_set_axis_attr( "nav_lev", paxis = gdept_1d )
#if defined key_si3
      CALL xios_add_child(axisgroup_hdl, axis_hdl, "numcat")
      CALL iom_set_axis_attr( "numcat", (/ (REAL(ji,wp), ji=1,jpl) /) )
#endif
      CALL xios_add_child(axisgroup_hdl, axis_hdl, "nav_hgt")
      CALL iom_set_axis_attr( "nav_hgt", (/ (REAL(ji,wp), ji=1,jpka) /) )
      CALL xios_get_handle("scalar_definition", scalargroup_hdl)
      CALL xios_add_child(scalargroup_hdl, scalar_hdl, "grid_scalar")
#endif
   END SUBROUTINE iom_set_rst_context


   SUBROUTINE set_xios_context(kdid, cdcont)
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE  iom_set_rst_context  ***
      !!
      !! ** Purpose : set correct XIOS context based on kdid
      !!
      !!---------------------------------------------------------------------
      INTEGER,           INTENT(IN)     :: kdid           ! Identifier of the file
      CHARACTER(LEN=lc), INTENT(OUT)    :: cdcont         ! name of the context for XIOS read/write

      cdcont = "NONE"

      IF(lrxios) THEN
         IF(kdid == numror) THEN
            cdcont = cr_ocerst_cxt
         ELSEIF(kdid == numrir) THEN
            cdcont = cr_icerst_cxt
         ELSEIF(kdid == numrar) THEN
            cdcont = cr_ablrst_cxt
         ELSEIF(kdid == numrtr) THEN
            cdcont = cr_toprst_cxt
         ELSEIF(kdid == numrsr) THEN
            cdcont = cr_sedrst_cxt
         ENDIF
      ENDIF

      IF(lwxios) THEN
         IF(kdid == numrow) THEN
            cdcont = cw_ocerst_cxt
         ELSEIF(kdid == numriw) THEN
            cdcont = cw_icerst_cxt
         ELSEIF(kdid == numraw) THEN
            cdcont = cw_ablrst_cxt
         ELSEIF(kdid == numrtw) THEN
            cdcont = cw_toprst_cxt
         ELSEIF(kdid == numrsw) THEN
            cdcont = cw_sedrst_cxt
         ENDIF
      ENDIF
   END SUBROUTINE set_xios_context


   SUBROUTINE iom_swap( cdname )
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE  iom_swap  ***
      !!
      !! ** Purpose :  swap context between different agrif grid for xmlio_server
      !!---------------------------------------------------------------------
      CHARACTER(len=*), INTENT(in) :: cdname
      CHARACTER(len=256)           :: clname, cltmpn
      INTEGER                      :: iln
Guillaume Samson's avatar
Guillaume Samson committed
#if defined key_xios
      TYPE(xios_context) :: nemo_hdl

      clname = TRIM(cdname)
      IF ( .NOT. Agrif_Root() ) THEN
         iln    = INDEX(clname,'/', BACK=.TRUE.)
         cltmpn = clname(1:iln)
         clname = clname(iln+1:LEN_TRIM(clname))
         clname = TRIM(cltmpn)//TRIM(Agrif_CFixed())//'_'//TRIM(clname)
Guillaume Samson's avatar
Guillaume Samson committed
      ENDIF
      !
Guillaume Samson's avatar
Guillaume Samson committed
      CALL xios_set_current_context(nemo_hdl)
#endif
      !
   END SUBROUTINE iom_swap


   SUBROUTINE iom_open( cdname, kiomid, ldwrt, ldstop, ldiof, kdlev, cdcomp )
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE  iom_open  ***
      !!
      !! ** Purpose :  open an input file (return 0 if not found)
      !!---------------------------------------------------------------------
      CHARACTER(len=*), INTENT(in   )           ::   cdname   ! File name
      INTEGER         , INTENT(  out)           ::   kiomid   ! iom identifier of the opened file
      LOGICAL         , INTENT(in   ), OPTIONAL ::   ldwrt    ! open in write modeb          (default = .FALSE.)
      LOGICAL         , INTENT(in   ), OPTIONAL ::   ldstop   ! stop if open to read a non-existing file (default = .TRUE.)
      LOGICAL         , INTENT(in   ), OPTIONAL ::   ldiof    ! Interp On the Fly, needed for AGRIF (default = .FALSE.)
      INTEGER         , INTENT(in   ), OPTIONAL ::   kdlev    ! number of vertical levels
      CHARACTER(len=3), INTENT(in   ), OPTIONAL ::   cdcomp   ! name of component calling iom_nf90_open
      !
      CHARACTER(LEN=256)    ::   clname    ! the name of the file based on cdname [[+clcpu]+clcpu]
      CHARACTER(LEN=256)    ::   cltmpn    ! tempory name to store clname (in writting mode)
      CHARACTER(LEN=10)     ::   clsuffix  ! ".nc"
      CHARACTER(LEN=15)     ::   clcpu     ! the cpu number (max jpmax_digits digits)
      CHARACTER(LEN=256)    ::   clinfo    ! info character
      LOGICAL               ::   llok      ! check the existence
      LOGICAL               ::   llwrt     ! local definition of ldwrt
      LOGICAL               ::   llstop    ! local definition of ldstop
      LOGICAL               ::   lliof     ! local definition of ldiof
      INTEGER               ::   icnt      ! counter for digits in clcpu (max = jpmax_digits)
      INTEGER               ::   iln, ils  ! lengths of character
      INTEGER               ::   istop     !
      ! local number of points for x,y dimensions
      ! position of first local point for x,y dimensions
      ! position of last local point for x,y dimensions
      ! start halo size for x,y dimensions
      ! end halo size for x,y dimensions
      !---------------------------------------------------------------------
      ! Initializations and control
      ! =============
      kiomid = -1
      clinfo = '                    iom_open ~~~  '
      istop = nstop
      ! if iom_open is called for the first time: initialize iom_file(:)%nfid to 0
      ! (could be done when defining iom_file in f95 but not in f90)
      IF( Agrif_Root() ) THEN
         IF( iom_open_init == 0 ) THEN
            iom_file(:)%nfid = 0
            iom_open_init = 1
         ENDIF
      ENDIF
      ! do we read or write the file?
      IF( PRESENT(ldwrt) ) THEN   ;   llwrt = ldwrt
      ELSE                        ;   llwrt = .FALSE.
      ENDIF
      ! do we call ctl_stop if we try to open a non-existing file in read mode?
      IF( PRESENT(ldstop) ) THEN   ;   llstop = ldstop
      ELSE                         ;   llstop = .TRUE.
      ENDIF
      ! are we using interpolation on the fly?
      IF( PRESENT(ldiof) ) THEN   ;   lliof = ldiof
      ELSE                        ;   lliof = .FALSE.
      ENDIF
      ! create the file name by added, if needed, TRIM(Agrif_CFixed()) and TRIM(clsuffix)
      ! =============
Guillaume Samson's avatar
Guillaume Samson committed
      IF ( .NOT. Agrif_Root() .AND. .NOT. lliof ) THEN
         iln    = INDEX(clname,'/', BACK=.TRUE.)
Guillaume Samson's avatar
Guillaume Samson committed
         cltmpn = clname(1:iln)
         clname = clname(iln+1:LEN_TRIM(clname))
         clname = TRIM(cltmpn)//TRIM(Agrif_CFixed())//'_'//TRIM(clname)
Guillaume Samson's avatar
Guillaume Samson committed
      ENDIF
      ! which suffix should we use?
      clsuffix = '.nc'
      ! Add the suffix if needed
      iln = LEN_TRIM(clname)
      ils = LEN_TRIM(clsuffix)
      IF( iln <= ils .OR. INDEX( TRIM(clname), TRIM(clsuffix), back = .TRUE. ) /= iln - ils + 1 )   &
         &   clname = TRIM(clname)//TRIM(clsuffix)
      cltmpn = clname   ! store this name
      ! try to find if the file to be opened already exist
      ! =============
      INQUIRE( FILE = clname, EXIST = llok )
      IF( .NOT.llok ) THEN
         ! we try to add the cpu number to the name
         WRITE(clcpu,*) narea-1

         clcpu  = TRIM(ADJUSTL(clcpu))
         iln = INDEX(clname,TRIM(clsuffix), back = .TRUE.)
         clname = clname(1:iln-1)//'_'//TRIM(clcpu)//TRIM(clsuffix)
         icnt = 0
         INQUIRE( FILE = clname, EXIST = llok )
         ! we try different formats for the cpu number by adding 0
         DO WHILE( .NOT.llok .AND. icnt < jpmax_digits )
Guillaume Samson's avatar
Guillaume Samson committed
            clname = clname(1:iln-1)//'_'//TRIM(clcpu)//TRIM(clsuffix)
            INQUIRE( FILE = clname, EXIST = llok )
            icnt = icnt + 1
         END DO
      ELSE
         lxios_sini = .TRUE.
      ENDIF
      ! Open the NetCDF file
      ! =============
      ! do we have some free file identifier?
      IF( MINVAL(iom_file(:)%nfid) /= 0 )   &
         &   CALL ctl_stop( TRIM(clinfo), 'No more free file identifier', 'increase jpmax_files in iom_def' )
      ! if no file was found...
      IF( .NOT. llok ) THEN
         IF( .NOT. llwrt ) THEN   ! we are in read mode
            IF( llstop ) THEN   ;   CALL ctl_stop( TRIM(clinfo), 'File '//TRIM(cltmpn)//'* not found' )
            ELSE                ;   istop = nstop + 1   ! make sure that istop /= nstop so we don't open the file
            ENDIF
         ELSE                     ! we are in write mode so we
            clname = cltmpn       ! get back the file name without the cpu number
         ENDIF
      ELSE
         IF( llwrt .AND. .NOT. ln_clobber ) THEN   ! we stop as we want to write in a new file
            CALL ctl_stop( TRIM(clinfo), 'We want to write in a new file but '//TRIM(clname)//' already exists...' )
            istop = nstop + 1                      ! make sure that istop /= nstop so we don't open the file
         ELSEIF( llwrt ) THEN     ! the file exists and we are in write mode with permission to
            clname = cltmpn       ! overwrite so get back the file name without the cpu number
         ENDIF
      ENDIF
      IF( istop == nstop ) THEN   ! no error within this routine
         CALL iom_nf90_open( clname, kiomid, llwrt, llok, kdlev = kdlev, cdcomp = cdcomp )
      ENDIF
      !
   END SUBROUTINE iom_open


   SUBROUTINE iom_close( kiomid )
      !!--------------------------------------------------------------------
      !!                   ***  SUBROUTINE  iom_close  ***
      !!
      !! ** Purpose : close an input file, or all files opened by iom
      !!--------------------------------------------------------------------
      INTEGER, INTENT(inout), OPTIONAL ::   kiomid   ! iom identifier of the file to be closed
      !                                              ! return 0 when file is properly closed
      !                                              ! No argument: all files opened by iom are closed

      INTEGER ::   jf         ! dummy loop indices
      INTEGER ::   i_s, i_e   ! temporary integer
      CHARACTER(LEN=100)    ::   clinfo    ! info character
      !---------------------------------------------------------------------
      !
      IF( iom_open_init == 0 )   RETURN   ! avoid to use iom_file(jf)%nfid that us not yet initialized
      !
      clinfo = '                    iom_close ~~~  '
      IF( PRESENT(kiomid) ) THEN
         i_s = kiomid
         i_e = kiomid
      ELSE
         i_s = 1
         i_e = jpmax_files
      ENDIF

      IF( i_s > 0 ) THEN
         DO jf = i_s, i_e
            IF( iom_file(jf)%nfid > 0 ) THEN
               CALL iom_nf90_close( jf )
               iom_file(jf)%nfid       = 0          ! free the id
               IF( PRESENT(kiomid) )   kiomid = 0   ! return 0 as id to specify that the file was closed
               IF(lwp) WRITE(numout,*) TRIM(clinfo)//' close file: '//TRIM(iom_file(jf)%name)//' ok'
            ELSEIF( PRESENT(kiomid) ) THEN
               WRITE(ctmp1,*) '--->',  kiomid
               CALL ctl_stop( TRIM(clinfo)//' Invalid file identifier', ctmp1 )
            ENDIF
         END DO
      ENDIF
      !
   END SUBROUTINE iom_close


   FUNCTION iom_varid ( kiomid, cdvar, kdimsz, kndims, lduld, ldstop )
      !!-----------------------------------------------------------------------
      !!                  ***  FUNCTION  iom_varid  ***
      !!
      !! ** Purpose : get the id of a variable in a file (return 0 if not found)
      !!-----------------------------------------------------------------------
      INTEGER              , INTENT(in   )           ::   kiomid   ! file Identifier
      CHARACTER(len=*)     , INTENT(in   )           ::   cdvar    ! name of the variable
      INTEGER, DIMENSION(:), INTENT(  out), OPTIONAL ::   kdimsz   ! size of each dimension
      INTEGER              , INTENT(  out), OPTIONAL ::   kndims   ! number of dimensions
      LOGICAL              , INTENT(  out), OPTIONAL ::   lduld    ! true if the last dimension is unlimited (time)
      LOGICAL              , INTENT(in   ), OPTIONAL ::   ldstop   ! stop if looking for non-existing variable (default = .TRUE.)
      !
      INTEGER                        ::   iom_varid, iiv, i_nvd
      LOGICAL                        ::   ll_fnd
      CHARACTER(LEN=100)             ::   clinfo                   ! info character
      LOGICAL                        ::   llstop                   ! local definition of ldstop
      !!-----------------------------------------------------------------------
      iom_varid = 0                         ! default definition
      ! do we call ctl_stop if we look for non-existing variable?
      IF( PRESENT(ldstop) ) THEN   ;   llstop = ldstop
      ELSE                         ;   llstop = .TRUE.
      ENDIF
      !
      IF( kiomid > 0 ) THEN
         clinfo = 'iom_varid, file: '//TRIM(iom_file(kiomid)%name)//', var: '//TRIM(cdvar)
Guillaume Samson's avatar
Guillaume Samson committed
         IF( iom_file(kiomid)%nfid == 0 ) THEN
            CALL ctl_stop( TRIM(clinfo), 'the file is not open' )
Guillaume Samson's avatar
Guillaume Samson committed
         ELSE
            ll_fnd  = .FALSE.
            iiv = 0
            !
            DO WHILE ( .NOT.ll_fnd .AND. iiv < iom_file(kiomid)%nvars )
               iiv = iiv + 1
               ll_fnd  = ( TRIM(cdvar) == TRIM(iom_file(kiomid)%cn_var(iiv)) )
            END DO
            !
            IF( .NOT.ll_fnd ) THEN
               iiv = iiv + 1
               IF( iiv <= jpmax_vars ) THEN
                  iom_varid = iom_nf90_varid( kiomid, cdvar, iiv, kdimsz, kndims, lduld )
               ELSE
                  CALL ctl_stop( TRIM(clinfo), 'Too many variables in the file '//iom_file(kiomid)%name,   &
Guillaume Samson's avatar
Guillaume Samson committed
                        &                      'increase the parameter jpmax_vars')
               ENDIF
               IF( llstop .AND. iom_varid == -1 )   CALL ctl_stop( TRIM(clinfo)//' not found' )
            ELSE
               iom_varid = iiv
               IF( PRESENT(kdimsz) ) THEN
                  i_nvd = iom_file(kiomid)%ndims(iiv)
                  IF( i_nvd <= size(kdimsz) ) THEN
                     kdimsz(1:i_nvd) = iom_file(kiomid)%dimsz(1:i_nvd,iiv)
                  ELSE
                     WRITE(ctmp1,*) i_nvd, size(kdimsz)
                     CALL ctl_stop( TRIM(clinfo), 'error in kdimsz size'//TRIM(ctmp1) )
Guillaume Samson's avatar
Guillaume Samson committed
                  ENDIF
               ENDIF
               IF( PRESENT(kndims) )  kndims = iom_file(kiomid)%ndims(iiv)
               IF( PRESENT( lduld) )  lduld  = iom_file(kiomid)%luld( iiv)
            ENDIF
         ENDIF
      ENDIF
      !
   END FUNCTION iom_varid


   !!----------------------------------------------------------------------
   !!                   INTERFACE iom_get
   !!----------------------------------------------------------------------
   SUBROUTINE iom_g0d_sp( kiomid, cdvar, pvar, ktime )
      INTEGER         , INTENT(in   )                 ::   kiomid    ! Identifier of the file
      CHARACTER(len=*), INTENT(in   )                 ::   cdvar     ! Name of the variable
      REAL(sp)        , INTENT(  out)                 ::   pvar      ! read field
      REAL(dp)                                        ::   ztmp_pvar ! tmp var to read field
      INTEGER         , INTENT(in   ),     OPTIONAL   ::   ktime     ! record number
      !
      INTEGER                                         ::   idvar     ! variable id
      INTEGER                                         ::   idmspc    ! number of spatial dimensions
      INTEGER         , DIMENSION(1)                  ::   itime     ! record number
      CHARACTER(LEN=100)                              ::   clinfo    ! info character
      CHARACTER(LEN=100)                              ::   clname    ! file name
      CHARACTER(LEN=1)                                ::   cldmspc   !
      CHARACTER(LEN=lc)                               ::   context
      !
      CALL set_xios_context(kiomid, context)

      IF(context == "NONE") THEN  ! read data using default library
         itime = 1
         IF( PRESENT(ktime) ) itime = ktime
         !
         clname = iom_file(kiomid)%name
         clinfo = '          iom_g0d, file: '//TRIM(clname)//', var: '//TRIM(cdvar)
Guillaume Samson's avatar
Guillaume Samson committed
         !
         IF( kiomid > 0 ) THEN
            idvar = iom_varid( kiomid, cdvar )
            IF( iom_file(kiomid)%nfid > 0 .AND. idvar > 0 ) THEN
               idmspc = iom_file ( kiomid )%ndims( idvar )
               IF( iom_file(kiomid)%luld(idvar) )  idmspc = idmspc - 1
               WRITE(cldmspc , fmt='(i1)') idmspc
               IF( idmspc > 0 )  CALL ctl_stop( TRIM(clinfo), 'When reading to a 0D array, we do not accept data', &
                                    &                         'with 1 or more spatial dimensions: '//cldmspc//' were found.' , &
                                    &                         'Use ncwa -a to suppress the unnecessary dimensions' )
               CALL iom_nf90_get( kiomid, idvar, ztmp_pvar, itime )
               pvar = ztmp_pvar
            ENDIF
         ENDIF
      ELSE
#if defined key_xios
         IF(lwp) WRITE(numout,*) 'XIOS RST READ (0D): ', TRIM(cdvar)
Guillaume Samson's avatar
Guillaume Samson committed
         CALL iom_swap(context)
         CALL xios_recv_field( TRIM(cdvar), pvar)
Guillaume Samson's avatar
Guillaume Samson committed
         CALL iom_swap(cxios_context)
#else
         WRITE(ctmp1,*) 'Can not use XIOS in iom_g0d, file: '//TRIM(clname)//', var:'//TRIM(cdvar)
Guillaume Samson's avatar
Guillaume Samson committed
         CALL ctl_stop( 'iom_g0d', ctmp1 )
#endif
      ENDIF
   END SUBROUTINE iom_g0d_sp

   SUBROUTINE iom_g0d_dp( kiomid, cdvar, pvar, ktime )
      INTEGER         , INTENT(in   )                 ::   kiomid    ! Identifier of the file
      CHARACTER(len=*), INTENT(in   )                 ::   cdvar     ! Name of the variable
      REAL(dp)        , INTENT(  out)                 ::   pvar      ! read field
      INTEGER         , INTENT(in   ),     OPTIONAL   ::   ktime     ! record number
      !