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 utils/tools/DOMAINcfg/src – NEMO

source: utils/tools/DOMAINcfg/src/iom.F90 @ 14243

Last change on this file since 14243 was 14243, checked in by smasson, 3 years ago

trunk: replace key_iomput by key_xios

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