New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
iom.F90 in NEMO/branches/2019/dev_r11756_SI3restart_XIOS/src/OCE/IOM – NEMO

source: NEMO/branches/2019/dev_r11756_SI3restart_XIOS/src/OCE/IOM/iom.F90 @ 11837

Last change on this file since 11837 was 11837, checked in by andmirek, 4 years ago

ticket #2323 read SI3 restart with XIOS

  • Property svn:keywords set to Id
File size: 140.8 KB
Line 
1MODULE iom
2   !!======================================================================
3   !!                    ***  MODULE  iom ***
4   !! Input/Output manager :  Library to read input files
5   !!======================================================================
6   !! History :  2.0  ! 2005-12  (J. Belier) Original code
7   !!            2.0  ! 2006-02  (S. Masson) Adaptation to NEMO
8   !!            3.0  ! 2007-07  (D. Storkey) Changes to iom_gettime
9   !!            3.4  ! 2012-12  (R. Bourdalle-Badie and G. Reffray)  add C1D case 
10   !!            3.6  ! 2014-15  DIMG format removed
11   !!            3.6  ! 2015-15  (J. Harle) Added procedure to read REAL attributes
12   !!            4.0  ! 2017-11  (M. Andrejczuk) Extend IOM interface to write any 3D fields
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   iom_open       : open a file read only
17   !!   iom_close      : close a file or all files opened by iom
18   !!   iom_get        : read a field (interfaced to several routines)
19   !!   iom_varid      : get the id of a variable in a file
20   !!   iom_rstput     : write a field in a restart file (interfaced to several routines)
21   !!----------------------------------------------------------------------
22   USE dom_oce         ! ocean space and time domain
23   USE c1d             ! 1D vertical configuration
24   USE flo_oce         ! floats module declarations
25   USE lbclnk          ! lateal boundary condition / mpp exchanges
26   USE iom_def         ! iom variables definitions
27   USE iom_nf90        ! NetCDF format with native NetCDF library
28   USE in_out_manager  ! I/O manager
29   USE lib_mpp           ! MPP library
30#if defined key_iomput
31   USE sbc_oce  , ONLY :   nn_fsbc         ! ocean space and time domain
32   USE trc_oce  , ONLY :   nn_dttrc        !  !: frequency of step on passive tracers
33   USE icb_oce  , ONLY :   nclasses, class_num       !  !: iceberg classes
34#if defined key_si3
35   USE ice      , ONLY :   jpl
36#endif
37   USE domngb          ! ocean space and time domain
38   USE phycst          ! physical constants
39   USE dianam          ! build name of file
40   USE xios
41# endif
42   USE ioipsl, ONLY :  ju2ymds    ! for calendar
43   USE crs             ! Grid coarsening
44#if defined key_top
45   USE trc, ONLY    :  profsed
46#endif
47   USE lib_fortran 
48   USE diurnal_bulk, ONLY : ln_diurnal_only, ln_diurnal
49#if defined key_si3
50   USE ice, ONLY: nlay_i, nlay_s
51#endif
52
53
54   IMPLICIT NONE
55   PUBLIC   !   must be public to be able to access iom_def through iom
56   
57#if defined key_iomput
58   LOGICAL, PUBLIC, PARAMETER ::   lk_iomput = .TRUE.        !: iom_put flag
59#else
60   LOGICAL, PUBLIC, PARAMETER ::   lk_iomput = .FALSE.       !: iom_put flag
61#endif
62   PUBLIC iom_init, iom_swap, iom_open, iom_close, iom_setkt, iom_varid, iom_get
63   PUBLIC iom_chkatt, iom_getatt, iom_putatt, iom_getszuld, iom_rstput, iom_delay_rst, iom_put
64   PUBLIC iom_use, iom_context_finalize, iom_miss_val
65
66   PRIVATE iom_rp0d, iom_rp1d, iom_rp2d, iom_rp3d
67   PRIVATE iom_g0d, iom_g1d, iom_g2d, iom_g3d, iom_g4d, iom_get_123d
68   PRIVATE iom_p1d, iom_p2d, iom_p3d
69#if defined key_iomput
70   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
71   PRIVATE set_grid, set_grid_bounds, set_scalar, set_xmlatt, set_mooring, iom_update_file_name, iom_sdate
72   PRIVATE iom_set_rst_context, iom_set_rstw_active, iom_set_rstr_active
73# endif
74   PUBLIC iom_set_rstw_var_active, iom_set_rstw_core, iom_set_rst_vars
75   PRIVATE iom_add
76   INTEGER, PRIVATE :: numrst
77
78   INTERFACE iom_get
79      MODULE PROCEDURE iom_g0d, iom_g1d, iom_g2d, iom_g3d, iom_g4d
80   END INTERFACE
81   INTERFACE iom_getatt
82      MODULE PROCEDURE iom_g0d_iatt, iom_g1d_iatt, iom_g0d_ratt, iom_g1d_ratt, iom_g0d_catt
83   END INTERFACE
84   INTERFACE iom_putatt
85      MODULE PROCEDURE iom_p0d_iatt, iom_p1d_iatt, iom_p0d_ratt, iom_p1d_ratt, iom_p0d_catt
86   END INTERFACE
87   INTERFACE iom_rstput
88      MODULE PROCEDURE iom_rp0d, iom_rp1d, iom_rp2d, iom_rp3d
89   END INTERFACE
90   INTERFACE iom_put
91      MODULE PROCEDURE iom_p0d, iom_p1d, iom_p2d, iom_p3d
92   END INTERFACE iom_put
93 
94   !!----------------------------------------------------------------------
95   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
96   !! $Id$
97   !! Software governed by the CeCILL license (see ./LICENSE)
98   !!----------------------------------------------------------------------
99CONTAINS
100
101   SUBROUTINE iom_init( cdname, fname, ld_tmppatch ) 
102      !!----------------------------------------------------------------------
103      !!                     ***  ROUTINE   ***
104      !!
105      !! ** Purpose :   
106      !!
107      !!----------------------------------------------------------------------
108      CHARACTER(len=*),           INTENT(in)  :: cdname
109      CHARACTER(len=*), OPTIONAL, INTENT(in)  :: fname
110      LOGICAL         , OPTIONAL, INTENT(in)  :: ld_tmppatch
111#if defined key_iomput
112      !
113      TYPE(xios_duration) :: dtime    = xios_duration(0, 0, 0, 0, 0, 0)
114      TYPE(xios_date)     :: start_date
115      CHARACTER(len=lc) :: clname
116      INTEGER           :: ji, jkmin
117      LOGICAL :: llrst_context              ! is context related to restart
118      !
119      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zt_bnds, zw_bnds
120      LOGICAL ::   ll_tmppatch = .TRUE.    !: seb: patch before we remove periodicity
121      INTEGER ::   nldi_save, nlei_save    !:      and close boundaries in output files
122      INTEGER ::   nldj_save, nlej_save    !:
123      !!----------------------------------------------------------------------
124      !
125      ! seb: patch before we remove periodicity and close boundaries in output files
126      IF( PRESENT(ld_tmppatch) ) THEN   ;   ll_tmppatch = ld_tmppatch
127      ELSE                              ;   ll_tmppatch = .TRUE.
128      ENDIF
129      IF ( ll_tmppatch ) THEN
130         nldi_save = nldi   ;   nlei_save = nlei
131         nldj_save = nldj   ;   nlej_save = nlej
132         IF( nimpp           ==      1 ) nldi = 1
133         IF( nimpp + jpi - 1 == jpiglo ) nlei = jpi
134         IF( njmpp           ==      1 ) nldj = 1
135         IF( njmpp + jpj - 1 == jpjglo ) nlej = jpj
136      ENDIF
137      !
138      ALLOCATE( zt_bnds(2,jpk), zw_bnds(2,jpk) )
139      !
140      clname = cdname
141      IF( TRIM(Agrif_CFixed()) /= '0' )   clname = TRIM(Agrif_CFixed())//"_"//TRIM(cdname)
142      CALL xios_context_initialize(TRIM(clname), mpi_comm_oce)
143      CALL iom_swap( cdname )
144      llrst_context =  (TRIM(cdname) == TRIM(crxios_context) .OR. TRIM(cdname) == TRIM(cwxios_context))
145
146      ! Calendar type is now defined in xml file
147      SELECT CASE ( nleapy )        ! Choose calendar for IOIPSL
148      CASE ( 1)   ; CALL xios_define_calendar( TYPE = "Gregorian", time_origin = xios_date(1900,01,01,00,00,00), &
149          &                                    start_date = xios_date(nyear,nmonth,nday,0,0,0) )
150      CASE ( 0)   ; CALL xios_define_calendar( TYPE = "NoLeap"   , time_origin = xios_date(1900,01,01,00,00,00), &
151          &                                    start_date = xios_date(nyear,nmonth,nday,0,0,0) )
152      CASE (30)   ; CALL xios_define_calendar( TYPE = "D360"     , time_origin = xios_date(1900,01,01,00,00,00), &
153          &                                    start_date = xios_date(nyear,nmonth,nday,0,0,0) )
154      END SELECT
155
156      ! horizontal grid definition
157      IF(.NOT.llrst_context) CALL set_scalar
158      !
159      IF( TRIM(cdname) == TRIM(cxios_context) ) THEN 
160         CALL set_grid( "T", glamt, gphit, .FALSE., .FALSE. ) 
161         CALL set_grid( "U", glamu, gphiu, .FALSE., .FALSE. )
162         CALL set_grid( "V", glamv, gphiv, .FALSE., .FALSE. )
163         CALL set_grid( "W", glamt, gphit, .FALSE., .FALSE. )
164         CALL set_grid_znl( gphit )
165         !
166         IF( ln_cfmeta ) THEN   ! Add additional grid metadata
167            CALL iom_set_domain_attr("grid_T", area = e1e2t(nldi:nlei, nldj:nlej))
168            CALL iom_set_domain_attr("grid_U", area = e1e2u(nldi:nlei, nldj:nlej))
169            CALL iom_set_domain_attr("grid_V", area = e1e2v(nldi:nlei, nldj:nlej))
170            CALL iom_set_domain_attr("grid_W", area = e1e2t(nldi:nlei, nldj:nlej))
171            CALL set_grid_bounds( "T", glamf, gphif, glamt, gphit )
172            CALL set_grid_bounds( "U", glamv, gphiv, glamu, gphiu )
173            CALL set_grid_bounds( "V", glamu, gphiu, glamv, gphiv )
174            CALL set_grid_bounds( "W", glamf, gphif, glamt, gphit )
175         ENDIF
176      ENDIF
177      !
178      IF( TRIM(cdname) == TRIM(cxios_context)//"_crs" ) THEN 
179         CALL dom_grid_crs   ! Save the parent grid information  & Switch to coarse grid domain
180         !
181         CALL set_grid( "T", glamt_crs, gphit_crs, .FALSE., .FALSE. ) 
182         CALL set_grid( "U", glamu_crs, gphiu_crs, .FALSE., .FALSE. ) 
183         CALL set_grid( "V", glamv_crs, gphiv_crs, .FALSE., .FALSE. ) 
184         CALL set_grid( "W", glamt_crs, gphit_crs, .FALSE., .FALSE. ) 
185         CALL set_grid_znl( gphit_crs )
186          !
187         CALL dom_grid_glo   ! Return to parent grid domain
188         !
189         IF( ln_cfmeta .AND. .NOT. llrst_context) THEN   ! Add additional grid metadata
190            CALL iom_set_domain_attr("grid_T", area = e1e2t_crs(nldi:nlei, nldj:nlej))
191            CALL iom_set_domain_attr("grid_U", area = e1u_crs(nldi:nlei, nldj:nlej) * e2u_crs(nldi:nlei, nldj:nlej))
192            CALL iom_set_domain_attr("grid_V", area = e1v_crs(nldi:nlei, nldj:nlej) * e2v_crs(nldi:nlei, nldj:nlej))
193            CALL iom_set_domain_attr("grid_W", area = e1e2t_crs(nldi:nlei, nldj:nlej))
194            CALL set_grid_bounds( "T", glamf_crs, gphif_crs, glamt_crs, gphit_crs )
195            CALL set_grid_bounds( "U", glamv_crs, gphiv_crs, glamu_crs, gphiu_crs )
196            CALL set_grid_bounds( "V", glamu_crs, gphiu_crs, glamv_crs, gphiv_crs )
197            CALL set_grid_bounds( "W", glamf_crs, gphif_crs, glamt_crs, gphit_crs )
198         ENDIF
199      ENDIF
200      !
201      ! vertical grid definition
202      IF(.NOT.llrst_context) THEN
203          CALL iom_set_axis_attr( "deptht",  paxis = gdept_1d )
204          CALL iom_set_axis_attr( "depthu",  paxis = gdept_1d )
205          CALL iom_set_axis_attr( "depthv",  paxis = gdept_1d )
206          CALL iom_set_axis_attr( "depthw",  paxis = gdepw_1d )
207
208          ! Add vertical grid bounds
209          jkmin = MIN(2,jpk)  ! in case jpk=1 (i.e. sas2D)
210          zt_bnds(2,:        ) = gdept_1d(:)
211          zt_bnds(1,jkmin:jpk) = gdept_1d(1:jpkm1)
212          zt_bnds(1,1        ) = gdept_1d(1) - e3w_1d(1)
213          zw_bnds(1,:        ) = gdepw_1d(:)
214          zw_bnds(2,1:jpkm1  ) = gdepw_1d(jkmin:jpk)
215          zw_bnds(2,jpk:     ) = gdepw_1d(jpk) + e3t_1d(jpk)
216          CALL iom_set_axis_attr( "deptht", bounds=zw_bnds )
217          CALL iom_set_axis_attr( "depthu", bounds=zw_bnds )
218          CALL iom_set_axis_attr( "depthv", bounds=zw_bnds )
219          CALL iom_set_axis_attr( "depthw", bounds=zt_bnds )
220          CALL iom_set_axis_attr( "nfloat", (/ (REAL(ji,wp), ji=1,jpnfl) /) )
221# if defined key_si3
222          CALL iom_set_axis_attr( "ncatice", (/ (REAL(ji,wp), ji=1,jpl) /) )
223          ! SIMIP diagnostics (4 main arctic straits)
224          CALL iom_set_axis_attr( "nstrait", (/ (REAL(ji,wp), ji=1,4) /) )
225# endif
226#if defined key_top
227          IF( ALLOCATED(profsed) ) CALL iom_set_axis_attr( "profsed", paxis = profsed )
228#endif
229          CALL iom_set_axis_attr( "icbcla", class_num )
230          CALL iom_set_axis_attr( "iax_20C", (/ REAL(20,wp) /) )   ! strange syntaxe and idea...
231          CALL iom_set_axis_attr( "iax_28C", (/ REAL(28,wp) /) )   ! strange syntaxe and idea...
232      ENDIF
233      !
234      ! automatic definitions of some of the xml attributs
235      IF( TRIM(cdname) == TRIM(crxios_context) .OR. TRIM(cdname) == TRIM(crixios_context) ) THEN
236          numrst = numror
237          IF(TRIM(cdname) == TRIM(crixios_context)) numrst = numrir 
238!set names of the fields in restart file IF using XIOS to read data
239          CALL iom_set_rst_context(.TRUE.)
240          CALL iom_set_rst_vars(rst_rfields)
241!set which fields are to be read from restart file
242          CALL iom_set_rstr_active(fname)
243      ELSE IF( TRIM(cdname) == TRIM(cwxios_context) ) THEN
244!set names of the fields in restart file IF using XIOS to write data
245          CALL iom_set_rst_context(.FALSE.)
246          CALL iom_set_rst_vars(rst_wfields)
247!set which fields are to be written to a restart file
248          CALL iom_set_rstw_active(rst_wfields, fname)
249      ELSE
250          CALL set_xmlatt
251      ENDIF
252      !
253      ! end file definition
254      dtime%second = rdt
255      CALL xios_set_timestep( dtime )
256      CALL xios_close_context_definition()
257      CALL xios_update_calendar( 0 )
258      !
259      DEALLOCATE( zt_bnds, zw_bnds )
260      !
261      IF ( ll_tmppatch ) THEN
262         nldi = nldi_save   ;   nlei = nlei_save
263         nldj = nldj_save   ;   nlej = nlej_save
264      ENDIF
265#endif
266      !
267   END SUBROUTINE iom_init
268
269   SUBROUTINE iom_set_rstw_var_active(drst_wfields, field)
270      !!---------------------------------------------------------------------
271      !!                   ***  SUBROUTINE  iom_set_rstw_var_active  ***
272      !!
273      !! ** Purpose :  enable variable in restart file when writing with XIOS
274      !!---------------------------------------------------------------------
275   CHARACTER(len = *), INTENT(IN) :: field
276   TYPE(RST_FIELD),    INTENT(INOUT) :: drst_wfields(max_rst_fields)
277   INTEGER :: i
278   LOGICAL :: llis_set
279   CHARACTER(LEN=256) :: clinfo    ! info character
280
281#if defined key_iomput
282   llis_set = .FALSE.
283
284   DO i = 1, max_rst_fields
285       IF(TRIM(drst_wfields(i)%vname) == field) THEN
286          drst_wfields(i)%active = .TRUE.
287          llis_set = .TRUE.
288          EXIT
289       ENDIF
290   ENDDO
291!Warn if variable is not in defined in drst_wfields
292   IF(.NOT.llis_set) THEN
293      WRITE(ctmp1,*) 'iom_set_rstw_var_active: variable ', field ,' is available for writing but not defined' 
294      CALL ctl_stop( 'iom_set_rstw_var_active:', ctmp1 )
295   ENDIF
296#else
297        clinfo = 'iom_set_rstw_var_active: key_iomput is needed to use XIOS restart read/write functionality'
298        CALL ctl_stop('STOP', TRIM(clinfo))
299#endif
300
301   END SUBROUTINE iom_set_rstw_var_active
302
303   SUBROUTINE iom_set_rstr_active(rst_file)
304      !!---------------------------------------------------------------------
305      !!                   ***  SUBROUTINE  iom_set_rstr_active  ***
306      !!
307      !! ** Purpose :  define file name in XIOS context for reading restart file,
308      !!               enable variables present in restart file for reading with XIOS
309      !!---------------------------------------------------------------------
310
311!sets enabled = .TRUE. for each field in restart file
312   CHARACTER(len=256) :: rst_file
313
314#if defined key_iomput
315   TYPE(xios_field) :: field_hdl
316   TYPE(xios_file) :: file_hdl
317   TYPE(xios_filegroup) :: filegroup_hdl
318   INTEGER :: i
319
320!set name of the restart file and enable available fields
321        if(lwp) WRITE(numout,*) 'Setting restart filename (for XIOS) to: ',rst_file
322        CALL xios_get_handle("file_definition", filegroup_hdl )
323        CALL xios_add_child(filegroup_hdl, file_hdl, 'rrestart')
324        CALL xios_set_file_attr( "rrestart", name=trim(rst_file), type="one_file", &
325             par_access="collective", enabled=.TRUE., mode="read",                 &
326             output_freq=xios_timestep)
327!define variables for restart context
328        DO i = 1, max_rst_fields
329         IF( TRIM(rst_rfields(i)%vname) /= "NO_NAME") THEN
330           IF( iom_varid( numrst, TRIM(rst_rfields(i)%vname), ldstop = .FALSE. ) > 0 ) THEN
331                IF(TRIM(rst_rfields(i)%vname) == "e_s") lrsi3_nemo = .FALSE.       ! SI3 restart written by XIOS
332                CALL xios_add_child(file_hdl, field_hdl, TRIM(rst_rfields(i)%vname))
333                SELECT CASE (TRIM(rst_rfields(i)%grid))
334!                CASE ("grid_Ni_4D")
335!                   CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(rst_rfields(i)%vname), &
336!                       domain_ref="grid_N", axis_ref="nlay_i", axis_ref="numcat", operation = "instant")
337!                CASE ("grid_Ns_4D")
338!                   CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(rst_rfields(i)%vname), &
339!                       domain_ref="grid_N", axis_ref="nlay_s", axis_ref="numcat", operation = "instant")
340                 CASE ("grid_Ni_3D")
341                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(rst_rfields(i)%vname), &
342                        domain_ref="grid_N", axis_ref="numcat", operation = "instant")
343                 CASE ("grid_N_3D")
344                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(rst_rfields(i)%vname), &
345                        domain_ref="grid_N", axis_ref="nav_lev", operation = "instant")
346                 CASE ("grid_N")
347                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(rst_rfields(i)%vname), &
348                        domain_ref="grid_N", operation = "instant") 
349                CASE ("grid_vector")
350                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(rst_rfields(i)%vname), &
351                         axis_ref="nav_lev", operation = "instant")
352                 CASE ("grid_scalar")
353                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(rst_rfields(i)%vname), &
354                        scalar_ref = "grid_scalar", operation = "instant")
355                END SELECT
356                IF(lwp) WRITE(numout,*) 'XIOS read: ', TRIM(rst_rfields(i)%vname), ' enabled in ', TRIM(rst_file)
357           ENDIF
358         ENDIF
359        END DO
360#endif
361   END SUBROUTINE iom_set_rstr_active
362
363   SUBROUTINE iom_set_rstw_core(cdmdl)
364      !!---------------------------------------------------------------------
365      !!                   ***  SUBROUTINE  iom_set_rstw_core  ***
366      !!
367      !! ** Purpose :  set variables which are always in restart file
368      !!---------------------------------------------------------------------
369   CHARACTER (len=*), INTENT (IN) :: cdmdl ! model OPA or SAS
370   CHARACTER(LEN=256)             :: clinfo    ! info character
371#if defined key_iomput
372   IF(cdmdl == "OPA") THEN
373!from restart.F90
374   CALL iom_set_rstw_var_active(rst_wfields, "rdt")
375   IF ( .NOT. ln_diurnal_only ) THEN
376        CALL iom_set_rstw_var_active(rst_wfields, 'ub'  )
377        CALL iom_set_rstw_var_active(rst_wfields, 'vb'  )
378        CALL iom_set_rstw_var_active(rst_wfields, 'tb'  )
379        CALL iom_set_rstw_var_active(rst_wfields, 'sb'  )
380        CALL iom_set_rstw_var_active(rst_wfields, 'sshb')
381        !
382        CALL iom_set_rstw_var_active(rst_wfields, 'un'  )
383        CALL iom_set_rstw_var_active(rst_wfields, 'vn'  )
384        CALL iom_set_rstw_var_active(rst_wfields, 'tn'  )
385        CALL iom_set_rstw_var_active(rst_wfields, 'sn'  )
386        CALL iom_set_rstw_var_active(rst_wfields, 'sshn')
387        CALL iom_set_rstw_var_active(rst_wfields, 'rhop')
388     ! extra variable needed for the ice sheet coupling
389        IF ( ln_iscpl ) THEN
390             CALL iom_set_rstw_var_active(rst_wfields, 'tmask')
391             CALL iom_set_rstw_var_active(rst_wfields, 'umask')
392             CALL iom_set_rstw_var_active(rst_wfields, 'vmask')
393             CALL iom_set_rstw_var_active(rst_wfields, 'smask')
394             CALL iom_set_rstw_var_active(rst_wfields, 'e3t_n')
395             CALL iom_set_rstw_var_active(rst_wfields, 'e3u_n')
396             CALL iom_set_rstw_var_active(rst_wfields, 'e3v_n')
397             CALL iom_set_rstw_var_active(rst_wfields, 'gdepw_n')
398        END IF
399      ENDIF
400      IF(ln_diurnal) CALL iom_set_rstw_var_active(rst_wfields, 'Dsst')
401!from trasbc.F90
402         CALL iom_set_rstw_var_active(rst_wfields, 'sbc_hc_b')
403         CALL iom_set_rstw_var_active(rst_wfields, 'sbc_sc_b')
404   ENDIF
405#else
406        clinfo = 'iom_set_rstw_core: key_iomput is needed to use XIOS restart read/write functionality'
407        CALL ctl_stop('STOP', TRIM(clinfo))
408#endif
409   END SUBROUTINE iom_set_rstw_core
410
411   SUBROUTINE iom_set_rst_vars(fields)
412      !!---------------------------------------------------------------------
413      !!                   ***  SUBROUTINE iom_set_rst_vars   ***
414      !!
415      !! ** Purpose :  Fill array fields with the information about all
416      !!               possible variables and corresponding grids definition
417      !!               for reading/writing restart with XIOS
418      !!---------------------------------------------------------------------
419   TYPE(RST_FIELD), INTENT(INOUT) :: fields(max_rst_fields)
420   INTEGER :: i, jk
421   CHARACTER(len=2)  ::   zchar
422   CHARACTER(len=25) ::   znam
423
424        i = 0
425        i = iom_add(i); fields(i)%vname="rdt";            fields(i)%grid="grid_scalar"
426        i = iom_add(i); fields(i)%vname="un";             fields(i)%grid="grid_N_3D"
427        i = iom_add(i); fields(i)%vname="ub";             fields(i)%grid="grid_N_3D"
428        i = iom_add(i); fields(i)%vname="vn";             fields(i)%grid="grid_N_3D"
429        i = iom_add(i); fields(i)%vname="vb";             fields(i)%grid="grid_N_3D" 
430        i = iom_add(i); fields(i)%vname="tn";             fields(i)%grid="grid_N_3D"
431        i = iom_add(i); fields(i)%vname="tb";             fields(i)%grid="grid_N_3D"
432        i = iom_add(i); fields(i)%vname="sn";             fields(i)%grid="grid_N_3D"
433        i = iom_add(i); fields(i)%vname="sb";             fields(i)%grid="grid_N_3D"
434        i = iom_add(i); fields(i)%vname="sshn";           fields(i)%grid="grid_N"
435        i = iom_add(i); fields(i)%vname="sshb";           fields(i)%grid="grid_N"
436        i = iom_add(i); fields(i)%vname="rhop";           fields(i)%grid="grid_N_3D"
437        i = iom_add(i); fields(i)%vname="kt";             fields(i)%grid="grid_scalar"
438        i = iom_add(i); fields(i)%vname="ndastp";         fields(i)%grid="grid_scalar"
439        i = iom_add(i); fields(i)%vname="adatrj";         fields(i)%grid="grid_scalar"
440        i = iom_add(i); fields(i)%vname="utau_b";         fields(i)%grid="grid_N"
441        i = iom_add(i); fields(i)%vname="vtau_b";         fields(i)%grid="grid_N"
442        i = iom_add(i); fields(i)%vname="qns_b";          fields(i)%grid="grid_N"
443        i = iom_add(i); fields(i)%vname="emp_b";          fields(i)%grid="grid_N"
444        i = iom_add(i); fields(i)%vname="sfx_b";          fields(i)%grid="grid_N"
445        i = iom_add(i); fields(i)%vname="en" ;            fields(i)%grid="grid_N_3D" 
446        i = iom_add(i); fields(i)%vname="avt_k";            fields(i)%grid="grid_N_3D"
447        i = iom_add(i); fields(i)%vname="avm_k";            fields(i)%grid="grid_N_3D"
448        i = iom_add(i); fields(i)%vname="dissl";          fields(i)%grid="grid_N_3D"
449        i = iom_add(i); fields(i)%vname="sbc_hc_b";       fields(i)%grid="grid_N"
450        i = iom_add(i); fields(i)%vname="sbc_sc_b";       fields(i)%grid="grid_N"
451        i = iom_add(i); fields(i)%vname="qsr_hc_b";       fields(i)%grid="grid_N_3D"
452        i = iom_add(i); fields(i)%vname="fraqsr_1lev";    fields(i)%grid="grid_N"
453        i = iom_add(i); fields(i)%vname="greenland_icesheet_mass"
454                                               fields(i)%grid="grid_scalar"
455        i = iom_add(i); fields(i)%vname="greenland_icesheet_timelapsed"
456                                               fields(i)%grid="grid_scalar"
457        i = iom_add(i); fields(i)%vname="greenland_icesheet_mass_roc"
458                                               fields(i)%grid="grid_scalar"
459        i = iom_add(i); fields(i)%vname="antarctica_icesheet_mass"
460                                               fields(i)%grid="grid_scalar"
461        i = iom_add(i); fields(i)%vname="antarctica_icesheet_timelapsed"
462                                               fields(i)%grid="grid_scalar"
463        i = iom_add(i); fields(i)%vname="antarctica_icesheet_mass_roc"
464                                               fields(i)%grid="grid_scalar"
465        i = iom_add(i); fields(i)%vname="frc_v";          fields(i)%grid="grid_scalar"
466        i = iom_add(i); fields(i)%vname="frc_t";          fields(i)%grid="grid_scalar"
467        i = iom_add(i); fields(i)%vname="frc_s";          fields(i)%grid="grid_scalar"
468        i = iom_add(i); fields(i)%vname="frc_wn_t";       fields(i)%grid="grid_scalar"
469        i = iom_add(i); fields(i)%vname="frc_wn_s";       fields(i)%grid="grid_scalar"
470        i = iom_add(i); fields(i)%vname="ssh_ini";        fields(i)%grid="grid_N"
471        i = iom_add(i); fields(i)%vname="e3t_ini";        fields(i)%grid="grid_N_3D"
472        i = iom_add(i); fields(i)%vname="hc_loc_ini";     fields(i)%grid="grid_N_3D"
473        i = iom_add(i); fields(i)%vname="sc_loc_ini";     fields(i)%grid="grid_N_3D"
474        i = iom_add(i); fields(i)%vname="ssh_hc_loc_ini"; fields(i)%grid="grid_N"
475        i = iom_add(i); fields(i)%vname="ssh_sc_loc_ini"; fields(i)%grid="grid_N"
476        i = iom_add(i); fields(i)%vname="tilde_e3t_b";    fields(i)%grid="grid_N"
477        i = iom_add(i); fields(i)%vname="tilde_e3t_n";    fields(i)%grid="grid_N"
478        i = iom_add(i); fields(i)%vname="hdiv_lf";        fields(i)%grid="grid_N"
479        i = iom_add(i); fields(i)%vname="ub2_b";          fields(i)%grid="grid_N"
480        i = iom_add(i); fields(i)%vname="vb2_b";          fields(i)%grid="grid_N"
481        i = iom_add(i); fields(i)%vname="sshbb_e";        fields(i)%grid="grid_N"
482        i = iom_add(i); fields(i)%vname="ubb_e";          fields(i)%grid="grid_N"
483        i = iom_add(i); fields(i)%vname="vbb_e";          fields(i)%grid="grid_N"
484        i = iom_add(i); fields(i)%vname="sshb_e";         fields(i)%grid="grid_N"
485        i = iom_add(i); fields(i)%vname="ub_e";           fields(i)%grid="grid_N"
486        i = iom_add(i); fields(i)%vname="vb_e";           fields(i)%grid="grid_N"
487        i = iom_add(i); fields(i)%vname="fwf_isf_b";      fields(i)%grid="grid_N"
488        i = iom_add(i); fields(i)%vname="isf_sc_b";       fields(i)%grid="grid_N"
489        i = iom_add(i); fields(i)%vname="isf_hc_b";       fields(i)%grid="grid_N"
490        i = iom_add(i); fields(i)%vname="ssh_ibb";        fields(i)%grid="grid_N"
491        i = iom_add(i); fields(i)%vname="rnf_b";          fields(i)%grid="grid_N"
492        i = iom_add(i); fields(i)%vname="rnf_hc_b";       fields(i)%grid="grid_N"
493        i = iom_add(i); fields(i)%vname="rnf_sc_b";       fields(i)%grid="grid_N"
494        i = iom_add(i); fields(i)%vname="nn_fsbc";        fields(i)%grid="grid_scalar"
495        i = iom_add(i); fields(i)%vname="ssu_m";          fields(i)%grid="grid_N"
496        i = iom_add(i); fields(i)%vname="ssv_m";          fields(i)%grid="grid_N"
497        i = iom_add(i); fields(i)%vname="sst_m";          fields(i)%grid="grid_N"
498        i = iom_add(i); fields(i)%vname="sss_m";          fields(i)%grid="grid_N"
499        i = iom_add(i); fields(i)%vname="ssh_m";          fields(i)%grid="grid_N"
500        i = iom_add(i); fields(i)%vname="e3t_m";          fields(i)%grid="grid_N"
501        i = iom_add(i); fields(i)%vname="frq_m";          fields(i)%grid="grid_N"
502        i = iom_add(i); fields(i)%vname="avmb";           fields(i)%grid="grid_vector"
503        i = iom_add(i); fields(i)%vname="avtb";           fields(i)%grid="grid_vector"
504        i = iom_add(i); fields(i)%vname="ub2_i_b";        fields(i)%grid="grid_N"
505        i = iom_add(i); fields(i)%vname="vb2_i_b";        fields(i)%grid="grid_N"
506        i = iom_add(i); fields(i)%vname="ntime";          fields(i)%grid="grid_scalar"
507        i = iom_add(i); fields(i)%vname="Dsst";           fields(i)%grid="grid_scalar"
508        i = iom_add(i); fields(i)%vname="tmask";          fields(i)%grid="grid_N_3D"
509        i = iom_add(i); fields(i)%vname="umask";          fields(i)%grid="grid_N_3D"
510        i = iom_add(i); fields(i)%vname="vmask";          fields(i)%grid="grid_N_3D"
511        i = iom_add(i); fields(i)%vname="smask";          fields(i)%grid="grid_N_3D"
512        i = iom_add(i); fields(i)%vname="gdepw_n";        fields(i)%grid="grid_N_3D"
513        i = iom_add(i); fields(i)%vname="e3t_n";          fields(i)%grid="grid_N_3D"
514        i = iom_add(i); fields(i)%vname="e3u_n";          fields(i)%grid="grid_N_3D"
515        i = iom_add(i); fields(i)%vname="e3v_n";          fields(i)%grid="grid_N_3D"
516        i = iom_add(i); fields(i)%vname="surf_ini";       fields(i)%grid="grid_N"
517        i = iom_add(i); fields(i)%vname="e3t_b";          fields(i)%grid="grid_N_3D"
518        i = iom_add(i); fields(i)%vname="hmxl_n";         fields(i)%grid="grid_N_3D"
519        i = iom_add(i); fields(i)%vname="un_bf";          fields(i)%grid="grid_N"
520        i = iom_add(i); fields(i)%vname="vn_bf";          fields(i)%grid="grid_N"
521        i = iom_add(i); fields(i)%vname="hbl";            fields(i)%grid="grid_N"
522        i = iom_add(i); fields(i)%vname="hbli";           fields(i)%grid="grid_N"
523        i = iom_add(i); fields(i)%vname="wn";             fields(i)%grid="grid_N_3D"
524#ifdef key_si3
525!ICE/icedia.F90
526        i = iom_add(i); fields(i)%vname="kt_ice";         fields(i)%grid="grid_scalar"
527        i = iom_add(i); fields(i)%vname="frc_voltop";     fields(i)%grid="grid_scalar"
528        i = iom_add(i); fields(i)%vname="frc_volbot";     fields(i)%grid="grid_scalar"
529        i = iom_add(i); fields(i)%vname="frc_temtop";     fields(i)%grid="grid_scalar"
530        i = iom_add(i); fields(i)%vname="frc_tembot";     fields(i)%grid="grid_scalar" 
531        i = iom_add(i); fields(i)%vname="frc_sal";        fields(i)%grid="grid_scalar"
532        i = iom_add(i); fields(i)%vname="vol_loc_ini";    fields(i)%grid="grid_N"
533        i = iom_add(i); fields(i)%vname="tem_loc_ini";    fields(i)%grid="grid_N"
534        i = iom_add(i); fields(i)%vname="sal_loc_ini";    fields(i)%grid="grid_N"
535        i = iom_add(i); fields(i)%vname="sal_loc_ini";    fields(i)%grid="grid_N"
536!ICE/icedyn_adv_pra.F90 grid_Ni_3D - jpi,jpj,jpl
537        i = iom_add(i); fields(i)%vname="sxice";          fields(i)%grid="grid_Ni_3D"
538        i = iom_add(i); fields(i)%vname="syice";          fields(i)%grid="grid_Ni_3D"
539        i = iom_add(i); fields(i)%vname="sxxice";         fields(i)%grid="grid_Ni_3D"
540        i = iom_add(i); fields(i)%vname="syyice";         fields(i)%grid="grid_Ni_3D"
541        i = iom_add(i); fields(i)%vname="sxyice";         fields(i)%grid="grid_Ni_3D"
542        i = iom_add(i); fields(i)%vname="sxsn";           fields(i)%grid="grid_Ni_3D"
543        i = iom_add(i); fields(i)%vname="sysn";           fields(i)%grid="grid_Ni_3D"
544        i = iom_add(i); fields(i)%vname="sxxsn";          fields(i)%grid="grid_Ni_3D"
545        i = iom_add(i); fields(i)%vname="syysn";          fields(i)%grid="grid_Ni_3D"
546        i = iom_add(i); fields(i)%vname="sxysn";          fields(i)%grid="grid_Ni_3D"
547        i = iom_add(i); fields(i)%vname="sxa";            fields(i)%grid="grid_Ni_3D"
548        i = iom_add(i); fields(i)%vname="sya";            fields(i)%grid="grid_Ni_3D"
549        i = iom_add(i); fields(i)%vname="sxxa";           fields(i)%grid="grid_Ni_3D"
550        i = iom_add(i); fields(i)%vname="syya";           fields(i)%grid="grid_Ni_3D"
551        i = iom_add(i); fields(i)%vname="sxya";           fields(i)%grid="grid_Ni_3D"
552        i = iom_add(i); fields(i)%vname="sxsal";          fields(i)%grid="grid_Ni_3D"
553        i = iom_add(i); fields(i)%vname="sysal";          fields(i)%grid="grid_Ni_3D"
554        i = iom_add(i); fields(i)%vname="sxxsal";         fields(i)%grid="grid_Ni_3D"
555        i = iom_add(i); fields(i)%vname="syysal";         fields(i)%grid="grid_Ni_3D"
556        i = iom_add(i); fields(i)%vname="sxysal";         fields(i)%grid="grid_Ni_3D"
557        i = iom_add(i); fields(i)%vname="sxage";          fields(i)%grid="grid_Ni_3D"
558        i = iom_add(i); fields(i)%vname="syage";          fields(i)%grid="grid_Ni_3D"
559        i = iom_add(i); fields(i)%vname="sxxage";         fields(i)%grid="grid_Ni_3D"
560        i = iom_add(i); fields(i)%vname="syyage";         fields(i)%grid="grid_Ni_3D"
561        i = iom_add(i); fields(i)%vname="sxyage";         fields(i)%grid="grid_Ni_3D"
562!this part is in addition to 4D grids, to handle old format of SI3 resart
563!NEMO format restart
564        DO jk = 1, nlay_s
565           WRITE(zchar,'(I2.2)') jk
566           znam = 'sxc0'//'_l'//zchar
567           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
568           znam = 'syc0'//'_l'//zchar
569           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
570           znam = 'sxxc0'//'_l'//zchar
571           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
572           znam = 'syyc0'//'_l'//zchar
573           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
574           znam = 'sxyc0'//'_l'//zchar
575           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
576        END DO
577!XIOS format restart
578        i = iom_add(i); fields(i)%vname="sxc0";           fields(i)%grid="grid_Ns_4D"
579        i = iom_add(i); fields(i)%vname="syc0";           fields(i)%grid="grid_Ns_4D"
580        i = iom_add(i); fields(i)%vname="sxxc0";          fields(i)%grid="grid_Ns_4D"
581        i = iom_add(i); fields(i)%vname="syyc0";          fields(i)%grid="grid_Ns_4D"
582        i = iom_add(i); fields(i)%vname="sxxc0";          fields(i)%grid="grid_Ns_4D"
583        !
584!NEMO format restart
585        DO jk = 1, nlay_i
586           WRITE(zchar,'(I2.2)') jk
587           znam = 'sxe'//'_l'//zchar
588           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
589           znam = 'sye'//'_l'//zchar 
590           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
591           znam = 'sxxe'//'_l'//zchar
592           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
593           znam = 'syye'//'_l'//zchar
594           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
595           znam = 'sxye'//'_l'//zchar
596           i = iom_add(i); fields(i)%vname=znam;          fields(i)%grid="grid_Ni_3D"
597        ENDDO
598!XIOS format restart
599        i = iom_add(i); fields(i)%vname="sxe";           fields(i)%grid="grid_Ni_4D"
600        i = iom_add(i); fields(i)%vname="sye";           fields(i)%grid="grid_Ni_4D"
601        i = iom_add(i); fields(i)%vname="sxxe";          fields(i)%grid="grid_Ni_4D"
602        i = iom_add(i); fields(i)%vname="syye";          fields(i)%grid="grid_Ni_4D"
603        i = iom_add(i); fields(i)%vname="sxye";          fields(i)%grid="grid_Ni_4D"
604!common for both
605        i = iom_add(i); fields(i)%vname="sxap";          fields(i)%grid="grid_Ni_3D"
606        i = iom_add(i); fields(i)%vname="syap";          fields(i)%grid="grid_Ni_3D"
607        i = iom_add(i); fields(i)%vname="sxxap";         fields(i)%grid="grid_Ni_3D"
608        i = iom_add(i); fields(i)%vname="syyap";         fields(i)%grid="grid_Ni_3D"
609        i = iom_add(i); fields(i)%vname="sxyap";         fields(i)%grid="grid_Ni_3D"
610        i = iom_add(i); fields(i)%vname="sxvp";          fields(i)%grid="grid_Ni_3D"
611        i = iom_add(i); fields(i)%vname="syvp";          fields(i)%grid="grid_Ni_3D"
612        i = iom_add(i); fields(i)%vname="sxxvp";         fields(i)%grid="grid_Ni_3D"
613        i = iom_add(i); fields(i)%vname="syyvp";         fields(i)%grid="grid_Ni_3D"
614        i = iom_add(i); fields(i)%vname="sxyvp";         fields(i)%grid="grid_Ni_3D"
615        i = iom_add(i); fields(i)%vname="stress1_i";     fields(i)%grid="grid_N"
616        i = iom_add(i); fields(i)%vname="stress2_i";     fields(i)%grid="grid_N"
617        i = iom_add(i); fields(i)%vname="stress12_i";    fields(i)%grid="grid_N"
618        i = iom_add(i); fields(i)%vname="nn_fsbc";       fields(i)%grid="grid_scalar"
619        i = iom_add(i); fields(i)%vname="kt_ice";        fields(i)%grid="grid_scalar"
620        i = iom_add(i); fields(i)%vname="v_i";           fields(i)%grid="grid_Ni_3D"
621        i = iom_add(i); fields(i)%vname="v_s";           fields(i)%grid="grid_Ni_3D"
622        i = iom_add(i); fields(i)%vname="sv_i";          fields(i)%grid="grid_Ni_3D"
623        i = iom_add(i); fields(i)%vname="a_i";           fields(i)%grid="grid_Ni_3D"
624        i = iom_add(i); fields(i)%vname="t_su";          fields(i)%grid="grid_Ni_3D"
625        i = iom_add(i); fields(i)%vname="u_ice";         fields(i)%grid="grid_N"
626        i = iom_add(i); fields(i)%vname="v_ice";         fields(i)%grid="grid_N"
627!NEMO format restat
628        DO jk = 1, nlay_s
629           WRITE(zchar,'(I2.2)') jk
630           znam = 'e_s'//'_l'//zchar
631           i = iom_add(i); fields(i)%vname=znam;         fields(i)%grid="grid_Ni_3D"
632        END DO
633!XIOS format restart
634        i = iom_add(i); fields(i)%vname="e_s";           fields(i)%grid="grid_Ns_4D"
635!NEMO format restat
636        DO jk = 1, nlay_i
637           WRITE(zchar,'(I2.2)') jk
638           znam = 'e_i'//'_l'//zchar
639           i = iom_add(i); fields(i)%vname=znam;         fields(i)%grid="grid_Ni_3D"
640        END DO
641!XIOS format restart
642        i = iom_add(i); fields(i)%vname="e_i";           fields(i)%grid="grid_Ni_4D"
643!common for both
644        i = iom_add(i); fields(i)%vname="oa_i";          fields(i)%grid="grid_Ni_3D"
645        i = iom_add(i); fields(i)%vname="a_ip";          fields(i)%grid="grid_Ni_3D"
646        i = iom_add(i); fields(i)%vname="v_ip";          fields(i)%grid="grid_Ni_3D"
647        i = iom_add(i); fields(i)%vname="cnd_ice";       fields(i)%grid="grid_Ni_3D"
648        i = iom_add(i); fields(i)%vname="t1_ice";        fields(i)%grid="grid_Ni_3D"
649        i = iom_add(i); fields(i)%vname="snwice_mass";   fields(i)%grid="grid_N"
650        i = iom_add(i); fields(i)%vname="snwice_mass_b"; fields(i)%grid="grid_N"
651#endif
652   END SUBROUTINE iom_set_rst_vars
653
654
655   SUBROUTINE iom_set_rstw_active(drst_wfields, cdrst_file)
656      !!---------------------------------------------------------------------
657      !!                   ***  SUBROUTINE iom_set_rstw_active   ***
658      !!
659      !! ** Purpose :  define file name in XIOS context for writing restart
660      !!               enable variables present in restart file for writing
661      !!---------------------------------------------------------------------
662!sets enabled = .TRUE. for each field in restart file
663   CHARACTER(len=*), INTENT(IN) :: cdrst_file
664   TYPE(RST_FIELD),  INTENT(IN) :: drst_wfields(max_rst_fields)
665#if defined key_iomput
666   TYPE(xios_field) :: field_hdl
667   TYPE(xios_file) :: file_hdl
668   TYPE(xios_filegroup) :: filegroup_hdl
669   INTEGER :: i
670   CHARACTER(lc)  ::   clpath
671
672!set name of the restart file and enable available fields
673        IF(lwp) WRITE(numout,*) 'Setting restart filename (for XIOS write) to: ',cdrst_file
674        CALL xios_get_handle("file_definition", filegroup_hdl )
675        CALL xios_add_child(filegroup_hdl, file_hdl, 'wrestart')
676        IF(nxioso.eq.1) THEN
677           CALL xios_set_file_attr( "wrestart", type="one_file", enabled=.TRUE.,& 
678                                    mode="write", output_freq=xios_timestep) 
679           if(lwp) write(numout,*) 'OPEN ', trim(cdrst_file), ' in one_file mode' 
680        ELSE 
681           CALL xios_set_file_attr( "wrestart", type="multiple_file", enabled=.TRUE.,& 
682                                    mode="write", output_freq=xios_timestep) 
683           if(lwp) write(numout,*) 'OPEN ', trim(cdrst_file), ' in multiple_file mode' 
684        ENDIF
685        CALL xios_set_file_attr( "wrestart", name=trim(cdrst_file))
686!define fields for restart context
687        DO i = 1, max_rst_fields
688         IF( drst_wfields(i)%active ) THEN
689                CALL xios_add_child(file_hdl, field_hdl, TRIM(drst_wfields(i)%vname))
690                SELECT CASE (TRIM(drst_wfields(i)%grid))
691                 CASE ("grid_N_3D")
692                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(drst_wfields(i)%vname), &
693                        domain_ref="grid_N", axis_ref="nav_lev", prec = 8, operation = "instant")
694                 CASE ("grid_N")
695                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(drst_wfields(i)%vname), &
696                        domain_ref="grid_N", prec = 8, operation = "instant") 
697                 CASE ("grid_vector")
698                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(drst_wfields(i)%vname), &
699                         axis_ref="nav_lev", prec = 8, operation = "instant")
700                 CASE ("grid_scalar")
701                    CALL xios_set_attr (field_hdl, enabled = .TRUE., name = TRIM(drst_wfields(i)%vname), &
702                        scalar_ref = "grid_scalar", prec = 8, operation = "instant")
703                END SELECT
704         ENDIF
705        END DO
706#endif
707   END SUBROUTINE iom_set_rstw_active
708
709   SUBROUTINE iom_set_rst_context(ld_rstr) 
710     !!---------------------------------------------------------------------
711      !!                   ***  SUBROUTINE  iom_set_rst_context  ***
712      !!
713      !! ** Purpose : Define domain, axis and grid for restart (read/write)
714      !!              context
715      !!               
716      !!---------------------------------------------------------------------
717   LOGICAL, INTENT(IN)               :: ld_rstr
718!ld_rstr is true for restart context. There is no need to define grid for
719!restart read, because it's read from file
720#if defined key_iomput
721   TYPE(xios_domaingroup)            :: domaingroup_hdl 
722   TYPE(xios_domain)                 :: domain_hdl 
723   TYPE(xios_axisgroup)              :: axisgroup_hdl 
724   TYPE(xios_axis)                   :: axis_hdl 
725   TYPE(xios_scalar)                 :: scalar_hdl 
726   TYPE(xios_scalargroup)            :: scalargroup_hdl 
727   INTEGER                           :: ji
728
729     CALL xios_get_handle("domain_definition",domaingroup_hdl) 
730     CALL xios_add_child(domaingroup_hdl, domain_hdl, "grid_N") 
731     CALL set_grid("N", glamt, gphit, .TRUE., ld_rstr) 
732 
733     CALL xios_get_handle("axis_definition",axisgroup_hdl) 
734     CALL xios_add_child(axisgroup_hdl, axis_hdl, "nav_lev") 
735!AGRIF fails to compile when unit= is in call to xios_set_axis_attr
736!    CALL xios_set_axis_attr( "nav_lev", long_name="Vertical levels",  unit="m", positive="down")
737     CALL xios_set_axis_attr( "nav_lev", long_name="Vertical levels in meters", positive="down")
738     CALL iom_set_axis_attr( "nav_lev", paxis = gdept_1d ) 
739
740     CALL xios_get_handle("scalar_definition", scalargroup_hdl) 
741     CALL xios_add_child(scalargroup_hdl, scalar_hdl, "grid_scalar") 
742#if defined key_si3
743     CALL xios_get_handle("axis_definition",axisgroup_hdl) 
744     CALL xios_add_child(axisgroup_hdl, axis_hdl, "numcat")
745     CALL iom_set_axis_attr( "numcat", (/ (REAL(ji,wp), ji=1,jpl) /) )
746     CALL xios_add_child(axisgroup_hdl, axis_hdl, "nlay_i")
747     CALL iom_set_axis_attr( "nlay_i", (/ (REAL(ji,wp), ji=1,nlay_i) /) )
748     CALL xios_add_child(axisgroup_hdl, axis_hdl, "nlay_s")
749     CALL iom_set_axis_attr( "nlay_s", (/ (REAL(ji,wp), ji=1,nlay_s) /) )
750#endif
751#endif
752   END SUBROUTINE iom_set_rst_context
753
754   SUBROUTINE iom_swap( cdname )
755      !!---------------------------------------------------------------------
756      !!                   ***  SUBROUTINE  iom_swap  ***
757      !!
758      !! ** Purpose :  swap context between different agrif grid for xmlio_server
759      !!---------------------------------------------------------------------
760      CHARACTER(len=*), INTENT(in) :: cdname
761#if defined key_iomput
762      TYPE(xios_context) :: nemo_hdl
763     
764      IF( TRIM(Agrif_CFixed()) == '0' ) THEN
765        CALL xios_get_handle(TRIM(cdname),nemo_hdl)
766      ELSE
767        CALL xios_get_handle(TRIM(Agrif_CFixed())//"_"//TRIM(cdname),nemo_hdl)
768      ENDIF
769      !
770      CALL xios_set_current_context(nemo_hdl)
771#endif
772      !
773   END SUBROUTINE iom_swap
774
775
776   SUBROUTINE iom_open( cdname, kiomid, ldwrt, kdom, ldstop, ldiof, kdlev )
777      !!---------------------------------------------------------------------
778      !!                   ***  SUBROUTINE  iom_open  ***
779      !!
780      !! ** Purpose :  open an input file (return 0 if not found)
781      !!---------------------------------------------------------------------
782      CHARACTER(len=*), INTENT(in   )           ::   cdname   ! File name
783      INTEGER         , INTENT(  out)           ::   kiomid   ! iom identifier of the opened file
784      LOGICAL         , INTENT(in   ), OPTIONAL ::   ldwrt    ! open in write modeb          (default = .FALSE.)
785      INTEGER         , INTENT(in   ), OPTIONAL ::   kdom     ! Type of domain to be written (default = jpdom_local_noovlap)
786      LOGICAL         , INTENT(in   ), OPTIONAL ::   ldstop   ! stop if open to read a non-existing file (default = .TRUE.)
787      LOGICAL         , INTENT(in   ), OPTIONAL ::   ldiof    ! Interp On the Fly, needed for AGRIF (default = .FALSE.)
788      INTEGER         , INTENT(in   ), OPTIONAL ::   kdlev    ! number of vertical levels
789      !
790      CHARACTER(LEN=256)    ::   clname    ! the name of the file based on cdname [[+clcpu]+clcpu]
791      CHARACTER(LEN=256)    ::   cltmpn    ! tempory name to store clname (in writting mode)
792      CHARACTER(LEN=10)     ::   clsuffix  ! ".nc"
793      CHARACTER(LEN=15)     ::   clcpu     ! the cpu number (max jpmax_digits digits)
794      CHARACTER(LEN=256)    ::   clinfo    ! info character
795      LOGICAL               ::   llok      ! check the existence
796      LOGICAL               ::   llwrt     ! local definition of ldwrt
797      LOGICAL               ::   llnoov    ! local definition to read overlap
798      LOGICAL               ::   llstop    ! local definition of ldstop
799      LOGICAL               ::   lliof     ! local definition of ldiof
800      INTEGER               ::   icnt      ! counter for digits in clcpu (max = jpmax_digits)
801      INTEGER               ::   iln, ils  ! lengths of character
802      INTEGER               ::   idom      ! type of domain
803      INTEGER               ::   istop     !
804      INTEGER, DIMENSION(2,5) ::   idompar ! domain parameters:
805      ! local number of points for x,y dimensions
806      ! position of first local point for x,y dimensions
807      ! position of last local point for x,y dimensions
808      ! start halo size for x,y dimensions
809      ! end halo size for x,y dimensions
810      !---------------------------------------------------------------------
811      ! Initializations and control
812      ! =============
813      kiomid = -1
814      clinfo = '                    iom_open ~~~  '
815      istop = nstop
816      ! if iom_open is called for the first time: initialize iom_file(:)%nfid to 0
817      ! (could be done when defining iom_file in f95 but not in f90)
818      IF( Agrif_Root() ) THEN
819         IF( iom_open_init == 0 ) THEN
820            iom_file(:)%nfid = 0
821            iom_open_init = 1
822         ENDIF
823      ENDIF
824      ! do we read or write the file?
825      IF( PRESENT(ldwrt) ) THEN   ;   llwrt = ldwrt
826      ELSE                        ;   llwrt = .FALSE.
827      ENDIF
828      ! do we call ctl_stop if we try to open a non-existing file in read mode?
829      IF( PRESENT(ldstop) ) THEN   ;   llstop = ldstop
830      ELSE                         ;   llstop = .TRUE.
831      ENDIF
832      ! are we using interpolation on the fly?
833      IF( PRESENT(ldiof) ) THEN   ;   lliof = ldiof
834      ELSE                        ;   lliof = .FALSE.
835      ENDIF
836      ! do we read the overlap
837      ! ugly patch SM+JMM+RB to overwrite global definition in some cases
838      llnoov = (jpni * jpnj ) == jpnij .AND. .NOT. lk_agrif
839      ! create the file name by added, if needed, TRIM(Agrif_CFixed()) and TRIM(clsuffix)
840      ! =============
841      clname   = trim(cdname)
842      IF ( .NOT. Agrif_Root() .AND. .NOT. lliof ) THEN
843!FUS         iln    = INDEX(clname,'/')
844         iln    = INDEX(clname,'/',BACK=.true.)  ! FUS: to insert the nest index at the right location within the string, the last / has to be found (search from the right to left)
845         cltmpn = clname(1:iln)
846         clname = clname(iln+1:LEN_TRIM(clname))
847         clname=TRIM(cltmpn)//TRIM(Agrif_CFixed())//'_'//TRIM(clname)
848      ENDIF
849      ! which suffix should we use?
850      clsuffix = '.nc'
851      ! Add the suffix if needed
852      iln = LEN_TRIM(clname)
853      ils = LEN_TRIM(clsuffix)
854      IF( iln <= ils .OR. INDEX( TRIM(clname), TRIM(clsuffix), back = .TRUE. ) /= iln - ils + 1 )   &
855         &   clname = TRIM(clname)//TRIM(clsuffix)
856      cltmpn = clname   ! store this name
857      ! try to find if the file to be opened already exist
858      ! =============
859      INQUIRE( FILE = clname, EXIST = llok )
860      IF( .NOT.llok ) THEN
861         ! we try to add the cpu number to the name
862         WRITE(clcpu,*) narea-1
863
864         clcpu  = TRIM(ADJUSTL(clcpu))
865         iln = INDEX(clname,TRIM(clsuffix), back = .TRUE.)
866         clname = clname(1:iln-1)//'_'//TRIM(clcpu)//TRIM(clsuffix)
867         icnt = 0
868         INQUIRE( FILE = clname, EXIST = llok ) 
869         ! we try different formats for the cpu number by adding 0
870         DO WHILE( .NOT.llok .AND. icnt < jpmax_digits )
871            clcpu  = "0"//trim(clcpu)
872            clname = clname(1:iln-1)//'_'//TRIM(clcpu)//TRIM(clsuffix)
873            INQUIRE( FILE = clname, EXIST = llok )
874            icnt = icnt + 1
875         END DO
876      ELSE
877         lxios_sini = .TRUE.
878      ENDIF
879      IF( llwrt ) THEN
880         ! check the domain definition
881! JMM + SM: ugly patch before getting the new version of lib_mpp)
882!         idom = jpdom_local_noovlap   ! default definition
883         IF( llnoov ) THEN   ;   idom = jpdom_local_noovlap   ! default definition
884         ELSE                ;   idom = jpdom_local_full      ! default definition
885         ENDIF
886         IF( PRESENT(kdom) )   idom = kdom
887         ! create the domain informations
888         ! =============
889         SELECT CASE (idom)
890         CASE (jpdom_local_full)
891            idompar(:,1) = (/ jpi             , jpj              /)
892            idompar(:,2) = (/ nimpp           , njmpp            /)
893            idompar(:,3) = (/ nimpp + jpi - 1 , njmpp + jpj - 1  /)
894            idompar(:,4) = (/ nldi - 1        , nldj - 1         /)
895            idompar(:,5) = (/ jpi - nlei      , jpj - nlej       /)
896         CASE (jpdom_local_noextra)
897            idompar(:,1) = (/ nlci            , nlcj             /)
898            idompar(:,2) = (/ nimpp           , njmpp            /)
899            idompar(:,3) = (/ nimpp + nlci - 1, njmpp + nlcj - 1 /)
900            idompar(:,4) = (/ nldi - 1        , nldj - 1         /)
901            idompar(:,5) = (/ nlci - nlei     , nlcj - nlej      /)
902         CASE (jpdom_local_noovlap)
903            idompar(:,1) = (/ nlei  - nldi + 1, nlej  - nldj + 1 /)
904            idompar(:,2) = (/ nimpp + nldi - 1, njmpp + nldj - 1 /)
905            idompar(:,3) = (/ nimpp + nlei - 1, njmpp + nlej - 1 /)
906            idompar(:,4) = (/ 0               , 0                /)
907            idompar(:,5) = (/ 0               , 0                /)
908         CASE DEFAULT
909            CALL ctl_stop( TRIM(clinfo), 'wrong value of kdom, only jpdom_local* cases are accepted' )
910         END SELECT
911      ENDIF
912      ! Open the NetCDF file
913      ! =============
914      ! do we have some free file identifier?
915      IF( MINVAL(iom_file(:)%nfid) /= 0 )   &
916         &   CALL ctl_stop( TRIM(clinfo), 'No more free file identifier', 'increase jpmax_files in iom_def' )
917      ! if no file was found...
918      IF( .NOT. llok ) THEN
919         IF( .NOT. llwrt ) THEN   ! we are in read mode
920            IF( llstop ) THEN   ;   CALL ctl_stop( TRIM(clinfo), 'File '//TRIM(cltmpn)//'* not found' )
921            ELSE                ;   istop = nstop + 1   ! make sure that istop /= nstop so we don't open the file
922            ENDIF
923         ELSE                     ! we are in write mode so we
924            clname = cltmpn       ! get back the file name without the cpu number
925         ENDIF
926      ELSE
927         IF( llwrt .AND. .NOT. ln_clobber ) THEN   ! we stop as we want to write in a new file
928            CALL ctl_stop( TRIM(clinfo), 'We want to write in a new file but '//TRIM(clname)//' already exists...' )
929            istop = nstop + 1                      ! make sure that istop /= nstop so we don't open the file
930         ELSEIF( llwrt ) THEN     ! the file exists and we are in write mode with permission to
931            clname = cltmpn       ! overwrite so get back the file name without the cpu number
932         ENDIF
933      ENDIF
934      IF( istop == nstop ) THEN   ! no error within this routine
935         CALL iom_nf90_open( clname, kiomid, llwrt, llok, idompar, kdlev = kdlev )
936      ENDIF
937      !
938   END SUBROUTINE iom_open
939
940
941   SUBROUTINE iom_close( kiomid )
942      !!--------------------------------------------------------------------
943      !!                   ***  SUBROUTINE  iom_close  ***
944      !!
945      !! ** Purpose : close an input file, or all files opened by iom
946      !!--------------------------------------------------------------------
947      INTEGER, INTENT(inout), OPTIONAL ::   kiomid   ! iom identifier of the file to be closed
948      !                                              ! return 0 when file is properly closed
949      !                                              ! No argument: all files opened by iom are closed
950
951      INTEGER ::   jf         ! dummy loop indices
952      INTEGER ::   i_s, i_e   ! temporary integer
953      CHARACTER(LEN=100)    ::   clinfo    ! info character
954      !---------------------------------------------------------------------
955      !
956      clinfo = '                    iom_close ~~~  '
957      IF( PRESENT(kiomid) ) THEN
958         i_s = kiomid
959         i_e = kiomid
960      ELSE
961         i_s = 1
962         i_e = jpmax_files
963      ENDIF
964
965      IF( i_s > 0 ) THEN
966         DO jf = i_s, i_e
967            IF( iom_file(jf)%nfid > 0 ) THEN
968               CALL iom_nf90_close( jf )
969               iom_file(jf)%nfid       = 0          ! free the id
970               IF( PRESENT(kiomid) )   kiomid = 0   ! return 0 as id to specify that the file was closed
971               IF(lwp) WRITE(numout,*) TRIM(clinfo)//' close file: '//TRIM(iom_file(jf)%name)//' ok'
972            ELSEIF( PRESENT(kiomid) ) THEN
973               WRITE(ctmp1,*) '--->',  kiomid
974               CALL ctl_stop( TRIM(clinfo)//' Invalid file identifier', ctmp1 )
975            ENDIF
976         END DO
977      ENDIF
978      !   
979   END SUBROUTINE iom_close
980
981
982   FUNCTION iom_varid ( kiomid, cdvar, kdimsz, kndims, lduld, ldstop ) 
983      !!-----------------------------------------------------------------------
984      !!                  ***  FUNCTION  iom_varid  ***
985      !!
986      !! ** Purpose : get the id of a variable in a file (return 0 if not found)
987      !!-----------------------------------------------------------------------
988      INTEGER              , INTENT(in   )           ::   kiomid   ! file Identifier
989      CHARACTER(len=*)     , INTENT(in   )           ::   cdvar    ! name of the variable
990      INTEGER, DIMENSION(:), INTENT(  out), OPTIONAL ::   kdimsz   ! size of each dimension
991      INTEGER              , INTENT(  out), OPTIONAL ::   kndims   ! number of dimensions
992      LOGICAL              , INTENT(  out), OPTIONAL ::   lduld    ! true if the last dimension is unlimited (time)
993      LOGICAL              , INTENT(in   ), OPTIONAL ::   ldstop   ! stop if looking for non-existing variable (default = .TRUE.)
994      !
995      INTEGER                        ::   iom_varid, iiv, i_nvd
996      LOGICAL                        ::   ll_fnd
997      CHARACTER(LEN=100)             ::   clinfo                   ! info character
998      LOGICAL                        ::   llstop                   ! local definition of ldstop
999      !!-----------------------------------------------------------------------
1000      iom_varid = 0                         ! default definition
1001      ! do we call ctl_stop if we look for non-existing variable?
1002      IF( PRESENT(ldstop) ) THEN   ;   llstop = ldstop
1003      ELSE                         ;   llstop = .TRUE.
1004      ENDIF
1005      !
1006      IF( kiomid > 0 ) THEN
1007         clinfo = 'iom_varid, file: '//trim(iom_file(kiomid)%name)//', var: '//trim(cdvar)
1008         IF( iom_file(kiomid)%nfid == 0 ) THEN
1009            CALL ctl_stop( trim(clinfo), 'the file is not open' )
1010         ELSE
1011            ll_fnd  = .FALSE.
1012            iiv = 0
1013            !
1014            DO WHILE ( .NOT.ll_fnd .AND. iiv < iom_file(kiomid)%nvars )
1015               iiv = iiv + 1
1016               ll_fnd  = ( TRIM(cdvar) == TRIM(iom_file(kiomid)%cn_var(iiv)) )
1017            END DO
1018            !
1019            IF( .NOT.ll_fnd ) THEN
1020               iiv = iiv + 1
1021               IF( iiv <= jpmax_vars ) THEN
1022                  iom_varid = iom_nf90_varid( kiomid, cdvar, iiv, kdimsz, kndims, lduld )
1023               ELSE
1024                  CALL ctl_stop( trim(clinfo), 'Too many variables in the file '//iom_file(kiomid)%name,   &
1025                        &                      'increase the parameter jpmax_vars')
1026               ENDIF
1027               IF( llstop .AND. iom_varid == -1 )   CALL ctl_stop( TRIM(clinfo)//' not found' ) 
1028            ELSE
1029               iom_varid = iiv
1030               IF( PRESENT(kdimsz) ) THEN
1031                  i_nvd = iom_file(kiomid)%ndims(iiv)
1032                  IF( i_nvd <= size(kdimsz) ) THEN
1033                     kdimsz(1:i_nvd) = iom_file(kiomid)%dimsz(1:i_nvd,iiv)
1034                  ELSE
1035                     WRITE(ctmp1,*) i_nvd, size(kdimsz)
1036                     CALL ctl_stop( trim(clinfo), 'error in kdimsz size'//trim(ctmp1) )
1037                  ENDIF
1038               ENDIF
1039               IF( PRESENT(kndims) )  kndims = iom_file(kiomid)%ndims(iiv)
1040               IF( PRESENT( lduld) )  lduld  = iom_file(kiomid)%luld( iiv)
1041            ENDIF
1042         ENDIF
1043      ENDIF
1044      !
1045   END FUNCTION iom_varid
1046
1047
1048   !!----------------------------------------------------------------------
1049   !!                   INTERFACE iom_get
1050   !!----------------------------------------------------------------------
1051   SUBROUTINE iom_g0d( kiomid, cdvar, pvar, ktime, ldxios )
1052      INTEGER         , INTENT(in   )                 ::   kiomid    ! Identifier of the file
1053      CHARACTER(len=*), INTENT(in   )                 ::   cdvar     ! Name of the variable
1054      REAL(wp)        , INTENT(  out)                 ::   pvar      ! read field
1055      INTEGER         , INTENT(in   ),     OPTIONAL   ::   ktime     ! record number
1056      LOGICAL         , INTENT(in   ),     OPTIONAL   ::   ldxios    ! use xios to read restart
1057      !
1058      INTEGER                                         ::   idvar     ! variable id
1059      INTEGER                                         ::   idmspc    ! number of spatial dimensions
1060      INTEGER         , DIMENSION(1)                  ::   itime     ! record number
1061      CHARACTER(LEN=100)                              ::   clinfo    ! info character
1062      CHARACTER(LEN=100)                              ::   clname    ! file name
1063      CHARACTER(LEN=1)                                ::   cldmspc   !
1064      LOGICAL                                         ::   llxios
1065      !
1066      llxios = .FALSE.
1067      IF( PRESENT(ldxios) ) llxios = ldxios
1068
1069      IF(.NOT.llxios) THEN  ! read data using default library
1070         itime = 1
1071         IF( PRESENT(ktime) ) itime = ktime
1072         !
1073         clname = iom_file(kiomid)%name
1074         clinfo = '          iom_g0d, file: '//trim(clname)//', var: '//trim(cdvar)
1075         !
1076         IF( kiomid > 0 ) THEN
1077            idvar = iom_varid( kiomid, cdvar )
1078            IF( iom_file(kiomid)%nfid > 0 .AND. idvar > 0 ) THEN
1079               idmspc = iom_file ( kiomid )%ndims( idvar )
1080               IF( iom_file(kiomid)%luld(idvar) )  idmspc = idmspc - 1
1081               WRITE(cldmspc , fmt='(i1)') idmspc
1082               IF( idmspc > 0 )  CALL ctl_stop( TRIM(clinfo), 'When reading to a 0D array, we do not accept data', &
1083                                    &                         'with 1 or more spatial dimensions: '//cldmspc//' were found.' , &
1084                                    &                         'Use ncwa -a to suppress the unnecessary dimensions' )
1085               CALL iom_nf90_get( kiomid, idvar, pvar, itime )
1086            ENDIF
1087         ENDIF
1088      ELSE
1089#if defined key_iomput
1090         IF(lwp) WRITE(numout,*) 'XIOS RST READ (0D): ', trim(cdvar)
1091         CALL xios_recv_field( trim(cdvar), pvar)
1092#else
1093         WRITE(ctmp1,*) 'Can not use XIOS in iom_g0d, file: '//trim(clname)//', var:'//trim(cdvar)
1094         CALL ctl_stop( 'iom_g0d', ctmp1 )
1095#endif
1096      ENDIF
1097   END SUBROUTINE iom_g0d
1098
1099   SUBROUTINE iom_g1d( kiomid, kdom, cdvar, pvar, ktime, kstart, kcount, ldxios )
1100      INTEGER         , INTENT(in   )                         ::   kiomid    ! Identifier of the file
1101      INTEGER         , INTENT(in   )                         ::   kdom      ! Type of domain to be read
1102      CHARACTER(len=*), INTENT(in   )                         ::   cdvar     ! Name of the variable
1103      REAL(wp)        , INTENT(  out), DIMENSION(:)           ::   pvar      ! read field
1104      INTEGER         , INTENT(in   )              , OPTIONAL ::   ktime     ! record number
1105      INTEGER         , INTENT(in   ), DIMENSION(1), OPTIONAL ::   kstart    ! start axis position of the reading
1106      INTEGER         , INTENT(in   ), DIMENSION(1), OPTIONAL ::   kcount    ! number of points in each axis
1107      LOGICAL         , INTENT(in   ),               OPTIONAL ::   ldxios    ! read data using XIOS
1108      !
1109      IF( kiomid > 0 ) THEN
1110         IF( iom_file(kiomid)%nfid > 0 ) CALL iom_get_123d( kiomid, kdom       , cdvar        , pv_r1d=pvar,   &
1111              &                                                     ktime=ktime, kstart=kstart, kcount=kcount, &
1112              &                                                     ldxios=ldxios )
1113      ENDIF
1114   END SUBROUTINE iom_g1d
1115
1116   SUBROUTINE iom_g2d( kiomid, kdom, cdvar, pvar, ktime, kstart, kcount, lrowattr, ldxios)
1117      INTEGER         , INTENT(in   )                           ::   kiomid    ! Identifier of the file
1118      INTEGER         , INTENT(in   )                           ::   kdom      ! Type of domain to be read
1119      CHARACTER(len=*), INTENT(in   )                           ::   cdvar     ! Name of the variable
1120      REAL(wp)        , INTENT(  out), DIMENSION(:,:)           ::   pvar      ! read field
1121      INTEGER         , INTENT(in   )                , OPTIONAL ::   ktime     ! record number
1122      INTEGER         , INTENT(in   ), DIMENSION(2)  , OPTIONAL ::   kstart    ! start axis position of the reading
1123      INTEGER         , INTENT(in   ), DIMENSION(2)  , OPTIONAL ::   kcount    ! number of points in each axis
1124      LOGICAL         , INTENT(in   )                , OPTIONAL ::   lrowattr  ! logical flag telling iom_get to
1125                                                                               ! look for and use a file attribute
1126                                                                               ! called open_ocean_jstart to set the start
1127                                                                               ! value for the 2nd dimension (netcdf only)
1128      LOGICAL         , INTENT(in   ),               OPTIONAL ::   ldxios      ! read data using XIOS
1129      !
1130      IF( kiomid > 0 ) THEN
1131         IF( iom_file(kiomid)%nfid > 0 ) CALL iom_get_123d( kiomid, kdom       , cdvar        , pv_r2d=pvar,   &
1132              &                                                     ktime=ktime, kstart=kstart, kcount=kcount, &
1133              &                                                     lrowattr=lrowattr,  ldxios=ldxios)
1134      ENDIF
1135   END SUBROUTINE iom_g2d
1136
1137   SUBROUTINE iom_g3d( kiomid, kdom, cdvar, pvar, ktime, kstart, kcount, lrowattr, ldxios )
1138      INTEGER         , INTENT(in   )                             ::   kiomid    ! Identifier of the file
1139      INTEGER         , INTENT(in   )                             ::   kdom      ! Type of domain to be read
1140      CHARACTER(len=*), INTENT(in   )                             ::   cdvar     ! Name of the variable
1141      REAL(wp)        , INTENT(  out), DIMENSION(:,:,:)           ::   pvar      ! read field
1142      INTEGER         , INTENT(in   )                  , OPTIONAL ::   ktime     ! record number
1143      INTEGER         , INTENT(in   ), DIMENSION(3)    , OPTIONAL ::   kstart    ! start axis position of the reading
1144      INTEGER         , INTENT(in   ), DIMENSION(3)    , OPTIONAL ::   kcount    ! number of points in each axis
1145      LOGICAL         , INTENT(in   )                  , OPTIONAL ::   lrowattr  ! logical flag telling iom_get to
1146                                                                                 ! look for and use a file attribute
1147                                                                                 ! called open_ocean_jstart to set the start
1148                                                                                 ! value for the 2nd dimension (netcdf only)
1149      LOGICAL         , INTENT(in   ),               OPTIONAL ::   ldxios        ! read data using XIOS
1150      !
1151      IF( kiomid > 0 ) THEN
1152         IF( iom_file(kiomid)%nfid > 0 ) CALL iom_get_123d( kiomid, kdom       , cdvar        , pv_r3d=pvar,   &
1153              &                                                     ktime=ktime, kstart=kstart, kcount=kcount, &
1154              &                                                     lrowattr=lrowattr, ldxios=ldxios )
1155      ENDIF
1156   END SUBROUTINE iom_g3d
1157
1158   SUBROUTINE iom_g4d( kiomid, kdom, cdvar, pvar, ktime, kstart, kcount, lrowattr, ldxios )
1159      INTEGER         , INTENT(in   )                             ::   kiomid    ! Identifier of the file
1160      INTEGER         , INTENT(in   )                             ::   kdom      ! Type of domain to be read
1161      CHARACTER(len=*), INTENT(in   )                             ::   cdvar     ! Name of the variable
1162      REAL(wp)        , INTENT(  out), DIMENSION(:,:,:,:)         ::   pvar      ! read field
1163      INTEGER         , INTENT(in   )                  , OPTIONAL ::   ktime     ! record number
1164      INTEGER         , INTENT(in   ), DIMENSION(4)    , OPTIONAL ::   kstart    ! start axis position of the reading
1165      INTEGER         , INTENT(in   ), DIMENSION(4)    , OPTIONAL ::   kcount    ! number of points in each axis
1166      LOGICAL         , INTENT(in   )                  , OPTIONAL ::   lrowattr  ! logical flag telling iom_get to
1167                                                                                 ! look for and use a file attribute
1168                                                                                 ! called open_ocean_jstart to set the start
1169                                                                                 ! value for the 2nd dimension (netcdf only)
1170      LOGICAL         , INTENT(in   ),               OPTIONAL ::   ldxios        ! read data using XIOS
1171      !
1172      IF( kiomid > 0 ) THEN
1173         IF( iom_file(kiomid)%nfid > 0 ) CALL iom_get_123d( kiomid, kdom       , cdvar        , pv_r4d=pvar,   &
1174              &                                                     ktime=ktime, kstart=kstart, kcount=kcount, &
1175              &                                                     lrowattr=lrowattr, ldxios=ldxios )
1176      ENDIF
1177   END SUBROUTINE iom_g4d
1178   !!----------------------------------------------------------------------
1179
1180   SUBROUTINE iom_get_123d( kiomid, kdom  , cdvar ,   &
1181         &                  pv_r1d, pv_r2d, pv_r3d, pv_r4d,  &
1182         &                  ktime , kstart, kcount,   &
1183         &                  lrowattr, ldxios        )
1184      !!-----------------------------------------------------------------------
1185      !!                  ***  ROUTINE  iom_get_123d  ***
1186      !!
1187      !! ** Purpose : read a 1D/2D/3D variable
1188      !!
1189      !! ** Method : read ONE record at each CALL
1190      !!-----------------------------------------------------------------------
1191      INTEGER                    , INTENT(in   )           ::   kiomid     ! Identifier of the file
1192      INTEGER                    , INTENT(in   )           ::   kdom       ! Type of domain to be read
1193      CHARACTER(len=*)           , INTENT(in   )           ::   cdvar      ! Name of the variable
1194      REAL(wp), DIMENSION(:)     , INTENT(  out), OPTIONAL ::   pv_r1d     ! read field (1D case)
1195      REAL(wp), DIMENSION(:,:)   , INTENT(  out), OPTIONAL ::   pv_r2d     ! read field (2D case)
1196      REAL(wp), DIMENSION(:,:,:) , INTENT(  out), OPTIONAL ::   pv_r3d     ! read field (3D case)
1197      REAL(wp), DIMENSION(:,:,:,:),INTENT(  out), OPTIONAL ::   pv_r4d     ! read field (4D case)
1198      INTEGER                    , INTENT(in   ), OPTIONAL ::   ktime      ! record number
1199      INTEGER , DIMENSION(:)     , INTENT(in   ), OPTIONAL ::   kstart     ! start position of the reading in each axis
1200      INTEGER , DIMENSION(:)     , INTENT(in   ), OPTIONAL ::   kcount     ! number of points to be read in each axis
1201      LOGICAL                    , INTENT(in   ), OPTIONAL ::   lrowattr   ! logical flag telling iom_get to
1202                                                                           ! look for and use a file attribute
1203                                                                           ! called open_ocean_jstart to set the start
1204                                                                           ! value for the 2nd dimension (netcdf only)
1205      LOGICAL                    , INTENT(in   ), OPTIONAL ::   ldxios     ! use XIOS to read restart
1206      !
1207      LOGICAL                        ::   llxios       ! local definition for XIOS read
1208      LOGICAL                        ::   llnoov      ! local definition to read overlap
1209      LOGICAL                        ::   luse_jattr  ! local definition to read open_ocean_jstart file attribute
1210      INTEGER                        ::   jstartrow   ! start point for 2nd dimension optionally set by file attribute
1211      INTEGER                        ::   jl          ! loop on number of dimension
1212      INTEGER                        ::   idom        ! type of domain
1213      INTEGER                        ::   idvar       ! id of the variable
1214      INTEGER                        ::   inbdim      ! number of dimensions of the variable
1215      INTEGER                        ::   idmspc      ! number of spatial dimensions
1216      INTEGER                        ::   itime       ! record number
1217      INTEGER                        ::   istop       ! temporary value of nstop
1218      INTEGER                        ::   ix1, ix2, iy1, iy2   ! subdomain indexes
1219      INTEGER                        ::   ji, jj      ! loop counters
1220      INTEGER                        ::   irankpv     !
1221      INTEGER                        ::   ind1, ind2  ! substring index
1222      INTEGER, DIMENSION(jpmax_dims) ::   istart      ! starting point to read for each axis
1223      INTEGER, DIMENSION(jpmax_dims) ::   icnt        ! number of value to read along each axis
1224      INTEGER, DIMENSION(jpmax_dims) ::   idimsz      ! size of the dimensions of the variable
1225      INTEGER, DIMENSION(jpmax_dims) ::   ishape      ! size of the dimensions of the variable
1226      REAL(wp)                       ::   zscf, zofs  ! sacle_factor and add_offset
1227      INTEGER                        ::   itmp        ! temporary integer
1228      CHARACTER(LEN=256)             ::   clinfo      ! info character
1229      CHARACTER(LEN=256)             ::   clname      ! file name
1230      CHARACTER(LEN=1)               ::   clrankpv, cldmspc      !
1231      LOGICAL                        ::   ll_depth_spec ! T => if kstart, kcount present then *only* use values for 3rd spatial dimension.
1232      INTEGER                        ::   inlev       ! number of levels for 3D data
1233      REAL(wp)                       ::   gma, gmi
1234      !---------------------------------------------------------------------
1235      !
1236      inlev = -1
1237      IF( PRESENT(pv_r3d) )   inlev = SIZE(pv_r3d, 3)
1238      !
1239      llxios = .FALSE.
1240      if(PRESENT(ldxios)) llxios = ldxios
1241      idvar = iom_varid( kiomid, cdvar ) 
1242      idom = kdom
1243      !
1244      IF(.NOT.llxios) THEN
1245         clname = iom_file(kiomid)%name   !   esier to read
1246         clinfo = '          iom_get_123d, file: '//trim(clname)//', var: '//trim(cdvar)
1247         ! local definition of the domain ?
1248         ! do we read the overlap
1249         ! ugly patch SM+JMM+RB to overwrite global definition in some cases
1250         llnoov = (jpni * jpnj ) == jpnij .AND. .NOT. lk_agrif 
1251         ! check kcount and kstart optionals parameters...
1252         IF( PRESENT(kcount) .AND. (.NOT. PRESENT(kstart)) ) CALL ctl_stop(trim(clinfo), 'kcount present needs kstart present')
1253         IF( PRESENT(kstart) .AND. (.NOT. PRESENT(kcount)) ) CALL ctl_stop(trim(clinfo), 'kstart present needs kcount present')
1254         IF( PRESENT(kstart) .AND. idom /= jpdom_unknown .AND.  idom /= jpdom_autoglo_xy  ) &
1255     &          CALL ctl_stop(trim(clinfo), 'kstart present needs kdom = jpdom_unknown or kdom = jpdom_autoglo_xy')
1256
1257         luse_jattr = .false.
1258         IF( PRESENT(lrowattr) ) THEN
1259            IF( lrowattr .AND. idom /= jpdom_data   ) CALL ctl_stop(trim(clinfo), 'lrowattr present and true needs kdom = jpdom_data')
1260            IF( lrowattr .AND. idom == jpdom_data   ) luse_jattr = .true.
1261         ENDIF
1262
1263         ! Search for the variable in the data base (eventually actualize data)
1264         istop = nstop
1265         !
1266         IF( idvar > 0 ) THEN
1267            ! to write iom_file(kiomid)%dimsz in a shorter way !
1268            idimsz(:) = iom_file(kiomid)%dimsz(:, idvar) 
1269            inbdim = iom_file(kiomid)%ndims(idvar)            ! number of dimensions in the file
1270            idmspc = inbdim                                   ! number of spatial dimensions in the file
1271            IF( iom_file(kiomid)%luld(idvar) )   idmspc = inbdim - 1
1272            IF( idmspc > 3 )   CALL ctl_stop(trim(clinfo), 'the file has more than 3 spatial dimensions this case is not coded...') 
1273            !
1274            ! update idom definition...
1275            ! Identify the domain in case of jpdom_auto(glo/dta) definition
1276            IF( idom == jpdom_autoglo_xy ) THEN
1277               ll_depth_spec = .TRUE.
1278               idom = jpdom_autoglo
1279            ELSE
1280               ll_depth_spec = .FALSE.
1281            ENDIF
1282            IF( idom == jpdom_autoglo .OR. idom == jpdom_autodta ) THEN           
1283               IF( idom == jpdom_autoglo ) THEN   ;   idom = jpdom_global 
1284               ELSE                               ;   idom = jpdom_data
1285               ENDIF
1286               ind1 = INDEX( clname, '_', back = .TRUE. ) + 1
1287               ind2 = INDEX( clname, '.', back = .TRUE. ) - 1
1288               IF( ind2 > ind1 ) THEN   ;   IF( VERIFY( clname(ind1:ind2), '0123456789' ) == 0 )   idom = jpdom_local   ;   ENDIF
1289            ENDIF
1290            ! Identify the domain in case of jpdom_local definition
1291            IF( idom == jpdom_local ) THEN
1292               IF(     idimsz(1) == jpi               .AND. idimsz(2) == jpj               ) THEN   ;   idom = jpdom_local_full
1293               ELSEIF( idimsz(1) == nlci              .AND. idimsz(2) == nlcj              ) THEN   ;   idom = jpdom_local_noextra
1294               ELSEIF( idimsz(1) == (nlei - nldi + 1) .AND. idimsz(2) == (nlej - nldj + 1) ) THEN   ;   idom = jpdom_local_noovlap
1295               ELSE   ;   CALL ctl_stop( trim(clinfo), 'impossible to identify the local domain' )
1296               ENDIF
1297            ENDIF
1298            !
1299            ! check the consistency between input array and data rank in the file
1300            !
1301            ! initializations
1302            itime = 1
1303            IF( PRESENT(ktime) ) itime = ktime
1304            !
1305            irankpv = 1 * COUNT( (/PRESENT(pv_r1d)/) ) + 2 * COUNT( (/PRESENT(pv_r2d)/) ) + 3 * COUNT( (/PRESENT(pv_r3d)/) )
1306            WRITE(clrankpv, fmt='(i1)') irankpv
1307            WRITE(cldmspc , fmt='(i1)') idmspc
1308            !
1309            IF(     idmspc <  irankpv ) THEN
1310               CALL ctl_stop( TRIM(clinfo), 'The file has only '//cldmspc//' spatial dimension',   &
1311                  &                         'it is impossible to read a '//clrankpv//'D array from this file...' )
1312            ELSEIF( idmspc == irankpv ) THEN
1313               IF( PRESENT(pv_r1d) .AND. idom /= jpdom_unknown )   &
1314                  &   CALL ctl_stop( TRIM(clinfo), 'case not coded...You must use jpdom_unknown' )
1315            ELSEIF( idmspc >  irankpv ) THEN
1316                  IF( PRESENT(pv_r2d) .AND. itime == 1 .AND. idimsz(3) == 1 .AND. idmspc == 3 ) THEN
1317                     CALL ctl_warn( trim(clinfo), '2D array but 3 spatial dimensions for the data...'              ,   &
1318                           &         'As the size of the z dimension is 1 and as we try to read the first record, ',   &
1319                           &         'we accept this case, even if there is a possible mix-up between z and time dimension' )   
1320                     idmspc = idmspc - 1
1321                  ELSE
1322                     CALL ctl_stop( TRIM(clinfo), 'To keep iom lisibility, when reading a '//clrankpv//'D array,'         ,   &
1323                        &                         'we do not accept data with '//cldmspc//' spatial dimensions',   &
1324                        &                         'Use ncwa -a to suppress the unnecessary dimensions' )
1325                  ENDIF
1326            ENDIF
1327            !
1328            ! definition of istart and icnt
1329            !
1330            icnt  (:) = 1
1331            istart(:) = 1
1332            istart(idmspc+1) = itime
1333   
1334            IF( PRESENT(kstart) .AND. .NOT. ll_depth_spec ) THEN
1335               istart(1:idmspc) = kstart(1:idmspc) 
1336               icnt  (1:idmspc) = kcount(1:idmspc)
1337            ELSE
1338               IF(idom == jpdom_unknown ) THEN
1339                  icnt(1:idmspc) = idimsz(1:idmspc)
1340               ELSE
1341                  IF( .NOT. PRESENT(pv_r1d) ) THEN   !   not a 1D array
1342                     IF(     idom == jpdom_data    ) THEN
1343                        jstartrow = 1
1344                        IF( luse_jattr ) THEN
1345                           CALL iom_getatt(kiomid, 'open_ocean_jstart', jstartrow ) ! -999 is returned if the attribute is not found
1346                           jstartrow = MAX(1,jstartrow)
1347                        ENDIF
1348                        istart(1:2) = (/ mig(1), mjg(1) + jstartrow - 1 /)  ! icnt(1:2) done below
1349                     ELSEIF( idom == jpdom_global  ) THEN ; istart(1:2) = (/ nimpp , njmpp  /)  ! icnt(1:2) done below
1350                     ENDIF
1351                     ! we do not read the overlap                     -> we start to read at nldi, nldj
1352! JMM + SM: ugly patch before getting the new version of lib_mpp)
1353!                  IF( idom /= jpdom_local_noovlap )   istart(1:2) = istart(1:2) + (/ nldi - 1, nldj - 1 /)
1354                     IF( llnoov .AND. idom /= jpdom_local_noovlap ) istart(1:2) = istart(1:2) + (/ nldi - 1, nldj - 1 /)
1355                  ! we do not read the overlap and the extra-halos -> from nldi to nlei and from nldj to nlej
1356! JMM + SM: ugly patch before getting the new version of lib_mpp)
1357!                  icnt(1:2) = (/ nlei - nldi + 1, nlej - nldj + 1 /)
1358                     IF( llnoov ) THEN   ;   icnt(1:2) = (/ nlei - nldi + 1, nlej - nldj + 1 /)
1359                     ELSE                ;   icnt(1:2) = (/ nlci           , nlcj            /)
1360                     ENDIF
1361                     IF( PRESENT(pv_r3d) ) THEN
1362                        IF( idom == jpdom_data ) THEN                        ;                               icnt(3) = inlev
1363                        ELSEIF( ll_depth_spec .AND. PRESENT(kstart) ) THEN   ;   istart(3) = kstart(3)   ;   icnt(3) = kcount(3)
1364                        ELSE                                                 ;                               icnt(3) = inlev
1365                        ENDIF
1366                     ENDIF
1367                  ENDIF
1368               ENDIF
1369            ENDIF
1370
1371            ! check that istart and icnt can be used with this file
1372            !-
1373            DO jl = 1, jpmax_dims
1374               itmp = istart(jl)+icnt(jl)-1
1375               IF( itmp > idimsz(jl) .AND. idimsz(jl) /= 0 ) THEN
1376                  WRITE( ctmp1, FMT="('(istart(', i1, ') + icnt(', i1, ') - 1) = ', i5)" ) jl, jl, itmp
1377                  WRITE( ctmp2, FMT="(' is larger than idimsz(', i1,') = ', i5)"         ) jl, idimsz(jl)
1378                  CALL ctl_stop( trim(clinfo), 'start and count too big regarding to the size of the data, ', ctmp1, ctmp2 )     
1379               ENDIF
1380            END DO
1381
1382            ! check that icnt matches the input array
1383            !-     
1384            IF( idom == jpdom_unknown ) THEN
1385               IF( irankpv == 1 )        ishape(1:1) = SHAPE(pv_r1d)
1386               IF( irankpv == 2 )        ishape(1:2) = SHAPE(pv_r2d)
1387               IF( irankpv == 3 )        ishape(1:3) = SHAPE(pv_r3d)
1388               ctmp1 = 'd'
1389            ELSE
1390               IF( irankpv == 2 ) THEN
1391! JMM + SM: ugly patch before getting the new version of lib_mpp)
1392!               ishape(1:2) = SHAPE(pv_r2d(nldi:nlei,nldj:nlej  ))   ;   ctmp1 = 'd(nldi:nlei,nldj:nlej)'
1393                  IF( llnoov ) THEN ; ishape(1:2)=SHAPE(pv_r2d(nldi:nlei,nldj:nlej  )) ; ctmp1='d(nldi:nlei,nldj:nlej)'
1394                  ELSE              ; ishape(1:2)=SHAPE(pv_r2d(1   :nlci,1   :nlcj  )) ; ctmp1='d(1:nlci,1:nlcj)'
1395                  ENDIF
1396               ENDIF
1397               IF( irankpv == 3 ) THEN 
1398! JMM + SM: ugly patch before getting the new version of lib_mpp)
1399!               ishape(1:3) = SHAPE(pv_r3d(nldi:nlei,nldj:nlej,:))   ;   ctmp1 = 'd(nldi:nlei,nldj:nlej,:)'
1400                  IF( llnoov ) THEN ; ishape(1:3)=SHAPE(pv_r3d(nldi:nlei,nldj:nlej,:)) ; ctmp1='d(nldi:nlei,nldj:nlej,:)'
1401                  ELSE              ; ishape(1:3)=SHAPE(pv_r3d(1   :nlci,1   :nlcj,:)) ; ctmp1='d(1:nlci,1:nlcj,:)'
1402                  ENDIF
1403               ENDIF
1404            ENDIF
1405         
1406            DO jl = 1, irankpv
1407               WRITE( ctmp2, FMT="(', ', i1,'): ', i5,' /= icnt(', i1,'):', i5)" ) jl, ishape(jl), jl, icnt(jl)
1408               IF( ishape(jl) /= icnt(jl) )   CALL ctl_stop( TRIM(clinfo), 'size(pv_r'//clrankpv//TRIM(ctmp1)//TRIM(ctmp2) )
1409            END DO
1410
1411         ENDIF
1412
1413         ! read the data
1414         !-     
1415         IF( idvar > 0 .AND. istop == nstop ) THEN   ! no additional errors until this point...
1416            !
1417         ! find the right index of the array to be read
1418! JMM + SM: ugly patch before getting the new version of lib_mpp)
1419!         IF( idom /= jpdom_unknown ) THEN   ;   ix1 = nldi   ;   ix2 = nlei      ;   iy1 = nldj   ;   iy2 = nlej
1420!         ELSE                               ;   ix1 = 1      ;   ix2 = icnt(1)   ;   iy1 = 1      ;   iy2 = icnt(2)
1421!         ENDIF
1422            IF( llnoov ) THEN
1423               IF( idom /= jpdom_unknown ) THEN   ;   ix1 = nldi   ;   ix2 = nlei      ;   iy1 = nldj   ;   iy2 = nlej
1424               ELSE                               ;   ix1 = 1      ;   ix2 = icnt(1)   ;   iy1 = 1      ;   iy2 = icnt(2)
1425               ENDIF
1426            ELSE
1427               IF( idom /= jpdom_unknown ) THEN   ;   ix1 = 1      ;   ix2 = nlci      ;   iy1 = 1      ;   iy2 = nlcj
1428               ELSE                               ;   ix1 = 1      ;   ix2 = icnt(1)   ;   iy1 = 1      ;   iy2 = icnt(2)
1429               ENDIF
1430            ENDIF
1431     
1432            CALL iom_nf90_get( kiomid, idvar, inbdim, istart, icnt, ix1, ix2, iy1, iy2, pv_r1d, pv_r2d, pv_r3d )
1433
1434            IF( istop == nstop ) THEN   ! no additional errors until this point...
1435               IF(lwp) WRITE(numout,"(10x,' read ',a,' (rec: ',i6,') in ',a,' ok')") TRIM(cdvar), itime, TRIM(iom_file(kiomid)%name)
1436             
1437               !--- overlap areas and extra hallows (mpp)
1438               IF(     PRESENT(pv_r2d) .AND. idom /= jpdom_unknown ) THEN
1439                  CALL lbc_lnk( 'iom', pv_r2d,'Z', -999., kfillmode = jpfillnothing )
1440               ELSEIF( PRESENT(pv_r3d) .AND. idom /= jpdom_unknown ) THEN
1441                  ! this if could be simplified with the new lbc_lnk that works with any size of the 3rd dimension
1442                  IF( icnt(3) == inlev ) THEN
1443                     CALL lbc_lnk( 'iom', pv_r3d,'Z', -999., kfillmode = jpfillnothing )
1444                  ELSE   ! put some arbitrary value (a call to lbc_lnk will be done later...)
1445                     DO jj = nlcj+1, jpj   ;   pv_r3d(1:nlci, jj, :) = pv_r3d(1:nlci, nlej, :)   ;   END DO
1446                     DO ji = nlci+1, jpi   ;   pv_r3d(ji    , : , :) = pv_r3d(nlei  , :   , :)   ;   END DO
1447                  ENDIF
1448               ENDIF
1449               !
1450            ELSE
1451               ! return if istop == nstop is false
1452               RETURN
1453            ENDIF
1454         ELSE
1455            ! return if statment idvar > 0 .AND. istop == nstop is false
1456            RETURN
1457         ENDIF
1458         !
1459      ELSE        ! read using XIOS. Only if KEY_IOMPUT is defined
1460#if defined key_iomput
1461!would be good to be able to check which context is active and swap only if current is not restart
1462         IF( PRESENT(pv_r4d) ) THEN
1463            pv_r4d(:, :, :, :) = 0.
1464            if(lwp) write(numout,*) 'XIOS RST READ (4D): ',trim(cdvar)
1465            CALL xios_recv_field( trim(cdvar), pv_r4d)
1466            IF(idom /= jpdom_unknown ) then
1467                CALL lbc_lnk( 'iom', pv_r4d,'Z', -999., kfillmode = jpfillnothing)
1468            ENDIF
1469         ELSEIF( PRESENT(pv_r3d) ) THEN
1470            pv_r3d(:, :, :) = 0.
1471            if(lwp) write(numout,*) 'XIOS RST READ (3D): ',trim(cdvar)
1472            CALL xios_recv_field( trim(cdvar), pv_r3d)
1473            IF(idom /= jpdom_unknown ) then
1474                CALL lbc_lnk( 'iom', pv_r3d,'Z', -999., kfillmode = jpfillnothing)
1475            ENDIF
1476         ELSEIF( PRESENT(pv_r2d) ) THEN
1477            pv_r2d(:, :) = 0.
1478            if(lwp) write(numout,*) 'XIOS RST READ (2D): ', trim(cdvar)
1479            CALL xios_recv_field( trim(cdvar), pv_r2d)
1480            IF(idom /= jpdom_unknown ) THEN
1481                CALL lbc_lnk('iom', pv_r2d,'Z',-999., kfillmode = jpfillnothing)
1482            ENDIF
1483         ELSEIF( PRESENT(pv_r1d) ) THEN
1484            pv_r1d(:) = 0.
1485            if(lwp) write(numout,*) 'XIOS RST READ (1D): ', trim(cdvar)
1486            CALL xios_recv_field( trim(cdvar), pv_r1d)
1487         ENDIF
1488#else
1489         istop = istop + 1 
1490         clinfo = 'Can not use XIOS in iom_get_123d, file: '//trim(clname)//', var:'//trim(cdvar)
1491#endif
1492      ENDIF
1493!some final adjustments
1494      ! C1D case : always call lbc_lnk to replicate the central value over the whole 3X3 domain
1495      IF( lk_c1d .AND. PRESENT(pv_r2d) )   CALL lbc_lnk( 'iom', pv_r2d,'Z',1. )
1496      IF( lk_c1d .AND. PRESENT(pv_r3d) )   CALL lbc_lnk( 'iom', pv_r3d,'Z',1. )
1497
1498      !--- Apply scale_factor and offset
1499      zscf = iom_file(kiomid)%scf(idvar)      ! scale factor
1500      zofs = iom_file(kiomid)%ofs(idvar)      ! offset
1501      IF(     PRESENT(pv_r1d) ) THEN
1502         IF( zscf /= 1. )   pv_r1d(:) = pv_r1d(:) * zscf 
1503         IF( zofs /= 0. )   pv_r1d(:) = pv_r1d(:) + zofs
1504      ELSEIF( PRESENT(pv_r2d) ) THEN
1505         IF( zscf /= 1.)   pv_r2d(:,:) = pv_r2d(:,:) * zscf
1506         IF( zofs /= 0.)   pv_r2d(:,:) = pv_r2d(:,:) + zofs
1507      ELSEIF( PRESENT(pv_r3d) ) THEN
1508         IF( zscf /= 1.)   pv_r3d(:,:,:) = pv_r3d(:,:,:) * zscf
1509         IF( zofs /= 0.)   pv_r3d(:,:,:) = pv_r3d(:,:,:) + zofs
1510      ELSEIF( PRESENT(pv_r4d) ) THEN
1511         IF( zscf /= 1.)   pv_r4d(:,:,:,:) = pv_r4d(:,:,:,:) * zscf
1512         IF( zofs /= 0.)   pv_r4d(:,:,:,:) = pv_r4d(:,:,:,:) + zofs
1513      ENDIF
1514      !
1515   END SUBROUTINE iom_get_123d
1516
1517
1518   FUNCTION iom_getszuld ( kiomid ) 
1519      !!-----------------------------------------------------------------------
1520      !!                  ***  FUNCTION  iom_getszuld  ***
1521      !!
1522      !! ** Purpose : get the size of the unlimited dimension in a file
1523      !!              (return -1 if not found)
1524      !!-----------------------------------------------------------------------
1525      INTEGER, INTENT(in   ) ::   kiomid   ! file Identifier
1526      !
1527      INTEGER                ::   iom_getszuld
1528      !!-----------------------------------------------------------------------
1529      iom_getszuld = -1
1530      IF( kiomid > 0 ) THEN
1531         IF( iom_file(kiomid)%iduld > 0 )   iom_getszuld = iom_file(kiomid)%lenuld
1532      ENDIF
1533   END FUNCTION iom_getszuld
1534   
1535
1536   !!----------------------------------------------------------------------
1537   !!                   INTERFACE iom_chkatt
1538   !!----------------------------------------------------------------------
1539   SUBROUTINE iom_chkatt( kiomid, cdatt, llok, ksize, cdvar )
1540      INTEGER         , INTENT(in   )                 ::   kiomid    ! Identifier of the file
1541      CHARACTER(len=*), INTENT(in   )                 ::   cdatt     ! Name of the attribute
1542      LOGICAL         , INTENT(  out)                 ::   llok      ! Error code
1543      INTEGER         , INTENT(  out), OPTIONAL       ::   ksize     ! Size of the attribute array
1544      CHARACTER(len=*), INTENT(in   ), OPTIONAL       ::   cdvar     ! Name of the variable
1545      !
1546      IF( kiomid > 0 ) THEN
1547         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_chkatt( kiomid, cdatt, llok, ksize=ksize, cdvar=cdvar )
1548      ENDIF
1549      !
1550   END SUBROUTINE iom_chkatt
1551
1552   !!----------------------------------------------------------------------
1553   !!                   INTERFACE iom_getatt
1554   !!----------------------------------------------------------------------
1555   SUBROUTINE iom_g0d_iatt( kiomid, cdatt, katt0d, cdvar )
1556      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1557      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1558      INTEGER               , INTENT(  out)           ::   katt0d    ! read field
1559      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1560      !
1561      IF( kiomid > 0 ) THEN
1562         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  katt0d =  katt0d, cdvar=cdvar )
1563      ENDIF
1564   END SUBROUTINE iom_g0d_iatt
1565
1566   SUBROUTINE iom_g1d_iatt( kiomid, cdatt, katt1d, cdvar )
1567      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1568      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1569      INTEGER, DIMENSION(:) , INTENT(  out)           ::   katt1d    ! read field
1570      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1571      !
1572      IF( kiomid > 0 ) THEN
1573         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  katt1d =  katt1d, cdvar=cdvar )
1574      ENDIF
1575   END SUBROUTINE iom_g1d_iatt
1576
1577   SUBROUTINE iom_g0d_ratt( kiomid, cdatt, patt0d, cdvar )
1578      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1579      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1580      REAL(wp)              , INTENT(  out)           ::   patt0d    ! read field
1581      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1582      !
1583      IF( kiomid > 0 ) THEN
1584         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  patt0d =  patt0d, cdvar=cdvar )
1585      ENDIF
1586   END SUBROUTINE iom_g0d_ratt
1587
1588   SUBROUTINE iom_g1d_ratt( kiomid, cdatt, patt1d, cdvar )
1589      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1590      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1591      REAL(wp), DIMENSION(:), INTENT(  out)           ::   patt1d    ! read field
1592      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1593      !
1594      IF( kiomid > 0 ) THEN
1595         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  patt1d =  patt1d, cdvar=cdvar )
1596      ENDIF
1597   END SUBROUTINE iom_g1d_ratt
1598   
1599   SUBROUTINE iom_g0d_catt( kiomid, cdatt, cdatt0d, cdvar )
1600      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1601      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1602      CHARACTER(len=*)      , INTENT(  out)           ::   cdatt0d   ! read field
1603      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1604      !
1605      IF( kiomid > 0 ) THEN
1606         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt, cdatt0d = cdatt0d, cdvar=cdvar )
1607      ENDIF
1608   END SUBROUTINE iom_g0d_catt
1609
1610
1611   !!----------------------------------------------------------------------
1612   !!                   INTERFACE iom_putatt
1613   !!----------------------------------------------------------------------
1614   SUBROUTINE iom_p0d_iatt( kiomid, cdatt, katt0d, cdvar )
1615      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1616      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1617      INTEGER               , INTENT(in   )           ::   katt0d    ! written field
1618      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1619      !
1620      IF( kiomid > 0 ) THEN
1621         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  katt0d =  katt0d, cdvar=cdvar )
1622      ENDIF
1623   END SUBROUTINE iom_p0d_iatt
1624
1625   SUBROUTINE iom_p1d_iatt( kiomid, cdatt, katt1d, cdvar )
1626      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1627      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1628      INTEGER, DIMENSION(:) , INTENT(in   )           ::   katt1d    ! written field
1629      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1630      !
1631      IF( kiomid > 0 ) THEN
1632         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  katt1d =  katt1d, cdvar=cdvar )
1633      ENDIF
1634   END SUBROUTINE iom_p1d_iatt
1635
1636   SUBROUTINE iom_p0d_ratt( kiomid, cdatt, patt0d, cdvar )
1637      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1638      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1639      REAL(wp)              , INTENT(in   )           ::   patt0d    ! written field
1640      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1641      !
1642      IF( kiomid > 0 ) THEN
1643         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  patt0d =  patt0d, cdvar=cdvar )
1644      ENDIF
1645   END SUBROUTINE iom_p0d_ratt
1646
1647   SUBROUTINE iom_p1d_ratt( kiomid, cdatt, patt1d, cdvar )
1648      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1649      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1650      REAL(wp), DIMENSION(:), INTENT(in   )           ::   patt1d    ! written field
1651      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1652      !
1653      IF( kiomid > 0 ) THEN
1654         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  patt1d =  patt1d, cdvar=cdvar )
1655      ENDIF
1656   END SUBROUTINE iom_p1d_ratt
1657   
1658   SUBROUTINE iom_p0d_catt( kiomid, cdatt, cdatt0d, cdvar )
1659      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1660      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1661      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt0d   ! written field
1662      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1663      !
1664      IF( kiomid > 0 ) THEN
1665         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt, cdatt0d = cdatt0d, cdvar=cdvar )
1666      ENDIF
1667   END SUBROUTINE iom_p0d_catt
1668
1669
1670   !!----------------------------------------------------------------------
1671   !!                   INTERFACE iom_rstput
1672   !!----------------------------------------------------------------------
1673   SUBROUTINE iom_rp0d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1674      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1675      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1676      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1677      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1678      REAL(wp)        , INTENT(in)                         ::   pvar     ! written field
1679      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1680      LOGICAL, OPTIONAL :: ldxios   ! xios write flag
1681      LOGICAL :: llx                ! local xios write flag
1682      INTEGER :: ivid   ! variable id
1683
1684      llx = .FALSE.
1685      IF(PRESENT(ldxios)) llx = ldxios
1686      IF( llx ) THEN
1687#ifdef key_iomput
1688      IF( kt == kwrite ) THEN
1689          IF(lwp) write(numout,*) 'RESTART: write (XIOS 0D) ',trim(cdvar)
1690          CALL xios_send_field(trim(cdvar), pvar)
1691      ENDIF
1692#endif
1693      ELSE
1694         IF( kiomid > 0 ) THEN
1695            IF( iom_file(kiomid)%nfid > 0 ) THEN
1696               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1697               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r0d = pvar )
1698            ENDIF
1699         ENDIF
1700      ENDIF
1701   END SUBROUTINE iom_rp0d
1702
1703   SUBROUTINE iom_rp1d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1704      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1705      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1706      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1707      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1708      REAL(wp)        , INTENT(in), DIMENSION(          :) ::   pvar     ! written field
1709      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1710      LOGICAL, OPTIONAL                                    ::   ldxios   ! xios write flag
1711      LOGICAL :: llx                ! local xios write flag
1712      INTEGER :: ivid   ! variable id
1713
1714      llx = .FALSE.
1715      IF(PRESENT(ldxios)) llx = ldxios
1716      IF( llx ) THEN
1717#ifdef key_iomput
1718      IF( kt == kwrite ) THEN
1719         IF(lwp) write(numout,*) 'RESTART: write (XIOS 1D) ',trim(cdvar)
1720         CALL xios_send_field(trim(cdvar), pvar)
1721      ENDIF
1722#endif
1723      ELSE
1724         IF( kiomid > 0 ) THEN
1725            IF( iom_file(kiomid)%nfid > 0 ) THEN
1726               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1727               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r1d = pvar )
1728            ENDIF
1729         ENDIF
1730      ENDIF
1731   END SUBROUTINE iom_rp1d
1732
1733   SUBROUTINE iom_rp2d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1734      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1735      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1736      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1737      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1738      REAL(wp)        , INTENT(in), DIMENSION(:,    :    ) ::   pvar     ! written field
1739      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1740      LOGICAL, OPTIONAL :: ldxios   ! xios write flag
1741      LOGICAL :: llx
1742      INTEGER :: ivid   ! variable id
1743
1744      llx = .FALSE.
1745      IF(PRESENT(ldxios)) llx = ldxios
1746      IF( llx ) THEN
1747#ifdef key_iomput
1748      IF( kt == kwrite ) THEN
1749         IF(lwp) write(numout,*) 'RESTART: write (XIOS 2D) ',trim(cdvar)
1750         CALL xios_send_field(trim(cdvar), pvar)
1751      ENDIF
1752#endif
1753      ELSE
1754         IF( kiomid > 0 ) THEN
1755            IF( iom_file(kiomid)%nfid > 0 ) THEN
1756               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1757               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r2d = pvar )
1758            ENDIF
1759         ENDIF
1760      ENDIF
1761   END SUBROUTINE iom_rp2d
1762
1763   SUBROUTINE iom_rp3d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1764      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1765      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1766      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1767      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1768      REAL(wp)        , INTENT(in),       DIMENSION(:,:,:) ::   pvar     ! written field
1769      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1770      LOGICAL, OPTIONAL :: ldxios   ! xios write flag
1771      LOGICAL :: llx                 ! local xios write flag
1772      INTEGER :: ivid   ! variable id
1773
1774      llx = .FALSE.
1775      IF(PRESENT(ldxios)) llx = ldxios
1776      IF( llx ) THEN
1777#ifdef key_iomput
1778      IF( kt == kwrite ) THEN
1779         IF(lwp) write(numout,*) 'RESTART: write (XIOS 3D) ',trim(cdvar)
1780         CALL xios_send_field(trim(cdvar), pvar)
1781      ENDIF
1782#endif
1783      ELSE
1784         IF( kiomid > 0 ) THEN
1785            IF( iom_file(kiomid)%nfid > 0 ) THEN
1786               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1787               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r3d = pvar )
1788            ENDIF
1789         ENDIF
1790      ENDIF
1791   END SUBROUTINE iom_rp3d
1792
1793
1794  SUBROUTINE iom_delay_rst( cdaction, cdcpnt, kncid )
1795      !!---------------------------------------------------------------------
1796      !!   Routine iom_delay_rst: used read/write restart related to mpp_delay
1797      !!
1798      !!---------------------------------------------------------------------
1799      CHARACTER(len=*), INTENT(in   ) ::   cdaction        !
1800      CHARACTER(len=*), INTENT(in   ) ::   cdcpnt
1801      INTEGER         , INTENT(in   ) ::   kncid
1802      !
1803      INTEGER  :: ji
1804      INTEGER  :: indim
1805      LOGICAL  :: llattexist
1806      REAL(wp), ALLOCATABLE, DIMENSION(:) ::   zreal1d
1807      !!---------------------------------------------------------------------
1808      !
1809      !                                      ===================================
1810      IF( TRIM(cdaction) == 'READ' ) THEN   ! read restart related to mpp_delay !
1811         !                                   ===================================
1812         DO ji = 1, nbdelay
1813            IF ( c_delaycpnt(ji) == cdcpnt ) THEN
1814               CALL iom_chkatt( kncid, 'DELAY_'//c_delaylist(ji), llattexist, indim )
1815               IF( llattexist )  THEN
1816                  ALLOCATE( todelay(ji)%z1d(indim) )
1817                  CALL iom_getatt( kncid, 'DELAY_'//c_delaylist(ji), todelay(ji)%z1d(:) )
1818                  ndelayid(ji) = 0   ! set to 0 to specify that the value was read in the restart
1819               ENDIF
1820           ENDIF
1821         END DO
1822         !                                   ====================================
1823      ELSE                                  ! write restart related to mpp_delay !
1824         !                                   ====================================
1825         DO ji = 1, nbdelay   ! save only ocean delayed global communication variables
1826            IF ( c_delaycpnt(ji) == cdcpnt ) THEN
1827               IF( ASSOCIATED(todelay(ji)%z1d) ) THEN
1828                  CALL mpp_delay_rcv(ji)   ! make sure %z1d is received
1829                  CALL iom_putatt( kncid, 'DELAY_'//c_delaylist(ji), todelay(ji)%z1d(:) )
1830               ENDIF
1831            ENDIF
1832         END DO
1833         !
1834      ENDIF
1835     
1836   END SUBROUTINE iom_delay_rst
1837 
1838   
1839
1840   !!----------------------------------------------------------------------
1841   !!                   INTERFACE iom_put
1842   !!----------------------------------------------------------------------
1843   SUBROUTINE iom_p0d( cdname, pfield0d )
1844      CHARACTER(LEN=*), INTENT(in) ::   cdname
1845      REAL(wp)        , INTENT(in) ::   pfield0d
1846!!      REAL(wp)        , DIMENSION(jpi,jpj) ::   zz     ! masson
1847#if defined key_iomput
1848!!clem      zz(:,:)=pfield0d
1849!!clem      CALL xios_send_field(cdname, zz)
1850      CALL xios_send_field(cdname, (/pfield0d/)) 
1851#else
1852      IF( .FALSE. )   WRITE(numout,*) cdname, pfield0d   ! useless test to avoid compilation warnings
1853#endif
1854   END SUBROUTINE iom_p0d
1855
1856   SUBROUTINE iom_p1d( cdname, pfield1d )
1857      CHARACTER(LEN=*)          , INTENT(in) ::   cdname
1858      REAL(wp),     DIMENSION(:), INTENT(in) ::   pfield1d
1859#if defined key_iomput
1860      CALL xios_send_field( cdname, RESHAPE( (/pfield1d/), (/1,1,SIZE(pfield1d)/) ) )
1861#else
1862      IF( .FALSE. )   WRITE(numout,*) cdname, pfield1d   ! useless test to avoid compilation warnings
1863#endif
1864   END SUBROUTINE iom_p1d
1865
1866   SUBROUTINE iom_p2d( cdname, pfield2d )
1867      CHARACTER(LEN=*)            , INTENT(in) ::   cdname
1868      REAL(wp),     DIMENSION(:,:), INTENT(in) ::   pfield2d
1869#if defined key_iomput
1870      CALL xios_send_field(cdname, pfield2d)
1871#else
1872      IF( .FALSE. )   WRITE(numout,*) cdname, pfield2d   ! useless test to avoid compilation warnings
1873#endif
1874   END SUBROUTINE iom_p2d
1875
1876   SUBROUTINE iom_p3d( cdname, pfield3d )
1877      CHARACTER(LEN=*)                , INTENT(in) ::   cdname
1878      REAL(wp),       DIMENSION(:,:,:), INTENT(in) ::   pfield3d
1879#if defined key_iomput
1880      CALL xios_send_field( cdname, pfield3d )
1881#else
1882      IF( .FALSE. )   WRITE(numout,*) cdname, pfield3d   ! useless test to avoid compilation warnings
1883#endif
1884   END SUBROUTINE iom_p3d
1885
1886#if defined key_iomput
1887   !!----------------------------------------------------------------------
1888   !!   'key_iomput'                                         XIOS interface
1889   !!----------------------------------------------------------------------
1890
1891   SUBROUTINE iom_set_domain_attr( cdid, ni_glo, nj_glo, ibegin, jbegin, ni, nj,                                               &
1892      &                                    data_dim, data_ibegin, data_ni, data_jbegin, data_nj, lonvalue, latvalue, mask,     &
1893      &                                    nvertex, bounds_lon, bounds_lat, area )
1894      !!----------------------------------------------------------------------
1895      !!----------------------------------------------------------------------
1896      CHARACTER(LEN=*)                  , INTENT(in) ::   cdid
1897      INTEGER                 , OPTIONAL, INTENT(in) ::   ni_glo, nj_glo, ibegin, jbegin, ni, nj
1898      INTEGER                 , OPTIONAL, INTENT(in) ::   data_dim, data_ibegin, data_ni, data_jbegin, data_nj
1899      INTEGER                 , OPTIONAL, INTENT(in) ::   nvertex
1900      REAL(wp), DIMENSION(:)  , OPTIONAL, INTENT(in) ::   lonvalue, latvalue
1901      REAL(wp), DIMENSION(:,:), OPTIONAL, INTENT(in) ::   bounds_lon, bounds_lat, area
1902      LOGICAL , DIMENSION(:)  , OPTIONAL, INTENT(in) ::   mask
1903      !!----------------------------------------------------------------------
1904      !
1905      IF( xios_is_valid_domain     (cdid) ) THEN
1906         CALL xios_set_domain_attr     ( cdid, ni_glo=ni_glo, nj_glo=nj_glo, ibegin=ibegin, jbegin=jbegin, ni=ni, nj=nj,   &
1907            &    data_dim=data_dim, data_ibegin=data_ibegin, data_ni=data_ni, data_jbegin=data_jbegin, data_nj=data_nj ,   &
1908            &    lonvalue_1D=lonvalue, latvalue_1D=latvalue, mask_1D=mask, nvertex=nvertex, bounds_lon_1D=bounds_lon,      &
1909            &    bounds_lat_1D=bounds_lat, area=area, type='curvilinear')
1910      ENDIF
1911      IF( xios_is_valid_domaingroup(cdid) ) THEN
1912         CALL xios_set_domaingroup_attr( cdid, ni_glo=ni_glo, nj_glo=nj_glo, ibegin=ibegin, jbegin=jbegin, ni=ni, nj=nj,   &
1913            &    data_dim=data_dim, data_ibegin=data_ibegin, data_ni=data_ni, data_jbegin=data_jbegin, data_nj=data_nj ,   &
1914            &    lonvalue_1D=lonvalue, latvalue_1D=latvalue, mask_1D=mask, nvertex=nvertex, bounds_lon_1D=bounds_lon,      &
1915            &    bounds_lat_1D=bounds_lat, area=area, type='curvilinear' )
1916      ENDIF
1917      !
1918      CALL xios_solve_inheritance()
1919      !
1920   END SUBROUTINE iom_set_domain_attr
1921
1922
1923   SUBROUTINE iom_set_zoom_domain_attr( cdid, ibegin, jbegin, ni, nj )
1924      !!----------------------------------------------------------------------
1925      !!----------------------------------------------------------------------
1926      CHARACTER(LEN=*), INTENT(in) ::   cdid
1927      INTEGER         , INTENT(in) ::   ibegin, jbegin, ni, nj
1928      !
1929      TYPE(xios_gridgroup) :: gridgroup_hdl
1930      TYPE(xios_grid)      :: grid_hdl
1931      TYPE(xios_domain)    :: domain_hdl 
1932      TYPE(xios_axis)      :: axis_hdl 
1933      CHARACTER(LEN=64)    :: cldomrefid   ! domain_ref name
1934      CHARACTER(len=1)     :: cl1          ! last character of this name
1935      !!----------------------------------------------------------------------
1936      !
1937      IF( xios_is_valid_zoom_domain(cdid) ) THEN
1938         ! define the zoom_domain attributs
1939         CALL xios_set_zoom_domain_attr( cdid, ibegin=ibegin, jbegin=jbegin, ni=ni, nj=nj )
1940         ! define a new 2D grid with this new domain
1941         CALL xios_get_handle("grid_definition", gridgroup_hdl )
1942         CALL xios_add_child(gridgroup_hdl, grid_hdl, TRIM(cdid)//'_2D' )   ! add a new 2D grid to grid_definition
1943         CALL xios_add_child(grid_hdl, domain_hdl, TRIM(cdid) )             ! add its domain
1944         ! define a new 3D grid with this new domain
1945         CALL xios_add_child(gridgroup_hdl, grid_hdl, TRIM(cdid)//'_3D' )   ! add a new 3D grid to grid_definition
1946         CALL xios_add_child(grid_hdl, domain_hdl, TRIM(cdid) )             ! add its domain
1947         ! vertical axis
1948         cl1 = cdid(LEN_TRIM(cdid):)                                        ! last letter of cdid
1949         cl1 = CHAR(ICHAR(cl1)+32)                                          ! from upper to lower case
1950         CALL xios_add_child(grid_hdl, axis_hdl, 'depth'//cl1)              ! add its axis
1951      ENDIF
1952      !     
1953   END SUBROUTINE iom_set_zoom_domain_attr
1954
1955
1956   SUBROUTINE iom_set_axis_attr( cdid, paxis, bounds )
1957      !!----------------------------------------------------------------------
1958      !!----------------------------------------------------------------------
1959      CHARACTER(LEN=*)      , INTENT(in) ::   cdid
1960      REAL(wp), DIMENSION(:)  , OPTIONAL, INTENT(in) ::   paxis
1961      REAL(wp), DIMENSION(:,:), OPTIONAL, INTENT(in) ::   bounds
1962      !!----------------------------------------------------------------------
1963      IF( PRESENT(paxis) ) THEN
1964         IF( xios_is_valid_axis     (cdid) )   CALL xios_set_axis_attr     ( cdid, n_glo=SIZE(paxis), value=paxis )
1965         IF( xios_is_valid_axisgroup(cdid) )   CALL xios_set_axisgroup_attr( cdid, n_glo=SIZE(paxis), value=paxis )
1966      ENDIF
1967      IF( xios_is_valid_axis     (cdid) )   CALL xios_set_axis_attr     ( cdid, bounds=bounds )
1968      IF( xios_is_valid_axisgroup(cdid) )   CALL xios_set_axisgroup_attr( cdid, bounds=bounds )
1969      CALL xios_solve_inheritance()
1970   END SUBROUTINE iom_set_axis_attr
1971
1972
1973   SUBROUTINE iom_set_field_attr( cdid, freq_op, freq_offset )
1974      !!----------------------------------------------------------------------
1975      !!----------------------------------------------------------------------
1976      CHARACTER(LEN=*)             , INTENT(in) ::   cdid
1977      TYPE(xios_duration), OPTIONAL, INTENT(in) ::   freq_op
1978      TYPE(xios_duration), OPTIONAL, INTENT(in) ::   freq_offset
1979      !!----------------------------------------------------------------------
1980      IF( xios_is_valid_field     (cdid) )   CALL xios_set_field_attr     ( cdid, freq_op=freq_op, freq_offset=freq_offset )
1981      IF( xios_is_valid_fieldgroup(cdid) )   CALL xios_set_fieldgroup_attr( cdid, freq_op=freq_op, freq_offset=freq_offset )
1982      CALL xios_solve_inheritance()
1983   END SUBROUTINE iom_set_field_attr
1984
1985
1986   SUBROUTINE iom_set_file_attr( cdid, name, name_suffix )
1987      !!----------------------------------------------------------------------
1988      !!----------------------------------------------------------------------
1989      CHARACTER(LEN=*)          , INTENT(in) ::   cdid
1990      CHARACTER(LEN=*),OPTIONAL , INTENT(in) ::   name, name_suffix
1991      !!----------------------------------------------------------------------
1992      IF( xios_is_valid_file     (cdid) )   CALL xios_set_file_attr     ( cdid, name=name, name_suffix=name_suffix )
1993      IF( xios_is_valid_filegroup(cdid) )   CALL xios_set_filegroup_attr( cdid, name=name, name_suffix=name_suffix )
1994      CALL xios_solve_inheritance()
1995   END SUBROUTINE iom_set_file_attr
1996
1997
1998   SUBROUTINE iom_get_file_attr( cdid, name, name_suffix, output_freq )
1999      !!----------------------------------------------------------------------
2000      !!----------------------------------------------------------------------
2001      CHARACTER(LEN=*)          , INTENT(in ) ::   cdid
2002      CHARACTER(LEN=*),OPTIONAL , INTENT(out) ::   name, name_suffix
2003      TYPE(xios_duration), OPTIONAL , INTENT(out) :: output_freq
2004      LOGICAL                                 ::   llexist1,llexist2,llexist3
2005      !---------------------------------------------------------------------
2006      IF( PRESENT( name        ) )   name = ''          ! default values
2007      IF( PRESENT( name_suffix ) )   name_suffix = ''
2008      IF( PRESENT( output_freq ) )   output_freq = xios_duration(0,0,0,0,0,0)
2009      IF( xios_is_valid_file     (cdid) ) THEN
2010         CALL xios_solve_inheritance()
2011         CALL xios_is_defined_file_attr     ( cdid, name = llexist1, name_suffix = llexist2, output_freq = llexist3)
2012         IF(llexist1)   CALL xios_get_file_attr     ( cdid, name = name )
2013         IF(llexist2)   CALL xios_get_file_attr     ( cdid, name_suffix = name_suffix )
2014         IF(llexist3)   CALL xios_get_file_attr     ( cdid, output_freq = output_freq )
2015      ENDIF
2016      IF( xios_is_valid_filegroup(cdid) ) THEN
2017         CALL xios_solve_inheritance()
2018         CALL xios_is_defined_filegroup_attr( cdid, name = llexist1, name_suffix = llexist2, output_freq = llexist3)
2019         IF(llexist1)   CALL xios_get_filegroup_attr( cdid, name = name )
2020         IF(llexist2)   CALL xios_get_filegroup_attr( cdid, name_suffix = name_suffix )
2021         IF(llexist3)   CALL xios_get_filegroup_attr( cdid, output_freq = output_freq )
2022      ENDIF
2023   END SUBROUTINE iom_get_file_attr
2024
2025
2026   SUBROUTINE iom_set_grid_attr( cdid, mask )
2027      !!----------------------------------------------------------------------
2028      !!----------------------------------------------------------------------
2029      CHARACTER(LEN=*)                   , INTENT(in) ::   cdid
2030      LOGICAL, DIMENSION(:,:,:), OPTIONAL, INTENT(in) ::   mask
2031      !!----------------------------------------------------------------------
2032      IF( xios_is_valid_grid     (cdid) )   CALL xios_set_grid_attr     ( cdid, mask_3D=mask )
2033      IF( xios_is_valid_gridgroup(cdid) )   CALL xios_set_gridgroup_attr( cdid, mask_3D=mask )
2034      CALL xios_solve_inheritance()
2035   END SUBROUTINE iom_set_grid_attr
2036
2037   SUBROUTINE iom_setkt( kt, cdname )
2038      !!----------------------------------------------------------------------
2039      !!----------------------------------------------------------------------
2040      INTEGER         , INTENT(in) ::   kt 
2041      CHARACTER(LEN=*), INTENT(in) ::   cdname
2042      !!----------------------------------------------------------------------
2043      CALL iom_swap( cdname )   ! swap to cdname context
2044      CALL xios_update_calendar(kt)
2045      IF( cdname /= TRIM(cxios_context) )   CALL iom_swap( TRIM(cxios_context) )   ! return back to nemo context
2046   END SUBROUTINE iom_setkt
2047
2048   SUBROUTINE iom_context_finalize( cdname )
2049      !!----------------------------------------------------------------------
2050      !!----------------------------------------------------------------------
2051      CHARACTER(LEN=*), INTENT(in) :: cdname
2052      CHARACTER(LEN=120)           :: clname
2053      !!----------------------------------------------------------------------
2054      clname = cdname
2055      IF( TRIM(Agrif_CFixed()) .NE. '0' ) clname = TRIM(Agrif_CFixed())//"_"//clname 
2056      IF( xios_is_valid_context(clname) ) THEN
2057         CALL iom_swap( cdname )   ! swap to cdname context
2058         CALL xios_context_finalize() ! finalize the context
2059         IF( cdname /= TRIM(cxios_context) ) CALL iom_swap( TRIM(cxios_context) )   ! return back to nemo context
2060      ENDIF
2061      !
2062   END SUBROUTINE iom_context_finalize
2063
2064
2065   SUBROUTINE set_grid( cdgrd, plon, plat, ldxios, ldrxios )
2066      !!----------------------------------------------------------------------
2067      !!                     ***  ROUTINE set_grid  ***
2068      !!
2069      !! ** Purpose :   define horizontal grids
2070      !!----------------------------------------------------------------------
2071      CHARACTER(LEN=1)            , INTENT(in) ::   cdgrd
2072      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   plon
2073      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   plat
2074      !
2075      INTEGER  :: ni, nj
2076      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zmask
2077      LOGICAL, INTENT(IN) :: ldxios, ldrxios
2078      !!----------------------------------------------------------------------
2079      !
2080      ni = nlei-nldi+1
2081      nj = nlej-nldj+1
2082      !
2083      CALL iom_set_domain_attr("grid_"//cdgrd, ni_glo=jpiglo, nj_glo=jpjglo, ibegin=nimpp+nldi-2, jbegin=njmpp+nldj-2, ni=ni, nj=nj)
2084      CALL iom_set_domain_attr("grid_"//cdgrd, data_dim=2, data_ibegin = 1-nldi, data_ni = jpi, data_jbegin = 1-nldj, data_nj = jpj)
2085!don't define lon and lat for restart reading context.
2086      IF ( .NOT.ldrxios ) &
2087         CALL iom_set_domain_attr("grid_"//cdgrd, lonvalue = RESHAPE(plon(nldi:nlei, nldj:nlej),(/ ni*nj /)),   &
2088         &                                     latvalue = RESHAPE(plat(nldi:nlei, nldj:nlej),(/ ni*nj /))) 
2089      !
2090      IF ( ln_mskland .AND. (.NOT.ldxios) ) THEN
2091         ! mask land points, keep values on coast line -> specific mask for U, V and W points
2092         SELECT CASE ( cdgrd )
2093         CASE('T')   ;   zmask(:,:,:)       = tmask(:,:,:)
2094         CASE('U')   ;   zmask(2:jpim1,:,:) = tmask(2:jpim1,:,:) + tmask(3:jpi,:,:)   ;   CALL lbc_lnk( 'iom', zmask, 'U', 1. )
2095         CASE('V')   ;   zmask(:,2:jpjm1,:) = tmask(:,2:jpjm1,:) + tmask(:,3:jpj,:)   ;   CALL lbc_lnk( 'iom', zmask, 'V', 1. )
2096         CASE('W')   ;   zmask(:,:,2:jpk  ) = tmask(:,:,1:jpkm1) + tmask(:,:,2:jpk)   ;   zmask(:,:,1) = tmask(:,:,1)
2097         END SELECT
2098         !
2099         CALL iom_set_domain_attr( "grid_"//cdgrd       , mask = RESHAPE(zmask(nldi:nlei,nldj:nlej,1),(/ni*nj    /)) /= 0. )
2100         CALL iom_set_grid_attr  ( "grid_"//cdgrd//"_3D", mask = RESHAPE(zmask(nldi:nlei,nldj:nlej,:),(/ni,nj,jpk/)) /= 0. )
2101      ENDIF
2102      !
2103   END SUBROUTINE set_grid
2104
2105
2106   SUBROUTINE set_grid_bounds( cdgrd, plon_cnr, plat_cnr, plon_pnt, plat_pnt )
2107      !!----------------------------------------------------------------------
2108      !!                   ***  ROUTINE set_grid_bounds  ***
2109      !!
2110      !! ** Purpose :   define horizontal grid corners
2111      !!
2112      !!----------------------------------------------------------------------
2113      CHARACTER(LEN=1)                      , INTENT(in) :: cdgrd
2114      REAL(wp), DIMENSION(jpi,jpj)          , INTENT(in) :: plon_cnr, plat_cnr  ! Lat/lon coord. of a contiguous vertex of cell (i,j)
2115      REAL(wp), DIMENSION(jpi,jpj), OPTIONAL, INTENT(in) :: plon_pnt, plat_pnt  ! Lat/lon coord. of the point of cell (i,j)
2116      !
2117      INTEGER :: ji, jj, jn, ni, nj
2118      INTEGER :: icnr, jcnr                                    ! Offset such that the vertex coordinate (i+icnr,j+jcnr)
2119      !                                                        ! represents the bottom-left corner of cell (i,j)
2120      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: z_bnds      ! Lat/lon coordinates of the vertices of cell (i,j)
2121      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z_fld       ! Working array to determine where to rotate cells
2122      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z_rot       ! Lat/lon working array for rotation of cells
2123      !!----------------------------------------------------------------------
2124      !
2125      ALLOCATE( z_bnds(4,jpi,jpj,2), z_fld(jpi,jpj), z_rot(4,2)  )
2126      !
2127      ! Offset of coordinate representing bottom-left corner
2128      SELECT CASE ( TRIM(cdgrd) )
2129      CASE ('T', 'W')   ;   icnr = -1   ;   jcnr = -1
2130      CASE ('U')        ;   icnr =  0   ;   jcnr = -1
2131      CASE ('V')        ;   icnr = -1   ;   jcnr =  0
2132      END SELECT
2133      !
2134      ni = nlei-nldi+1   ! Dimensions of subdomain interior
2135      nj = nlej-nldj+1
2136      !
2137      z_fld(:,:) = 1._wp
2138      CALL lbc_lnk( 'iom', z_fld, cdgrd, -1. )    ! Working array for location of northfold
2139      !
2140      ! Cell vertices that can be defined
2141      DO jj = 2, jpjm1
2142         DO ji = 2, jpim1
2143            z_bnds(1,ji,jj,1) = plat_cnr(ji+icnr,  jj+jcnr  ) ! Bottom-left
2144            z_bnds(2,ji,jj,1) = plat_cnr(ji+icnr+1,jj+jcnr  ) ! Bottom-right
2145            z_bnds(3,ji,jj,1) = plat_cnr(ji+icnr+1,jj+jcnr+1) ! Top-right
2146            z_bnds(4,ji,jj,1) = plat_cnr(ji+icnr,  jj+jcnr+1) ! Top-left
2147            z_bnds(1,ji,jj,2) = plon_cnr(ji+icnr,  jj+jcnr  ) ! Bottom-left
2148            z_bnds(2,ji,jj,2) = plon_cnr(ji+icnr+1,jj+jcnr  ) ! Bottom-right
2149            z_bnds(3,ji,jj,2) = plon_cnr(ji+icnr+1,jj+jcnr+1) ! Top-right
2150            z_bnds(4,ji,jj,2) = plon_cnr(ji+icnr,  jj+jcnr+1) ! Top-left
2151         END DO
2152      END DO
2153      !
2154      ! Cell vertices on boundries
2155      DO jn = 1, 4
2156         CALL lbc_lnk( 'iom', z_bnds(jn,:,:,1), cdgrd, 1., pfillval=999._wp )
2157         CALL lbc_lnk( 'iom', z_bnds(jn,:,:,2), cdgrd, 1., pfillval=999._wp )
2158      END DO
2159      !
2160      ! Zero-size cells at closed boundaries if cell points provided,
2161      ! otherwise they are closed cells with unrealistic bounds
2162      IF( PRESENT(plon_pnt) .AND. PRESENT(plat_pnt) ) THEN
2163         IF( (nbondi == -1 .OR. nbondi == 2) .AND. .NOT. (jperio == 1 .OR. jperio == 4 .OR. jperio == 6) ) THEN
2164            DO jn = 1, 4        ! (West or jpni = 1), closed E-W
2165               z_bnds(jn,1,:,1) = plat_pnt(1,:)  ;  z_bnds(jn,1,:,2) = plon_pnt(1,:)
2166            END DO
2167         ENDIF
2168         IF( (nbondi == 1 .OR. nbondi == 2) .AND. .NOT. (jperio == 1 .OR. jperio == 4 .OR. jperio == 6) ) THEN
2169            DO jn = 1, 4        ! (East or jpni = 1), closed E-W
2170               z_bnds(jn,nlci,:,1) = plat_pnt(nlci,:)  ;  z_bnds(jn,nlci,:,2) = plon_pnt(nlci,:)
2171            END DO
2172         ENDIF
2173         IF( nbondj == -1 .OR. (nbondj == 2 .AND. jperio /= 2) ) THEN
2174            DO jn = 1, 4        ! South or (jpnj = 1, not symmetric)
2175               z_bnds(jn,:,1,1) = plat_pnt(:,1)  ;  z_bnds(jn,:,1,2) = plon_pnt(:,1)
2176            END DO
2177         ENDIF
2178         IF( (nbondj == 1 .OR. nbondj == 2) .AND. jperio  < 3 ) THEN
2179            DO jn = 1, 4        ! (North or jpnj = 1), no north fold
2180               z_bnds(jn,:,nlcj,1) = plat_pnt(:,nlcj)  ;  z_bnds(jn,:,nlcj,2) = plon_pnt(:,nlcj)
2181            END DO
2182         ENDIF
2183      ENDIF
2184      !
2185      IF( (nbondj == 1 .OR. nbondj == 2) .AND. jperio >= 3 ) THEN    ! Rotate cells at the north fold
2186         DO jj = 1, jpj
2187            DO ji = 1, jpi
2188               IF( z_fld(ji,jj) == -1. ) THEN
2189                  z_rot(1,:) = z_bnds(3,ji,jj,:) ; z_rot(2,:) = z_bnds(4,ji,jj,:)
2190                  z_rot(3,:) = z_bnds(1,ji,jj,:) ; z_rot(4,:) = z_bnds(2,ji,jj,:)
2191                  z_bnds(:,ji,jj,:) = z_rot(:,:)
2192               ENDIF
2193            END DO
2194         END DO
2195      ELSE IF( nbondj == 2 .AND. jperio == 2 ) THEN                  ! Invert cells at the symmetric equator
2196         DO ji = 1, jpi
2197            z_rot(1:2,:) = z_bnds(3:4,ji,1,:)
2198            z_rot(3:4,:) = z_bnds(1:2,ji,1,:)
2199            z_bnds(:,ji,1,:) = z_rot(:,:)
2200         END DO
2201      ENDIF
2202      !
2203      CALL iom_set_domain_attr("grid_"//cdgrd, bounds_lat = RESHAPE(z_bnds(:,nldi:nlei,nldj:nlej,1),(/ 4,ni*nj /)),           &
2204          &                                    bounds_lon = RESHAPE(z_bnds(:,nldi:nlei,nldj:nlej,2),(/ 4,ni*nj /)), nvertex=4 )
2205      !
2206      DEALLOCATE( z_bnds, z_fld, z_rot ) 
2207      !
2208   END SUBROUTINE set_grid_bounds
2209
2210
2211   SUBROUTINE set_grid_znl( plat )
2212      !!----------------------------------------------------------------------
2213      !!                     ***  ROUTINE set_grid_znl  ***
2214      !!
2215      !! ** Purpose :   define grids for zonal mean
2216      !!
2217      !!----------------------------------------------------------------------
2218      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   plat
2219      !
2220      INTEGER  :: ni, nj, ix, iy
2221      REAL(wp), DIMENSION(:), ALLOCATABLE  ::   zlon
2222      !!----------------------------------------------------------------------
2223      !
2224      ni=nlei-nldi+1       ! define zonal mean domain (jpj*jpk)
2225      nj=nlej-nldj+1
2226      ALLOCATE( zlon(ni*nj) )       ;       zlon(:) = 0._wp
2227      !
2228      CALL dom_ngb( -168.53, 65.03, ix, iy, 'T' ) !  i-line that passes through Bering Strait: Reference latitude (used in plots)
2229!      CALL dom_ngb( 180., 90., ix, iy, 'T' ) !  i-line that passes near the North Pole : Reference latitude (used in plots)
2230      CALL iom_set_domain_attr("gznl", ni_glo=jpiglo, nj_glo=jpjglo, ibegin=nimpp+nldi-2, jbegin=njmpp+nldj-2, ni=ni, nj=nj)
2231      CALL iom_set_domain_attr("gznl", data_dim=2, data_ibegin = 1-nldi, data_ni = jpi, data_jbegin = 1-nldj, data_nj = jpj)
2232      CALL iom_set_domain_attr("gznl", lonvalue = zlon,   &
2233         &                             latvalue = RESHAPE(plat(nldi:nlei, nldj:nlej),(/ ni*nj /))) 
2234      CALL iom_set_zoom_domain_attr("znl_T", ibegin=ix-1, jbegin=0, ni=1, nj=jpjglo)
2235      CALL iom_set_zoom_domain_attr("znl_W", ibegin=ix-1, jbegin=0, ni=1, nj=jpjglo)
2236      !
2237      CALL iom_update_file_name('ptr')
2238      !
2239   END SUBROUTINE set_grid_znl
2240
2241
2242   SUBROUTINE set_scalar
2243      !!----------------------------------------------------------------------
2244      !!                     ***  ROUTINE set_scalar  ***
2245      !!
2246      !! ** Purpose :   define fake grids for scalar point
2247      !!
2248      !!----------------------------------------------------------------------
2249      REAL(wp), DIMENSION(1)   ::   zz = 1.
2250      !!----------------------------------------------------------------------
2251      !
2252      CALL iom_set_domain_attr('scalarpoint', ni_glo=jpnij, nj_glo=1, ibegin=narea-1, jbegin=0, ni=1, nj=1)
2253      CALL iom_set_domain_attr('scalarpoint', data_dim=2, data_ibegin = 1, data_ni = 1, data_jbegin = 1, data_nj = 1)
2254      !
2255      zz = REAL( narea, wp )
2256      CALL iom_set_domain_attr('scalarpoint', lonvalue=zz, latvalue=zz)
2257      !
2258   END SUBROUTINE set_scalar
2259
2260
2261   SUBROUTINE set_xmlatt
2262      !!----------------------------------------------------------------------
2263      !!                     ***  ROUTINE set_xmlatt  ***
2264      !!
2265      !! ** Purpose :   automatic definitions of some of the xml attributs...
2266      !!
2267      !!----------------------------------------------------------------------
2268      CHARACTER(len=1),DIMENSION( 3) ::   clgrd                    ! suffix name
2269      CHARACTER(len=256)             ::   clsuff                   ! suffix name
2270      CHARACTER(len=1)               ::   cl1                      ! 1 character
2271      CHARACTER(len=2)               ::   cl2                      ! 2 characters
2272      CHARACTER(len=3)               ::   cl3                      ! 3 characters
2273      INTEGER                        ::   ji, jg                   ! loop counters
2274      INTEGER                        ::   ix, iy                   ! i-,j- index
2275      REAL(wp)        ,DIMENSION(11) ::   zlontao                  ! longitudes of tao    moorings
2276      REAL(wp)        ,DIMENSION( 7) ::   zlattao                  ! latitudes  of tao    moorings
2277      REAL(wp)        ,DIMENSION( 4) ::   zlonrama                 ! longitudes of rama   moorings
2278      REAL(wp)        ,DIMENSION(11) ::   zlatrama                 ! latitudes  of rama   moorings
2279      REAL(wp)        ,DIMENSION( 3) ::   zlonpira                 ! longitudes of pirata moorings
2280      REAL(wp)        ,DIMENSION( 9) ::   zlatpira                 ! latitudes  of pirata moorings
2281      TYPE(xios_duration)            ::   f_op, f_of
2282      !!----------------------------------------------------------------------
2283      !
2284      ! frequency of the call of iom_put (attribut: freq_op)
2285      f_op%timestep = 1        ;  f_of%timestep =  0  ; CALL iom_set_field_attr('field_definition', freq_op=f_op, freq_offset=f_of)
2286      f_op%timestep = 2        ;  f_of%timestep =  0  ; CALL iom_set_field_attr('trendT_even'     , freq_op=f_op, freq_offset=f_of)
2287      f_op%timestep = 2        ;  f_of%timestep = -1  ; CALL iom_set_field_attr('trendT_odd'      , freq_op=f_op, freq_offset=f_of)
2288      f_op%timestep = nn_fsbc  ;  f_of%timestep =  0  ; CALL iom_set_field_attr('SBC'             , freq_op=f_op, freq_offset=f_of)
2289      f_op%timestep = nn_fsbc  ;  f_of%timestep =  0  ; CALL iom_set_field_attr('SBC_scalar'      , freq_op=f_op, freq_offset=f_of)
2290      f_op%timestep = nn_dttrc ;  f_of%timestep =  0  ; CALL iom_set_field_attr('ptrc_T'          , freq_op=f_op, freq_offset=f_of)
2291      f_op%timestep = nn_dttrc ;  f_of%timestep =  0  ; CALL iom_set_field_attr('diad_T'          , freq_op=f_op, freq_offset=f_of)
2292
2293      ! output file names (attribut: name)
2294      DO ji = 1, 9
2295         WRITE(cl1,'(i1)') ji 
2296         CALL iom_update_file_name('file'//cl1)
2297      END DO
2298      DO ji = 1, 99
2299         WRITE(cl2,'(i2.2)') ji 
2300         CALL iom_update_file_name('file'//cl2)
2301      END DO
2302      DO ji = 1, 999
2303         WRITE(cl3,'(i3.3)') ji 
2304         CALL iom_update_file_name('file'//cl3)
2305      END DO
2306
2307      ! Zooms...
2308      clgrd = (/ 'T', 'U', 'W' /) 
2309      DO jg = 1, SIZE(clgrd)                                                                   ! grid type
2310         cl1 = clgrd(jg)
2311         ! Equatorial section (attributs: jbegin, ni, name_suffix)
2312         CALL dom_ngb( 0., 0., ix, iy, cl1 )
2313         CALL iom_set_zoom_domain_attr('Eq'//cl1, ibegin=0, jbegin=iy-1, ni=jpiglo, nj=1 )
2314         CALL iom_get_file_attr   ('Eq'//cl1, name_suffix = clsuff             )
2315         CALL iom_set_file_attr   ('Eq'//cl1, name_suffix = TRIM(clsuff)//'_Eq')
2316         CALL iom_update_file_name('Eq'//cl1)
2317      END DO
2318      ! TAO moorings (attributs: ibegin, jbegin, name_suffix)
2319      zlontao = (/ 137.0, 147.0, 156.0, 165.0, -180.0, -170.0, -155.0, -140.0, -125.0, -110.0, -95.0 /)
2320      zlattao = (/  -8.0,  -5.0,  -2.0,   0.0,    2.0,    5.0,    8.0 /)
2321      CALL set_mooring( zlontao, zlattao )
2322      ! RAMA moorings (attributs: ibegin, jbegin, name_suffix)
2323      zlonrama = (/  55.0,  67.0, 80.5, 90.0 /)
2324      zlatrama = (/ -16.0, -12.0, -8.0, -4.0, -1.5, 0.0, 1.5, 4.0, 8.0, 12.0, 15.0 /)
2325      CALL set_mooring( zlonrama, zlatrama )
2326      ! PIRATA moorings (attributs: ibegin, jbegin, name_suffix)
2327      zlonpira = (/ -38.0, -23.0, -10.0 /)
2328      zlatpira = (/ -19.0, -14.0,  -8.0, 0.0, 4.0, 8.0, 12.0, 15.0, 20.0 /)
2329      CALL set_mooring( zlonpira, zlatpira )
2330      !
2331   END SUBROUTINE set_xmlatt
2332
2333
2334   SUBROUTINE set_mooring( plon, plat )
2335      !!----------------------------------------------------------------------
2336      !!                     ***  ROUTINE set_mooring  ***
2337      !!
2338      !! ** Purpose :   automatic definitions of moorings xml attributs...
2339      !!
2340      !!----------------------------------------------------------------------
2341      REAL(wp), DIMENSION(:), INTENT(in) ::   plon, plat   ! longitudes/latitudes oft the mooring
2342      !
2343!!$      CHARACTER(len=1),DIMENSION(4) ::   clgrd = (/ 'T', 'U', 'V', 'W' /)   ! suffix name
2344      CHARACTER(len=1),DIMENSION(1) ::   clgrd = (/ 'T' /)        ! suffix name
2345      CHARACTER(len=256)            ::   clname                   ! file name
2346      CHARACTER(len=256)            ::   clsuff                   ! suffix name
2347      CHARACTER(len=1)              ::   cl1                      ! 1 character
2348      CHARACTER(len=6)              ::   clon,clat                ! name of longitude, latitude
2349      INTEGER                       ::   ji, jj, jg               ! loop counters
2350      INTEGER                       ::   ix, iy                   ! i-,j- index
2351      REAL(wp)                      ::   zlon, zlat
2352      !!----------------------------------------------------------------------
2353      DO jg = 1, SIZE(clgrd)
2354         cl1 = clgrd(jg)
2355         DO ji = 1, SIZE(plon)
2356            DO jj = 1, SIZE(plat)
2357               zlon = plon(ji)
2358               zlat = plat(jj)
2359               ! modifications for RAMA moorings
2360               IF( zlon ==  67. .AND. zlat ==  15. )   zlon =  65.
2361               IF( zlon ==  90. .AND. zlat <=  -4. )   zlon =  95.
2362               IF( zlon ==  95. .AND. zlat ==  -4. )   zlat =  -5.
2363               ! modifications for PIRATA moorings
2364               IF( zlon == -38. .AND. zlat == -19. )   zlon = -34.
2365               IF( zlon == -38. .AND. zlat == -14. )   zlon = -32.
2366               IF( zlon == -38. .AND. zlat ==  -8. )   zlon = -30.
2367               IF( zlon == -38. .AND. zlat ==   0. )   zlon = -35.
2368               IF( zlon == -23. .AND. zlat ==  20. )   zlat =  21.
2369               IF( zlon == -10. .AND. zlat == -14. )   zlat = -10.
2370               IF( zlon == -10. .AND. zlat ==  -8. )   zlat =  -6.
2371               IF( zlon == -10. .AND. zlat ==   4. ) THEN   ;   zlon = 0.   ;   zlat = 0.   ;   ENDIF
2372               CALL dom_ngb( zlon, zlat, ix, iy, cl1 )
2373               IF( zlon >= 0. ) THEN 
2374                  IF( zlon == REAL(NINT(zlon), wp) ) THEN   ;   WRITE(clon, '(i3,  a)') NINT( zlon), 'e'
2375                  ELSE                                      ;   WRITE(clon, '(f5.1,a)')       zlon , 'e'
2376                  ENDIF
2377               ELSE             
2378                  IF( zlon == REAL(NINT(zlon), wp) ) THEN   ;   WRITE(clon, '(i3,  a)') NINT(-zlon), 'w'
2379                  ELSE                                      ;   WRITE(clon, '(f5.1,a)')      -zlon , 'w'
2380                  ENDIF
2381               ENDIF
2382               IF( zlat >= 0. ) THEN 
2383                  IF( zlat == REAL(NINT(zlat), wp) ) THEN   ;   WRITE(clat, '(i2,  a)') NINT( zlat), 'n'
2384                  ELSE                                      ;   WRITE(clat, '(f4.1,a)')       zlat , 'n'
2385                  ENDIF
2386               ELSE             
2387                  IF( zlat == REAL(NINT(zlat), wp) ) THEN   ;   WRITE(clat, '(i2,  a)') NINT(-zlat), 's'
2388                  ELSE                                      ;   WRITE(clat, '(f4.1,a)')      -zlat , 's'
2389                  ENDIF
2390               ENDIF
2391               clname = TRIM(ADJUSTL(clat))//TRIM(ADJUSTL(clon))
2392               CALL iom_set_zoom_domain_attr(TRIM(clname)//cl1, ibegin= ix-1, jbegin= iy-1, ni=1, nj=1)
2393
2394               CALL iom_get_file_attr   (TRIM(clname)//cl1, name_suffix = clsuff                         )
2395               CALL iom_set_file_attr   (TRIM(clname)//cl1, name_suffix = TRIM(clsuff)//'_'//TRIM(clname))
2396               CALL iom_update_file_name(TRIM(clname)//cl1)
2397            END DO
2398         END DO
2399      END DO
2400     
2401   END SUBROUTINE set_mooring
2402
2403   
2404   SUBROUTINE iom_update_file_name( cdid )
2405      !!----------------------------------------------------------------------
2406      !!                     ***  ROUTINE iom_update_file_name  ***
2407      !!
2408      !! ** Purpose :   
2409      !!
2410      !!----------------------------------------------------------------------
2411      CHARACTER(LEN=*)          , INTENT(in) ::   cdid
2412      !
2413      CHARACTER(LEN=256) ::   clname
2414      CHARACTER(LEN=20)  ::   clfreq
2415      CHARACTER(LEN=20)  ::   cldate
2416      CHARACTER(LEN=256) ::   cltmpn                 !FUS needed for correct path with AGRIF
2417      INTEGER            ::   iln                    !FUS needed for correct path with AGRIF
2418      INTEGER            ::   idx
2419      INTEGER            ::   jn
2420      INTEGER            ::   itrlen
2421      INTEGER            ::   iyear, imonth, iday, isec
2422      REAL(wp)           ::   zsec
2423      LOGICAL            ::   llexist
2424      TYPE(xios_duration)   ::   output_freq 
2425      !!----------------------------------------------------------------------
2426      !
2427      DO jn = 1, 2
2428         !
2429         output_freq = xios_duration(0,0,0,0,0,0)
2430         IF( jn == 1 )   CALL iom_get_file_attr( cdid, name        = clname, output_freq = output_freq )
2431         IF( jn == 2 )   CALL iom_get_file_attr( cdid, name_suffix = clname )
2432         !
2433         IF ( TRIM(clname) /= '' ) THEN 
2434            !
2435            idx = INDEX(clname,'@expname@') + INDEX(clname,'@EXPNAME@')
2436            DO WHILE ( idx /= 0 ) 
2437               clname = clname(1:idx-1)//TRIM(cexper)//clname(idx+9:LEN_TRIM(clname))
2438               idx = INDEX(clname,'@expname@') + INDEX(clname,'@EXPNAME@')
2439            END DO
2440            !
2441            idx = INDEX(clname,'@freq@') + INDEX(clname,'@FREQ@')
2442            DO WHILE ( idx /= 0 ) 
2443              IF ( output_freq%timestep /= 0) THEN
2444                  WRITE(clfreq,'(I18,A2)')INT(output_freq%timestep),'ts' 
2445                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2446              ELSE IF ( output_freq%second /= 0 ) THEN
2447                  WRITE(clfreq,'(I19,A1)')INT(output_freq%second),'s' 
2448                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2449              ELSE IF ( output_freq%minute /= 0 ) THEN
2450                  WRITE(clfreq,'(I18,A2)')INT(output_freq%minute),'mi' 
2451                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2452              ELSE IF ( output_freq%hour /= 0 ) THEN
2453                  WRITE(clfreq,'(I19,A1)')INT(output_freq%hour),'h' 
2454                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2455              ELSE IF ( output_freq%day /= 0 ) THEN
2456                  WRITE(clfreq,'(I19,A1)')INT(output_freq%day),'d' 
2457                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2458              ELSE IF ( output_freq%month /= 0 ) THEN   
2459                  WRITE(clfreq,'(I19,A1)')INT(output_freq%month),'m' 
2460                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2461              ELSE IF ( output_freq%year /= 0 ) THEN   
2462                  WRITE(clfreq,'(I19,A1)')INT(output_freq%year),'y' 
2463                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2464              ELSE
2465                  CALL ctl_stop('error in the name of file id '//TRIM(cdid),   &
2466                     & ' attribute output_freq is undefined -> cannot replace @freq@ in '//TRIM(clname) )
2467              ENDIF
2468              clname = clname(1:idx-1)//TRIM(ADJUSTL(clfreq))//clname(idx+6:LEN_TRIM(clname))
2469              idx = INDEX(clname,'@freq@') + INDEX(clname,'@FREQ@')
2470            END DO
2471            !
2472            idx = INDEX(clname,'@startdate@') + INDEX(clname,'@STARTDATE@')
2473            DO WHILE ( idx /= 0 ) 
2474               cldate = iom_sdate( fjulday - rdt / rday )
2475               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+11:LEN_TRIM(clname))
2476               idx = INDEX(clname,'@startdate@') + INDEX(clname,'@STARTDATE@')
2477            END DO
2478            !
2479            idx = INDEX(clname,'@startdatefull@') + INDEX(clname,'@STARTDATEFULL@')
2480            DO WHILE ( idx /= 0 ) 
2481               cldate = iom_sdate( fjulday - rdt / rday, ldfull = .TRUE. )
2482               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+15:LEN_TRIM(clname))
2483               idx = INDEX(clname,'@startdatefull@') + INDEX(clname,'@STARTDATEFULL@')
2484            END DO
2485            !
2486            idx = INDEX(clname,'@enddate@') + INDEX(clname,'@ENDDATE@')
2487            DO WHILE ( idx /= 0 ) 
2488               cldate = iom_sdate( fjulday + rdt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE. )
2489               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+9:LEN_TRIM(clname))
2490               idx = INDEX(clname,'@enddate@') + INDEX(clname,'@ENDDATE@')
2491            END DO
2492            !
2493            idx = INDEX(clname,'@enddatefull@') + INDEX(clname,'@ENDDATEFULL@')
2494            DO WHILE ( idx /= 0 ) 
2495               cldate = iom_sdate( fjulday + rdt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE., ldfull = .TRUE. )
2496               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+13:LEN_TRIM(clname))
2497               idx = INDEX(clname,'@enddatefull@') + INDEX(clname,'@ENDDATEFULL@')
2498            END DO
2499            !
2500!FUS            IF( jn == 1 .AND. TRIM(Agrif_CFixed()) /= '0' )   clname = TRIM(Agrif_CFixed())//"_"//TRIM(clname)
2501!FUS see comment line 700
2502            IF( jn == 1 .AND. TRIM(Agrif_CFixed()) /= '0' ) THEN
2503             iln    = INDEX(clname,'/',BACK=.true.)
2504             cltmpn = clname(1:iln)
2505             clname = clname(iln+1:LEN_TRIM(clname))
2506             clname = TRIM(cltmpn)//TRIM(Agrif_CFixed())//'_'//TRIM(clname)
2507            ENDIF
2508!FUS
2509            IF( jn == 1 )   CALL iom_set_file_attr( cdid, name        = clname )
2510            IF( jn == 2 )   CALL iom_set_file_attr( cdid, name_suffix = clname )
2511            !
2512         ENDIF
2513         !
2514      END DO
2515      !
2516   END SUBROUTINE iom_update_file_name
2517
2518
2519   FUNCTION iom_sdate( pjday, ld24, ldfull )
2520      !!----------------------------------------------------------------------
2521      !!                     ***  ROUTINE iom_sdate  ***
2522      !!
2523      !! ** Purpose :   send back the date corresponding to the given julian day
2524      !!----------------------------------------------------------------------
2525      REAL(wp), INTENT(in   )           ::   pjday    ! julian day
2526      LOGICAL , INTENT(in   ), OPTIONAL ::   ld24     ! true to force 24:00 instead of 00:00
2527      LOGICAL , INTENT(in   ), OPTIONAL ::   ldfull   ! true to get the compleate date: yyyymmdd_hh:mm:ss
2528      !
2529      CHARACTER(LEN=20) ::   iom_sdate
2530      CHARACTER(LEN=50) ::   clfmt                         !  format used to write the date
2531      INTEGER           ::   iyear, imonth, iday, ihour, iminute, isec
2532      REAL(wp)          ::   zsec
2533      LOGICAL           ::   ll24, llfull
2534      !!----------------------------------------------------------------------
2535      !
2536      IF( PRESENT(ld24) ) THEN   ;   ll24 = ld24
2537      ELSE                       ;   ll24 = .FALSE.
2538      ENDIF
2539      !
2540      IF( PRESENT(ldfull) ) THEN   ;   llfull = ldfull
2541      ELSE                         ;   llfull = .FALSE.
2542      ENDIF
2543      !
2544      CALL ju2ymds( pjday, iyear, imonth, iday, zsec )
2545      isec = NINT(zsec)
2546      !
2547      IF ( ll24 .AND. isec == 0 ) THEN   ! 00:00 of the next day -> move to 24:00 of the current day
2548         CALL ju2ymds( pjday - 1., iyear, imonth, iday, zsec )
2549         isec = 86400
2550      ENDIF
2551      !
2552      IF( iyear < 10000 ) THEN   ;   clfmt = "i4.4,2i2.2"                ! format used to write the date
2553      ELSE                       ;   WRITE(clfmt, "('i',i1,',2i2.2')") INT(LOG10(REAL(iyear,wp))) + 1
2554      ENDIF
2555      !
2556!$AGRIF_DO_NOT_TREAT     
2557      ! needed in the conv
2558      IF( llfull ) THEN
2559         clfmt = TRIM(clfmt)//",'_',i2.2,':',i2.2,':',i2.2"
2560         ihour   = isec / 3600
2561         isec    = MOD(isec, 3600)
2562         iminute = isec / 60
2563         isec    = MOD(isec, 60)
2564         WRITE(iom_sdate, '('//TRIM(clfmt)//')') iyear, imonth, iday, ihour, iminute, isec    ! date of the end of run
2565      ELSE
2566         WRITE(iom_sdate, '('//TRIM(clfmt)//')') iyear, imonth, iday                          ! date of the end of run
2567      ENDIF
2568!$AGRIF_END_DO_NOT_TREAT     
2569      !
2570   END FUNCTION iom_sdate
2571
2572   INTEGER FUNCTION iom_add(idm)
2573   INTEGER :: idm
2574
2575   iom_add = idm + 1
2576   IF( iom_add > max_rst_fields) THEN
2577     WRITE(ctmp1,*) 'E R R O R : iom_set_rst_vars SIZE of RST_FIELD array is too small'
2578     CALL ctl_stop( 'iom_add:', ctmp1 )
2579   ENDIF
2580
2581   END FUNCTION iom_add
2582
2583#else
2584   !!----------------------------------------------------------------------
2585   !!   NOT 'key_iomput'                               a few dummy routines
2586   !!----------------------------------------------------------------------
2587   SUBROUTINE iom_setkt( kt, cdname )
2588      INTEGER         , INTENT(in)::   kt 
2589      CHARACTER(LEN=*), INTENT(in) ::   cdname
2590      IF( .FALSE. )   WRITE(numout,*) kt, cdname   ! useless test to avoid compilation warnings
2591   END SUBROUTINE iom_setkt
2592
2593   SUBROUTINE iom_context_finalize( cdname )
2594      CHARACTER(LEN=*), INTENT(in) ::   cdname
2595      IF( .FALSE. )   WRITE(numout,*)  cdname   ! useless test to avoid compilation warnings
2596   END SUBROUTINE iom_context_finalize
2597
2598#endif
2599
2600   LOGICAL FUNCTION iom_use( cdname )
2601      CHARACTER(LEN=*), INTENT(in) ::   cdname
2602#if defined key_iomput
2603      iom_use = xios_field_is_active( cdname )
2604#else
2605      iom_use = .FALSE.
2606#endif
2607   END FUNCTION iom_use
2608
2609   SUBROUTINE iom_miss_val( cdname, pmiss_val )
2610      CHARACTER(LEN=*), INTENT(in ) ::   cdname
2611      REAL(wp)        , INTENT(out) ::   pmiss_val   
2612#if defined key_iomput
2613      ! get missing value
2614      CALL xios_get_field_attr( cdname, default_value = pmiss_val )
2615#else
2616      IF( .FALSE. )   WRITE(numout,*) cdname, pmiss_val   ! useless test to avoid compilation warnings
2617#endif
2618   END SUBROUTINE iom_miss_val
2619 
2620   !!======================================================================
2621END MODULE iom
Note: See TracBrowser for help on using the repository browser.