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

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

source: NEMO/branches/2019/ENHANCE-03_domcfg/src/iom.F90 @ 11602

Last change on this file since 11602 was 11196, checked in by mathiot, 5 years ago

ENHANCE-03_domcfg: USE dianam useless in iom.F90 (ticket #2143)

File size: 124.6 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_iomput
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_iomput
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_iomput
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_iomput
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_iomput
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_iomput 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_iomput
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_iomput
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_iomput
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_iomput
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      ! Initializations and control
617      ! =============
618      kiomid = -1
619      clinfo = '                    iom_open ~~~  '
620      istop = nstop
621      ! if iom_open is called for the first time: initialize iom_file(:)%nfid to 0
622      ! (could be done when defining iom_file in f95 but not in f90)
623      IF( Agrif_Root() ) THEN
624         IF( iom_open_init == 0 ) THEN
625            iom_file(:)%nfid = 0
626            iom_open_init = 1
627         ENDIF
628      ENDIF
629      ! do we read or write the file?
630      IF( PRESENT(ldwrt) ) THEN   ;   llwrt = ldwrt
631      ELSE                        ;   llwrt = .FALSE.
632      ENDIF
633      ! do we call ctl_stop if we try to open a non-existing file in read mode?
634      IF( PRESENT(ldstop) ) THEN   ;   llstop = ldstop
635      ELSE                         ;   llstop = .TRUE.
636      ENDIF
637      ! do we add agrif suffix
638      IF( PRESENT(lagrif) ) THEN   ;   llagrif = lagrif
639      ELSE                         ;   llagrif = .TRUE.
640      ENDIF
641      ! are we using interpolation on the fly?
642      IF( PRESENT(ldiof) ) THEN   ;   lliof = ldiof
643      ELSE                        ;   lliof = .FALSE.
644      ENDIF
645      ! do we read the overlap
646      ! ugly patch SM+JMM+RB to overwrite global definition in some cases
647      llnoov = (jpni * jpnj ) == jpnij .AND. .NOT. lk_agrif
648      ! create the file name by added, if needed, TRIM(Agrif_CFixed()) and TRIM(clsuffix)
649      ! =============
650      clname   = trim(cdname)
651      IF ( .NOT. Agrif_Root() .AND. .NOT. lliof .AND. llagrif) THEN
652         iln    = INDEX(clname,'/') 
653         cltmpn = clname(1:iln)
654         clname = clname(iln+1:LEN_TRIM(clname))
655         clname=TRIM(cltmpn)//TRIM(Agrif_CFixed())//'_'//TRIM(clname)
656      ENDIF
657      ! which suffix should we use?
658      clsuffix = '.nc'
659      ! Add the suffix if needed
660      iln = LEN_TRIM(clname)
661      ils = LEN_TRIM(clsuffix)
662      IF( iln <= ils .OR. INDEX( TRIM(clname), TRIM(clsuffix), back = .TRUE. ) /= iln - ils + 1 )   &
663         &   clname = TRIM(clname)//TRIM(clsuffix)
664      cltmpn = clname   ! store this name
665      ! try to find if the file to be opened already exist
666      ! =============
667      INQUIRE( FILE = clname, EXIST = llok )
668      IF( .NOT.llok ) THEN
669         ! we try to add the cpu number to the name
670         WRITE(clcpu,*) narea-1
671
672         clcpu  = TRIM(ADJUSTL(clcpu))
673         iln = INDEX(clname,TRIM(clsuffix), back = .TRUE.)
674         clname = clname(1:iln-1)//'_'//TRIM(clcpu)//TRIM(clsuffix)
675         icnt = 0
676         INQUIRE( FILE = clname, EXIST = llok ) 
677         ! we try different formats for the cpu number by adding 0
678         DO WHILE( .NOT.llok .AND. icnt < jpmax_digits )
679            clcpu  = "0"//trim(clcpu)
680            clname = clname(1:iln-1)//'_'//TRIM(clcpu)//TRIM(clsuffix)
681            INQUIRE( FILE = clname, EXIST = llok )
682            icnt = icnt + 1
683         END DO
684      ELSE
685         lxios_sini = .TRUE.
686      ENDIF
687      IF( llwrt ) THEN
688         ! check the domain definition
689! JMM + SM: ugly patch before getting the new version of lib_mpp)
690!         idom = jpdom_local_noovlap   ! default definition
691         IF( llnoov ) THEN   ;   idom = jpdom_local_noovlap   ! default definition
692         ELSE                ;   idom = jpdom_local_full      ! default definition
693         ENDIF
694         IF( PRESENT(kdom) )   idom = kdom
695         ! create the domain informations
696         ! =============
697         SELECT CASE (idom)
698         CASE (jpdom_local_full)
699            idompar(:,1) = (/ jpi             , jpj              /)
700            idompar(:,2) = (/ nimpp           , njmpp            /)
701            idompar(:,3) = (/ nimpp + jpi - 1 , njmpp + jpj - 1  /)
702            idompar(:,4) = (/ nldi - 1        , nldj - 1         /)
703            idompar(:,5) = (/ jpi - nlei      , jpj - nlej       /)
704         CASE (jpdom_local_noextra)
705            idompar(:,1) = (/ nlci            , nlcj             /)
706            idompar(:,2) = (/ nimpp           , njmpp            /)
707            idompar(:,3) = (/ nimpp + nlci - 1, njmpp + nlcj - 1 /)
708            idompar(:,4) = (/ nldi - 1        , nldj - 1         /)
709            idompar(:,5) = (/ nlci - nlei     , nlcj - nlej      /)
710         CASE (jpdom_local_noovlap)
711            idompar(:,1) = (/ nlei  - nldi + 1, nlej  - nldj + 1 /)
712            idompar(:,2) = (/ nimpp + nldi - 1, njmpp + nldj - 1 /)
713            idompar(:,3) = (/ nimpp + nlei - 1, njmpp + nlej - 1 /)
714            idompar(:,4) = (/ 0               , 0                /)
715            idompar(:,5) = (/ 0               , 0                /)
716         CASE DEFAULT
717            CALL ctl_stop( TRIM(clinfo), 'wrong value of kdom, only jpdom_local* cases are accepted' )
718         END SELECT
719      ENDIF
720      ! Open the NetCDF file
721      ! =============
722      ! do we have some free file identifier?
723      IF( MINVAL(iom_file(:)%nfid) /= 0 )   &
724         &   CALL ctl_stop( TRIM(clinfo), 'No more free file identifier', 'increase jpmax_files in iom_def' )
725      ! if no file was found...
726      IF( .NOT. llok ) THEN
727         IF( .NOT. llwrt ) THEN   ! we are in read mode
728            IF( llstop ) THEN   ;   CALL ctl_stop( TRIM(clinfo), 'File '//TRIM(cltmpn)//'* not found' )
729            ELSE                ;   istop = nstop + 1   ! make sure that istop /= nstop so we don't open the file
730            ENDIF
731         ELSE                     ! we are in write mode so we
732            clname = cltmpn       ! get back the file name without the cpu number
733         ENDIF
734      ELSE
735         IF( llwrt .AND. .NOT. ln_clobber ) THEN   ! we stop as we want to write in a new file
736            CALL ctl_stop( TRIM(clinfo), 'We want to write in a new file but '//TRIM(clname)//' already exists...' )
737            istop = nstop + 1                      ! make sure that istop /= nstop so we don't open the file
738         ELSEIF( llwrt ) THEN     ! the file exists and we are in write mode with permission to
739            clname = cltmpn       ! overwrite so get back the file name without the cpu number
740         ENDIF
741      ENDIF
742      IF( istop == nstop ) THEN   ! no error within this routine
743         CALL iom_nf90_open( clname, kiomid, llwrt, llok, idompar, kdlev = kdlev )
744      ENDIF
745      !
746   END SUBROUTINE iom_open
747
748
749   SUBROUTINE iom_close( kiomid )
750      !!--------------------------------------------------------------------
751      !!                   ***  SUBROUTINE  iom_close  ***
752      !!
753      !! ** Purpose : close an input file, or all files opened by iom
754      !!--------------------------------------------------------------------
755      INTEGER, INTENT(inout), OPTIONAL ::   kiomid   ! iom identifier of the file to be closed
756      !                                              ! return 0 when file is properly closed
757      !                                              ! No argument: all files opened by iom are closed
758
759      INTEGER ::   jf         ! dummy loop indices
760      INTEGER ::   i_s, i_e   ! temporary integer
761      CHARACTER(LEN=100)    ::   clinfo    ! info character
762      !---------------------------------------------------------------------
763      !
764      clinfo = '                    iom_close ~~~  '
765      IF( PRESENT(kiomid) ) THEN
766         i_s = kiomid
767         i_e = kiomid
768      ELSE
769         i_s = 1
770         i_e = jpmax_files
771      ENDIF
772
773      IF( i_s > 0 ) THEN
774         DO jf = i_s, i_e
775            IF( iom_file(jf)%nfid > 0 ) THEN
776               CALL iom_nf90_close( jf )
777               iom_file(jf)%nfid       = 0          ! free the id
778               IF( PRESENT(kiomid) )   kiomid = 0   ! return 0 as id to specify that the file was closed
779               IF(lwp) WRITE(numout,*) TRIM(clinfo)//' close file: '//TRIM(iom_file(jf)%name)//' ok'
780            ELSEIF( PRESENT(kiomid) ) THEN
781               WRITE(ctmp1,*) '--->',  kiomid
782               CALL ctl_stop( TRIM(clinfo)//' Invalid file identifier', ctmp1 )
783            ENDIF
784         END DO
785      ENDIF
786      !   
787   END SUBROUTINE iom_close
788
789
790   FUNCTION iom_varid ( kiomid, cdvar, kdimsz, kndims, ldstop ) 
791      !!-----------------------------------------------------------------------
792      !!                  ***  FUNCTION  iom_varid  ***
793      !!
794      !! ** Purpose : get the id of a variable in a file (return 0 if not found)
795      !!-----------------------------------------------------------------------
796      INTEGER              , INTENT(in   )           ::   kiomid   ! file Identifier
797      CHARACTER(len=*)     , INTENT(in   )           ::   cdvar    ! name of the variable
798      INTEGER, DIMENSION(:), INTENT(  out), OPTIONAL ::   kdimsz   ! size of each dimension
799      INTEGER,               INTENT(  out), OPTIONAL ::   kndims   ! size of the dimensions
800      LOGICAL              , INTENT(in   ), OPTIONAL ::   ldstop   ! stop if looking for non-existing variable (default = .TRUE.)
801      !
802      INTEGER                        ::   iom_varid, iiv, i_nvd
803      LOGICAL                        ::   ll_fnd
804      CHARACTER(LEN=100)             ::   clinfo                   ! info character
805      LOGICAL                        ::   llstop                   ! local definition of ldstop
806      !!-----------------------------------------------------------------------
807      iom_varid = 0                         ! default definition
808      ! do we call ctl_stop if we look for non-existing variable?
809      IF( PRESENT(ldstop) ) THEN   ;   llstop = ldstop
810      ELSE                         ;   llstop = .TRUE.
811      ENDIF
812      !
813      IF( kiomid > 0 ) THEN
814         clinfo = 'iom_varid, file: '//trim(iom_file(kiomid)%name)//', var: '//trim(cdvar)
815         IF( iom_file(kiomid)%nfid == 0 ) THEN
816            CALL ctl_stop( trim(clinfo), 'the file is not open' )
817         ELSE
818            ll_fnd  = .FALSE.
819            iiv = 0
820            !
821            DO WHILE ( .NOT.ll_fnd .AND. iiv < iom_file(kiomid)%nvars )
822               iiv = iiv + 1
823               ll_fnd  = ( TRIM(cdvar) == TRIM(iom_file(kiomid)%cn_var(iiv)) )
824            END DO
825            !
826            IF( .NOT.ll_fnd ) THEN
827               iiv = iiv + 1
828               IF( iiv <= jpmax_vars ) THEN
829                  iom_varid = iom_nf90_varid( kiomid, cdvar, iiv, kdimsz, kndims )
830               ELSE
831                  CALL ctl_stop( trim(clinfo), 'Too many variables in the file '//iom_file(kiomid)%name,   &
832                        &                      'increase the parameter jpmax_vars')
833               ENDIF
834               IF( llstop .AND. iom_varid == -1 )   CALL ctl_stop( TRIM(clinfo)//' not found' ) 
835            ELSE
836               iom_varid = iiv
837               IF( PRESENT(kdimsz) ) THEN
838                  i_nvd = iom_file(kiomid)%ndims(iiv)
839                  IF( i_nvd <= size(kdimsz) ) THEN
840                     kdimsz(1:i_nvd) = iom_file(kiomid)%dimsz(1:i_nvd,iiv)
841                  ELSE
842                     WRITE(ctmp1,*) i_nvd, size(kdimsz)
843                     CALL ctl_stop( trim(clinfo), 'error in kdimsz size'//trim(ctmp1) )
844                  ENDIF
845               ENDIF
846               IF( PRESENT(kndims) )  kndims = iom_file(kiomid)%ndims(iiv)
847            ENDIF
848         ENDIF
849      ENDIF
850      !
851   END FUNCTION iom_varid
852
853
854   !!----------------------------------------------------------------------
855   !!                   INTERFACE iom_get
856   !!----------------------------------------------------------------------
857   SUBROUTINE iom_g0d( kiomid, cdvar, pvar, ktime, ldxios )
858      INTEGER         , INTENT(in   )                 ::   kiomid    ! Identifier of the file
859      CHARACTER(len=*), INTENT(in   )                 ::   cdvar     ! Name of the variable
860      REAL(wp)        , INTENT(  out)                 ::   pvar      ! read field
861      INTEGER         , INTENT(in   ),     OPTIONAL   ::   ktime     ! record number
862      LOGICAL         , INTENT(in   ),     OPTIONAL   ::   ldxios    ! use xios to read restart
863      !
864      INTEGER                                         ::   idvar     ! variable id
865      INTEGER                                         ::   idmspc    ! number of spatial dimensions
866      INTEGER         , DIMENSION(1)                  ::   itime     ! record number
867      CHARACTER(LEN=100)                              ::   clinfo    ! info character
868      CHARACTER(LEN=100)                              ::   clname    ! file name
869      CHARACTER(LEN=1)                                ::   cldmspc   !
870      LOGICAL                                         ::   llxios
871      !
872      llxios = .FALSE.
873      IF( PRESENT(ldxios) ) llxios = ldxios
874
875      IF(.NOT.llxios) THEN  ! read data using default library
876         itime = 1
877         IF( PRESENT(ktime) ) itime = ktime
878         !
879         clname = iom_file(kiomid)%name
880         clinfo = '          iom_g0d, file: '//trim(clname)//', var: '//trim(cdvar)
881         !
882         IF( kiomid > 0 ) THEN
883            idvar = iom_varid( kiomid, cdvar )
884            IF( iom_file(kiomid)%nfid > 0 .AND. idvar > 0 ) THEN
885               idmspc = iom_file ( kiomid )%ndims( idvar )
886               IF( iom_file(kiomid)%luld(idvar) )  idmspc = idmspc - 1
887               WRITE(cldmspc , fmt='(i1)') idmspc
888               IF( idmspc > 0 )  CALL ctl_stop( TRIM(clinfo), 'When reading to a 0D array, we do not accept data', &
889                                    &                         'with 1 or more spatial dimensions: '//cldmspc//' were found.' , &
890                                    &                         'Use ncwa -a to suppress the unnecessary dimensions' )
891               CALL iom_nf90_get( kiomid, idvar, pvar, itime )
892            ENDIF
893         ENDIF
894      ELSE
895#if defined key_iomput
896         IF(lwp) WRITE(numout,*) 'XIOS RST READ (0D): ', trim(cdvar)
897         CALL iom_swap( TRIM(crxios_context) )
898         CALL xios_recv_field( trim(cdvar), pvar)
899         CALL iom_swap( TRIM(cxios_context) )
900#else
901         WRITE(ctmp1,*) 'Can not use XIOS in iom_g0d, file: '//trim(clname)//', var:'//trim(cdvar)
902         CALL ctl_stop( 'iom_g0d', ctmp1 )
903#endif
904      ENDIF
905   END SUBROUTINE iom_g0d
906
907   SUBROUTINE iom_g1d( kiomid, kdom, cdvar, pvar, ktime, kstart, kcount, ldxios )
908      INTEGER         , INTENT(in   )                         ::   kiomid    ! Identifier of the file
909      INTEGER         , INTENT(in   )                         ::   kdom      ! Type of domain to be read
910      CHARACTER(len=*), INTENT(in   )                         ::   cdvar     ! Name of the variable
911      REAL(wp)        , INTENT(  out), DIMENSION(:)           ::   pvar      ! read field
912      INTEGER         , INTENT(in   )              , OPTIONAL ::   ktime     ! record number
913      INTEGER         , INTENT(in   ), DIMENSION(1), OPTIONAL ::   kstart    ! start axis position of the reading
914      INTEGER         , INTENT(in   ), DIMENSION(1), OPTIONAL ::   kcount    ! number of points in each axis
915      LOGICAL         , INTENT(in   ),               OPTIONAL ::   ldxios    ! read data using XIOS
916      !
917      IF( kiomid > 0 ) THEN
918         IF( iom_file(kiomid)%nfid > 0 ) CALL iom_get_123d( kiomid, kdom       , cdvar        , pv_r1d=pvar,   &
919              &                                                     ktime=ktime, kstart=kstart, kcount=kcount, &
920              &                                                     ldxios=ldxios )
921      ENDIF
922   END SUBROUTINE iom_g1d
923
924   SUBROUTINE iom_g2d( kiomid, kdom, cdvar, pvar, ktime, kstart, kcount, lrowattr, ldxios)
925      INTEGER         , INTENT(in   )                           ::   kiomid    ! Identifier of the file
926      INTEGER         , INTENT(in   )                           ::   kdom      ! Type of domain to be read
927      CHARACTER(len=*), INTENT(in   )                           ::   cdvar     ! Name of the variable
928      REAL(wp)        , INTENT(  out), DIMENSION(:,:)           ::   pvar      ! read field
929      INTEGER         , INTENT(in   )                , OPTIONAL ::   ktime     ! record number
930      INTEGER         , INTENT(in   ), DIMENSION(2)  , OPTIONAL ::   kstart    ! start axis position of the reading
931      INTEGER         , INTENT(in   ), DIMENSION(2)  , OPTIONAL ::   kcount    ! number of points in each axis
932      LOGICAL         , INTENT(in   )                , OPTIONAL ::   lrowattr  ! logical flag telling iom_get to
933                                                                               ! look for and use a file attribute
934                                                                               ! called open_ocean_jstart to set the start
935                                                                               ! value for the 2nd dimension (netcdf only)
936      LOGICAL         , INTENT(in   ),               OPTIONAL ::   ldxios      ! read data using XIOS
937      !
938      IF( kiomid > 0 ) THEN
939         IF( iom_file(kiomid)%nfid > 0 ) CALL iom_get_123d( kiomid, kdom       , cdvar        , pv_r2d=pvar,   &
940              &                                                     ktime=ktime, kstart=kstart, kcount=kcount, &
941              &                                                     lrowattr=lrowattr,  ldxios=ldxios)
942      ENDIF
943   END SUBROUTINE iom_g2d
944
945   SUBROUTINE iom_g3d( kiomid, kdom, cdvar, pvar, ktime, kstart, kcount, lrowattr, ldxios )
946      INTEGER         , INTENT(in   )                             ::   kiomid    ! Identifier of the file
947      INTEGER         , INTENT(in   )                             ::   kdom      ! Type of domain to be read
948      CHARACTER(len=*), INTENT(in   )                             ::   cdvar     ! Name of the variable
949      REAL(wp)        , INTENT(  out), DIMENSION(:,:,:)           ::   pvar      ! read field
950      INTEGER         , INTENT(in   )                  , OPTIONAL ::   ktime     ! record number
951      INTEGER         , INTENT(in   ), DIMENSION(3)    , OPTIONAL ::   kstart    ! start axis position of the reading
952      INTEGER         , INTENT(in   ), DIMENSION(3)    , OPTIONAL ::   kcount    ! number of points in each axis
953      LOGICAL         , INTENT(in   )                  , OPTIONAL ::   lrowattr  ! logical flag telling iom_get to
954                                                                                 ! look for and use a file attribute
955                                                                                 ! called open_ocean_jstart to set the start
956                                                                                 ! value for the 2nd dimension (netcdf only)
957      LOGICAL         , INTENT(in   ),               OPTIONAL ::   ldxios        ! read data using XIOS
958      !
959      IF( kiomid > 0 ) THEN
960         IF( iom_file(kiomid)%nfid > 0 ) CALL iom_get_123d( kiomid, kdom       , cdvar        , pv_r3d=pvar,   &
961              &                                                     ktime=ktime, kstart=kstart, kcount=kcount, &
962              &                                                     lrowattr=lrowattr, ldxios=ldxios )
963      ENDIF
964   END SUBROUTINE iom_g3d
965   !!----------------------------------------------------------------------
966
967   SUBROUTINE iom_get_123d( kiomid, kdom  , cdvar ,   &
968         &                  pv_r1d, pv_r2d, pv_r3d,   &
969         &                  ktime , kstart, kcount,   &
970         &                  lrowattr, ldxios        )
971      !!-----------------------------------------------------------------------
972      !!                  ***  ROUTINE  iom_get_123d  ***
973      !!
974      !! ** Purpose : read a 1D/2D/3D variable
975      !!
976      !! ** Method : read ONE record at each CALL
977      !!-----------------------------------------------------------------------
978      INTEGER                    , INTENT(in   )           ::   kiomid     ! Identifier of the file
979      INTEGER                    , INTENT(in   )           ::   kdom       ! Type of domain to be read
980      CHARACTER(len=*)           , INTENT(in   )           ::   cdvar      ! Name of the variable
981      REAL(wp), DIMENSION(:)     , INTENT(  out), OPTIONAL ::   pv_r1d     ! read field (1D case)
982      REAL(wp), DIMENSION(:,:)   , INTENT(  out), OPTIONAL ::   pv_r2d     ! read field (2D case)
983      REAL(wp), DIMENSION(:,:,:) , INTENT(  out), OPTIONAL ::   pv_r3d     ! read field (3D case)
984      INTEGER                    , INTENT(in   ), OPTIONAL ::   ktime      ! record number
985      INTEGER , DIMENSION(:)     , INTENT(in   ), OPTIONAL ::   kstart     ! start position of the reading in each axis
986      INTEGER , DIMENSION(:)     , INTENT(in   ), OPTIONAL ::   kcount     ! number of points to be read in each axis
987      LOGICAL                    , INTENT(in   ), OPTIONAL ::   lrowattr   ! logical flag telling iom_get to
988                                                                           ! look for and use a file attribute
989                                                                           ! called open_ocean_jstart to set the start
990                                                                           ! value for the 2nd dimension (netcdf only)
991      LOGICAL                    , INTENT(in   ), OPTIONAL ::   ldxios     ! use XIOS to read restart
992      !
993      LOGICAL                        ::   llxios       ! local definition for XIOS read
994      LOGICAL                        ::   llnoov      ! local definition to read overlap
995      LOGICAL                        ::   luse_jattr  ! local definition to read open_ocean_jstart file attribute
996      INTEGER                        ::   jstartrow   ! start point for 2nd dimension optionally set by file attribute
997      INTEGER                        ::   jl          ! loop on number of dimension
998      INTEGER                        ::   idom        ! type of domain
999      INTEGER                        ::   idvar       ! id of the variable
1000      INTEGER                        ::   inbdim      ! number of dimensions of the variable
1001      INTEGER                        ::   idmspc      ! number of spatial dimensions
1002      INTEGER                        ::   itime       ! record number
1003      INTEGER                        ::   istop       ! temporary value of nstop
1004      INTEGER                        ::   ix1, ix2, iy1, iy2   ! subdomain indexes
1005      INTEGER                        ::   ji, jj      ! loop counters
1006      INTEGER                        ::   irankpv     !
1007      INTEGER                        ::   ind1, ind2  ! substring index
1008      INTEGER, DIMENSION(jpmax_dims) ::   istart      ! starting point to read for each axis
1009      INTEGER, DIMENSION(jpmax_dims) ::   icnt        ! number of value to read along each axis
1010      INTEGER, DIMENSION(jpmax_dims) ::   idimsz      ! size of the dimensions of the variable
1011      INTEGER, DIMENSION(jpmax_dims) ::   ishape      ! size of the dimensions of the variable
1012      REAL(wp)                       ::   zscf, zofs  ! sacle_factor and add_offset
1013      INTEGER                        ::   itmp        ! temporary integer
1014      CHARACTER(LEN=256)             ::   clinfo      ! info character
1015      CHARACTER(LEN=256)             ::   clname      ! file name
1016      CHARACTER(LEN=1)               ::   clrankpv, cldmspc      !
1017      LOGICAL                        ::   ll_depth_spec ! T => if kstart, kcount present then *only* use values for 3rd spatial dimension.
1018      INTEGER                        ::   inlev       ! number of levels for 3D data
1019      REAL(wp)                       ::   gma, gmi
1020      !---------------------------------------------------------------------
1021      !
1022      inlev = -1
1023      IF( PRESENT(pv_r3d) )   inlev = SIZE(pv_r3d, 3)
1024      !
1025      llxios = .FALSE.
1026      if(PRESENT(ldxios)) llxios = ldxios
1027      idvar = iom_varid( kiomid, cdvar ) 
1028      idom = kdom
1029      !
1030      IF(.NOT.llxios) THEN
1031         clname = iom_file(kiomid)%name   !   esier to read
1032         clinfo = '          iom_get_123d, file: '//trim(clname)//', var: '//trim(cdvar)
1033         ! local definition of the domain ?
1034         ! do we read the overlap
1035         ! ugly patch SM+JMM+RB to overwrite global definition in some cases
1036         llnoov = (jpni * jpnj ) == jpnij .AND. .NOT. lk_agrif 
1037         ! check kcount and kstart optionals parameters...
1038         IF( PRESENT(kcount) .AND. (.NOT. PRESENT(kstart)) ) CALL ctl_stop(trim(clinfo), 'kcount present needs kstart present')
1039         IF( PRESENT(kstart) .AND. (.NOT. PRESENT(kcount)) ) CALL ctl_stop(trim(clinfo), 'kstart present needs kcount present')
1040         IF( PRESENT(kstart) .AND. idom /= jpdom_unknown .AND.  idom /= jpdom_autoglo_xy  ) &
1041     &          CALL ctl_stop(trim(clinfo), 'kstart present needs kdom = jpdom_unknown or kdom = jpdom_autoglo_xy')
1042
1043         luse_jattr = .false.
1044         IF( PRESENT(lrowattr) ) THEN
1045            IF( lrowattr .AND. idom /= jpdom_data   ) CALL ctl_stop(trim(clinfo), 'lrowattr present and true needs kdom = jpdom_data')
1046            IF( lrowattr .AND. idom == jpdom_data   ) luse_jattr = .true.
1047         ENDIF
1048
1049         ! Search for the variable in the data base (eventually actualize data)
1050         istop = nstop
1051         !
1052         IF( idvar > 0 ) THEN
1053            ! to write iom_file(kiomid)%dimsz in a shorter way !
1054            idimsz(:) = iom_file(kiomid)%dimsz(:, idvar) 
1055            inbdim = iom_file(kiomid)%ndims(idvar)            ! number of dimensions in the file
1056            idmspc = inbdim                                   ! number of spatial dimensions in the file
1057            IF( iom_file(kiomid)%luld(idvar) )   idmspc = inbdim - 1
1058            IF( idmspc > 3 )   CALL ctl_stop(trim(clinfo), 'the file has more than 3 spatial dimensions this case is not coded...') 
1059            !
1060            ! update idom definition...
1061            ! Identify the domain in case of jpdom_auto(glo/dta) definition
1062            IF( idom == jpdom_autoglo_xy ) THEN
1063               ll_depth_spec = .TRUE.
1064               idom = jpdom_autoglo
1065            ELSE
1066               ll_depth_spec = .FALSE.
1067            ENDIF
1068            IF( idom == jpdom_autoglo .OR. idom == jpdom_autodta ) THEN           
1069               IF( idom == jpdom_autoglo ) THEN   ;   idom = jpdom_global 
1070               ELSE                               ;   idom = jpdom_data
1071               ENDIF
1072               ind1 = INDEX( clname, '_', back = .TRUE. ) + 1
1073               ind2 = INDEX( clname, '.', back = .TRUE. ) - 1
1074               IF( ind2 > ind1 ) THEN   ;   IF( VERIFY( clname(ind1:ind2), '0123456789' ) == 0 )   idom = jpdom_local   ;   ENDIF
1075            ENDIF
1076            ! Identify the domain in case of jpdom_local definition
1077            IF( idom == jpdom_local ) THEN
1078               IF(     idimsz(1) == jpi               .AND. idimsz(2) == jpj               ) THEN   ;   idom = jpdom_local_full
1079               ELSEIF( idimsz(1) == nlci              .AND. idimsz(2) == nlcj              ) THEN   ;   idom = jpdom_local_noextra
1080               ELSEIF( idimsz(1) == (nlei - nldi + 1) .AND. idimsz(2) == (nlej - nldj + 1) ) THEN   ;   idom = jpdom_local_noovlap
1081               ELSE   ;   CALL ctl_stop( trim(clinfo), 'impossible to identify the local domain' )
1082               ENDIF
1083            ENDIF
1084            !
1085            ! check the consistency between input array and data rank in the file
1086            !
1087            ! initializations
1088            itime = 1
1089            IF( PRESENT(ktime) ) itime = ktime
1090            !
1091            irankpv = 1 * COUNT( (/PRESENT(pv_r1d)/) ) + 2 * COUNT( (/PRESENT(pv_r2d)/) ) + 3 * COUNT( (/PRESENT(pv_r3d)/) )
1092            WRITE(clrankpv, fmt='(i1)') irankpv
1093            WRITE(cldmspc , fmt='(i1)') idmspc
1094            !
1095            IF(     idmspc <  irankpv ) THEN
1096               CALL ctl_stop( TRIM(clinfo), 'The file has only '//cldmspc//' spatial dimension',   &
1097                  &                         'it is impossible to read a '//clrankpv//'D array from this file...' )
1098            ELSEIF( idmspc == irankpv ) THEN
1099               IF( PRESENT(pv_r1d) .AND. idom /= jpdom_unknown )   &
1100                  &   CALL ctl_stop( TRIM(clinfo), 'case not coded...You must use jpdom_unknown' )
1101            ELSEIF( idmspc >  irankpv ) THEN
1102                  IF( PRESENT(pv_r2d) .AND. itime == 1 .AND. idimsz(3) == 1 .AND. idmspc == 3 ) THEN
1103                     CALL ctl_warn( trim(clinfo), '2D array but 3 spatial dimensions for the data...'              ,   &
1104                           &         'As the size of the z dimension is 1 and as we try to read the first record, ',   &
1105                           &         'we accept this case, even if there is a possible mix-up between z and time dimension' )   
1106                     idmspc = idmspc - 1
1107                  ELSE
1108                     CALL ctl_stop( TRIM(clinfo), 'To keep iom lisibility, when reading a '//clrankpv//'D array,'         ,   &
1109                        &                         'we do not accept data with '//cldmspc//' spatial dimensions',   &
1110                        &                         'Use ncwa -a to suppress the unnecessary dimensions' )
1111                  ENDIF
1112            ENDIF
1113            !
1114            ! definition of istart and icnt
1115            !
1116            icnt  (:) = 1
1117            istart(:) = 1
1118            istart(idmspc+1) = itime
1119   
1120            IF( PRESENT(kstart) .AND. .NOT. ll_depth_spec ) THEN
1121               istart(1:idmspc) = kstart(1:idmspc) 
1122               icnt  (1:idmspc) = kcount(1:idmspc)
1123            ELSE
1124               IF(idom == jpdom_unknown ) THEN
1125                  icnt(1:idmspc) = idimsz(1:idmspc)
1126               ELSE
1127                  IF( .NOT. PRESENT(pv_r1d) ) THEN   !   not a 1D array
1128                     IF(     idom == jpdom_data    ) THEN
1129                        jstartrow = 1
1130                        IF( luse_jattr ) THEN
1131                           CALL iom_getatt(kiomid, 'open_ocean_jstart', jstartrow ) ! -999 is returned if the attribute is not found
1132                           jstartrow = MAX(1,jstartrow)
1133                        ENDIF
1134                        istart(1:2) = (/ mig(1), mjg(1) + jstartrow - 1 /)  ! icnt(1:2) done below
1135                     ELSEIF( idom == jpdom_global  ) THEN ; istart(1:2) = (/ nimpp , njmpp  /)  ! icnt(1:2) done below
1136                     ENDIF
1137                     ! we do not read the overlap                     -> we start to read at nldi, nldj
1138! JMM + SM: ugly patch before getting the new version of lib_mpp)
1139!                  IF( idom /= jpdom_local_noovlap )   istart(1:2) = istart(1:2) + (/ nldi - 1, nldj - 1 /)
1140                     IF( llnoov .AND. idom /= jpdom_local_noovlap ) istart(1:2) = istart(1:2) + (/ nldi - 1, nldj - 1 /)
1141                  ! we do not read the overlap and the extra-halos -> from nldi to nlei and from nldj to nlej
1142! JMM + SM: ugly patch before getting the new version of lib_mpp)
1143!                  icnt(1:2) = (/ nlei - nldi + 1, nlej - nldj + 1 /)
1144                     IF( llnoov ) THEN   ;   icnt(1:2) = (/ nlei - nldi + 1, nlej - nldj + 1 /)
1145                     ELSE                ;   icnt(1:2) = (/ nlci           , nlcj            /)
1146                     ENDIF
1147                     IF( PRESENT(pv_r3d) ) THEN
1148                        IF( idom == jpdom_data ) THEN                        ;                               icnt(3) = inlev
1149                        ELSEIF( ll_depth_spec .AND. PRESENT(kstart) ) THEN   ;   istart(3) = kstart(3)   ;   icnt(3) = kcount(3)
1150                        ELSE                                                 ;                               icnt(3) = inlev
1151                        ENDIF
1152                     ENDIF
1153                  ENDIF
1154               ENDIF
1155            ENDIF
1156
1157            ! check that istart and icnt can be used with this file
1158            !-
1159            DO jl = 1, jpmax_dims
1160               itmp = istart(jl)+icnt(jl)-1
1161               IF( itmp > idimsz(jl) .AND. idimsz(jl) /= 0 ) THEN
1162                  WRITE( ctmp1, FMT="('(istart(', i1, ') + icnt(', i1, ') - 1) = ', i5)" ) jl, jl, itmp
1163                  WRITE( ctmp2, FMT="(' is larger than idimsz(', i1,') = ', i5)"         ) jl, idimsz(jl)
1164                  CALL ctl_stop( trim(clinfo), 'start and count too big regarding to the size of the data, ', ctmp1, ctmp2 )     
1165               ENDIF
1166            END DO
1167
1168            ! check that icnt matches the input array
1169            !-     
1170            IF( idom == jpdom_unknown ) THEN
1171               IF( irankpv == 1 )        ishape(1:1) = SHAPE(pv_r1d)
1172               IF( irankpv == 2 )        ishape(1:2) = SHAPE(pv_r2d)
1173               IF( irankpv == 3 )        ishape(1:3) = SHAPE(pv_r3d)
1174               ctmp1 = 'd'
1175            ELSE
1176               IF( irankpv == 2 ) THEN
1177! JMM + SM: ugly patch before getting the new version of lib_mpp)
1178!               ishape(1:2) = SHAPE(pv_r2d(nldi:nlei,nldj:nlej  ))   ;   ctmp1 = 'd(nldi:nlei,nldj:nlej)'
1179                  IF( llnoov ) THEN ; ishape(1:2)=SHAPE(pv_r2d(nldi:nlei,nldj:nlej  )) ; ctmp1='d(nldi:nlei,nldj:nlej)'
1180                  ELSE              ; ishape(1:2)=SHAPE(pv_r2d(1   :nlci,1   :nlcj  )) ; ctmp1='d(1:nlci,1:nlcj)'
1181                  ENDIF
1182               ENDIF
1183               IF( irankpv == 3 ) THEN 
1184! JMM + SM: ugly patch before getting the new version of lib_mpp)
1185!               ishape(1:3) = SHAPE(pv_r3d(nldi:nlei,nldj:nlej,:))   ;   ctmp1 = 'd(nldi:nlei,nldj:nlej,:)'
1186                  IF( llnoov ) THEN ; ishape(1:3)=SHAPE(pv_r3d(nldi:nlei,nldj:nlej,:)) ; ctmp1='d(nldi:nlei,nldj:nlej,:)'
1187                  ELSE              ; ishape(1:3)=SHAPE(pv_r3d(1   :nlci,1   :nlcj,:)) ; ctmp1='d(1:nlci,1:nlcj,:)'
1188                  ENDIF
1189               ENDIF
1190            ENDIF
1191         
1192            DO jl = 1, irankpv
1193               WRITE( ctmp2, FMT="(', ', i1,'): ', i5,' /= icnt(', i1,'):', i5)" ) jl, ishape(jl), jl, icnt(jl)
1194               IF( ishape(jl) /= icnt(jl) )   CALL ctl_stop( TRIM(clinfo), 'size(pv_r'//clrankpv//TRIM(ctmp1)//TRIM(ctmp2) )
1195            END DO
1196
1197         ENDIF
1198
1199         ! read the data
1200         !-     
1201         IF( idvar > 0 .AND. istop == nstop ) THEN   ! no additional errors until this point...
1202            !
1203         ! find the right index of the array to be read
1204! JMM + SM: ugly patch before getting the new version of lib_mpp)
1205!         IF( idom /= jpdom_unknown ) THEN   ;   ix1 = nldi   ;   ix2 = nlei      ;   iy1 = nldj   ;   iy2 = nlej
1206!         ELSE                               ;   ix1 = 1      ;   ix2 = icnt(1)   ;   iy1 = 1      ;   iy2 = icnt(2)
1207!         ENDIF
1208            IF( llnoov ) THEN
1209               IF( idom /= jpdom_unknown ) THEN   ;   ix1 = nldi   ;   ix2 = nlei      ;   iy1 = nldj   ;   iy2 = nlej
1210               ELSE                               ;   ix1 = 1      ;   ix2 = icnt(1)   ;   iy1 = 1      ;   iy2 = icnt(2)
1211               ENDIF
1212            ELSE
1213               IF( idom /= jpdom_unknown ) THEN   ;   ix1 = 1      ;   ix2 = nlci      ;   iy1 = 1      ;   iy2 = nlcj
1214               ELSE                               ;   ix1 = 1      ;   ix2 = icnt(1)   ;   iy1 = 1      ;   iy2 = icnt(2)
1215               ENDIF
1216            ENDIF
1217     
1218            CALL iom_nf90_get( kiomid, idvar, inbdim, istart, icnt, ix1, ix2, iy1, iy2, pv_r1d, pv_r2d, pv_r3d )
1219
1220            IF( istop == nstop ) THEN   ! no additional errors until this point...
1221               IF(lwp) WRITE(numout,"(10x,' read ',a,' (rec: ',i6,') in ',a,' ok')") TRIM(cdvar), itime, TRIM(iom_file(kiomid)%name)
1222             
1223               !--- overlap areas and extra hallows (mpp)
1224               IF(     PRESENT(pv_r2d) .AND. idom /= jpdom_unknown ) THEN
1225                  CALL lbc_lnk( 'iom', pv_r2d,'Z',-999.,'no0' )
1226               ELSEIF( PRESENT(pv_r3d) .AND. idom /= jpdom_unknown ) THEN
1227                  ! this if could be simplified with the new lbc_lnk that works with any size of the 3rd dimension
1228                  IF( icnt(3) == inlev ) THEN
1229                     CALL lbc_lnk( 'iom', pv_r3d,'Z',-999.,'no0' )
1230                  ELSE   ! put some arbitrary value (a call to lbc_lnk will be done later...)
1231                     DO jj = nlcj+1, jpj   ;   pv_r3d(1:nlci, jj, :) = pv_r3d(1:nlci, nlej, :)   ;   END DO
1232                     DO ji = nlci+1, jpi   ;   pv_r3d(ji    , : , :) = pv_r3d(nlei  , :   , :)   ;   END DO
1233                  ENDIF
1234               ENDIF
1235               !
1236            ELSE
1237               ! return if istop == nstop is false
1238               RETURN
1239            ENDIF
1240         ELSE
1241            ! return if statment idvar > 0 .AND. istop == nstop is false
1242            RETURN
1243         ENDIF
1244         !
1245      ELSE        ! read using XIOS. Only if KEY_IOMPUT is defined
1246#if defined key_iomput
1247!would be good to be able to check which context is active and swap only if current is not restart
1248         CALL iom_swap( TRIM(crxios_context) ) 
1249         IF( PRESENT(pv_r3d) ) THEN
1250            pv_r3d(:, :, :) = 0.
1251            if(lwp) write(numout,*) 'XIOS RST READ (3D): ',trim(cdvar)
1252            CALL xios_recv_field( trim(cdvar), pv_r3d)
1253            IF(idom /= jpdom_unknown ) then
1254                CALL lbc_lnk( 'iom', pv_r3d,'Z',-999.,'no0' )
1255            ENDIF
1256         ELSEIF( PRESENT(pv_r2d) ) THEN
1257            pv_r2d(:, :) = 0.
1258            if(lwp) write(numout,*) 'XIOS RST READ (2D): ', trim(cdvar)
1259            CALL xios_recv_field( trim(cdvar), pv_r2d)
1260            IF(idom /= jpdom_unknown ) THEN
1261                CALL lbc_lnk('iom', pv_r2d,'Z',-999.,'no0')
1262            ENDIF
1263         ELSEIF( PRESENT(pv_r1d) ) THEN
1264            pv_r1d(:) = 0.
1265            if(lwp) write(numout,*) 'XIOS RST READ (1D): ', trim(cdvar)
1266            CALL xios_recv_field( trim(cdvar), pv_r1d)
1267         ENDIF
1268         CALL iom_swap( TRIM(cxios_context) )
1269#else
1270         istop = istop + 1 
1271         clinfo = 'Can not use XIOS in iom_get_123d, file: '//trim(clname)//', var:'//trim(cdvar)
1272#endif
1273      ENDIF
1274!some final adjustments
1275      ! C1D case : always call lbc_lnk to replicate the central value over the whole 3X3 domain
1276
1277      !--- Apply scale_factor and offset
1278      zscf = iom_file(kiomid)%scf(idvar)      ! scale factor
1279      zofs = iom_file(kiomid)%ofs(idvar)      ! offset
1280      IF(     PRESENT(pv_r1d) ) THEN
1281         IF( zscf /= 1. )   pv_r1d(:) = pv_r1d(:) * zscf 
1282         IF( zofs /= 0. )   pv_r1d(:) = pv_r1d(:) + zofs
1283      ELSEIF( PRESENT(pv_r2d) ) THEN
1284         IF( zscf /= 1.)   pv_r2d(:,:) = pv_r2d(:,:) * zscf
1285         IF( zofs /= 0.)   pv_r2d(:,:) = pv_r2d(:,:) + zofs
1286      ELSEIF( PRESENT(pv_r3d) ) THEN
1287         IF( zscf /= 1.)   pv_r3d(:,:,:) = pv_r3d(:,:,:) * zscf
1288         IF( zofs /= 0.)   pv_r3d(:,:,:) = pv_r3d(:,:,:) + zofs
1289      ENDIF
1290      !
1291   END SUBROUTINE iom_get_123d
1292
1293
1294   FUNCTION iom_getszuld ( kiomid ) 
1295      !!-----------------------------------------------------------------------
1296      !!                  ***  FUNCTION  iom_getszuld  ***
1297      !!
1298      !! ** Purpose : get the size of the unlimited dimension in a file
1299      !!              (return -1 if not found)
1300      !!-----------------------------------------------------------------------
1301      INTEGER, INTENT(in   ) ::   kiomid   ! file Identifier
1302      !
1303      INTEGER                ::   iom_getszuld
1304      !!-----------------------------------------------------------------------
1305      iom_getszuld = -1
1306      IF( kiomid > 0 ) THEN
1307         IF( iom_file(kiomid)%iduld > 0 )   iom_getszuld = iom_file(kiomid)%lenuld
1308      ENDIF
1309   END FUNCTION iom_getszuld
1310   
1311
1312   !!----------------------------------------------------------------------
1313   !!                   INTERFACE iom_chkatt
1314   !!----------------------------------------------------------------------
1315   SUBROUTINE iom_chkatt( kiomid, cdatt, llok, ksize, cdvar )
1316      INTEGER         , INTENT(in   )                 ::   kiomid    ! Identifier of the file
1317      CHARACTER(len=*), INTENT(in   )                 ::   cdatt     ! Name of the attribute
1318      LOGICAL         , INTENT(  out)                 ::   llok      ! Error code
1319      INTEGER         , INTENT(  out), OPTIONAL       ::   ksize     ! Size of the attribute array
1320      CHARACTER(len=*), INTENT(in   ), OPTIONAL       ::   cdvar     ! Name of the variable
1321      !
1322      IF( kiomid > 0 ) THEN
1323         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_chkatt( kiomid, cdatt, llok, ksize=ksize, cdvar=cdvar )
1324      ENDIF
1325      !
1326   END SUBROUTINE iom_chkatt
1327
1328   !!----------------------------------------------------------------------
1329   !!                   INTERFACE iom_getatt
1330   !!----------------------------------------------------------------------
1331   SUBROUTINE iom_g0d_iatt( kiomid, cdatt, katt0d, cdvar )
1332      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1333      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1334      INTEGER               , INTENT(  out)           ::   katt0d    ! read field
1335      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1336      !
1337      IF( kiomid > 0 ) THEN
1338         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  katt0d =  katt0d, cdvar=cdvar )
1339      ENDIF
1340   END SUBROUTINE iom_g0d_iatt
1341
1342   SUBROUTINE iom_g1d_iatt( kiomid, cdatt, katt1d, cdvar )
1343      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1344      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1345      INTEGER, DIMENSION(:) , INTENT(  out)           ::   katt1d    ! read field
1346      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1347      !
1348      IF( kiomid > 0 ) THEN
1349         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  katt1d =  katt1d, cdvar=cdvar )
1350      ENDIF
1351   END SUBROUTINE iom_g1d_iatt
1352
1353   SUBROUTINE iom_g0d_ratt( kiomid, cdatt, patt0d, cdvar )
1354      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1355      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1356      REAL(wp)              , INTENT(  out)           ::   patt0d    ! read field
1357      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1358      !
1359      IF( kiomid > 0 ) THEN
1360         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  patt0d =  patt0d, cdvar=cdvar )
1361      ENDIF
1362   END SUBROUTINE iom_g0d_ratt
1363
1364   SUBROUTINE iom_g1d_ratt( kiomid, cdatt, patt1d, cdvar )
1365      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1366      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1367      REAL(wp), DIMENSION(:), INTENT(  out)           ::   patt1d    ! read field
1368      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1369      !
1370      IF( kiomid > 0 ) THEN
1371         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt,  patt1d =  patt1d, cdvar=cdvar )
1372      ENDIF
1373   END SUBROUTINE iom_g1d_ratt
1374   
1375   SUBROUTINE iom_g0d_catt( kiomid, cdatt, cdatt0d, cdvar )
1376      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1377      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1378      CHARACTER(len=*)      , INTENT(  out)           ::   cdatt0d   ! read field
1379      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1380      !
1381      IF( kiomid > 0 ) THEN
1382         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_getatt( kiomid, cdatt, cdatt0d = cdatt0d, cdvar=cdvar )
1383      ENDIF
1384   END SUBROUTINE iom_g0d_catt
1385
1386
1387   !!----------------------------------------------------------------------
1388   !!                   INTERFACE iom_putatt
1389   !!----------------------------------------------------------------------
1390   SUBROUTINE iom_p0d_iatt( kiomid, cdatt, katt0d, cdvar )
1391      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1392      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1393      INTEGER               , INTENT(in   )           ::   katt0d    ! written field
1394      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1395      !
1396      IF( kiomid > 0 ) THEN
1397         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  katt0d =  katt0d, cdvar=cdvar )
1398      ENDIF
1399   END SUBROUTINE iom_p0d_iatt
1400
1401   SUBROUTINE iom_p1d_iatt( kiomid, cdatt, katt1d, cdvar )
1402      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1403      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1404      INTEGER, DIMENSION(:) , INTENT(in   )           ::   katt1d    ! written field
1405      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1406      !
1407      IF( kiomid > 0 ) THEN
1408         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  katt1d =  katt1d, cdvar=cdvar )
1409      ENDIF
1410   END SUBROUTINE iom_p1d_iatt
1411
1412   SUBROUTINE iom_p0d_ratt( kiomid, cdatt, patt0d, cdvar )
1413      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1414      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1415      REAL(wp)              , INTENT(in   )           ::   patt0d    ! written field
1416      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1417      !
1418      IF( kiomid > 0 ) THEN
1419         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  patt0d =  patt0d, cdvar=cdvar )
1420      ENDIF
1421   END SUBROUTINE iom_p0d_ratt
1422
1423   SUBROUTINE iom_p1d_ratt( kiomid, cdatt, patt1d, cdvar )
1424      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1425      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1426      REAL(wp), DIMENSION(:), INTENT(in   )           ::   patt1d    ! written field
1427      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1428      !
1429      IF( kiomid > 0 ) THEN
1430         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt,  patt1d =  patt1d, cdvar=cdvar )
1431      ENDIF
1432   END SUBROUTINE iom_p1d_ratt
1433   
1434   SUBROUTINE iom_p0d_catt( kiomid, cdatt, cdatt0d, cdvar )
1435      INTEGER               , INTENT(in   )           ::   kiomid    ! Identifier of the file
1436      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt     ! Name of the attribute
1437      CHARACTER(len=*)      , INTENT(in   )           ::   cdatt0d   ! written field
1438      CHARACTER(len=*)      , INTENT(in   ), OPTIONAL ::   cdvar     ! Name of the variable
1439      !
1440      IF( kiomid > 0 ) THEN
1441         IF( iom_file(kiomid)%nfid > 0 )   CALL iom_nf90_putatt( kiomid, cdatt, cdatt0d = cdatt0d, cdvar=cdvar )
1442      ENDIF
1443   END SUBROUTINE iom_p0d_catt
1444
1445
1446   !!----------------------------------------------------------------------
1447   !!                   INTERFACE iom_rstput
1448   !!----------------------------------------------------------------------
1449   SUBROUTINE iom_rp0d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1450      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1451      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1452      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1453      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1454      REAL(wp)        , INTENT(in)                         ::   pvar     ! written field
1455      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1456      LOGICAL, OPTIONAL :: ldxios   ! xios write flag
1457      LOGICAL :: llx                ! local xios write flag
1458      INTEGER :: ivid   ! variable id
1459
1460      llx = .FALSE.
1461      IF(PRESENT(ldxios)) llx = ldxios
1462      IF( llx ) THEN
1463#ifdef key_iomput
1464      IF( kt == kwrite ) THEN
1465          IF(lwp) write(numout,*) 'RESTART: write (XIOS 0D) ',trim(cdvar)
1466          CALL xios_send_field(trim(cdvar), pvar)
1467      ENDIF
1468#endif
1469      ELSE
1470         IF( kiomid > 0 ) THEN
1471            IF( iom_file(kiomid)%nfid > 0 ) THEN
1472               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1473               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r0d = pvar )
1474            ENDIF
1475         ENDIF
1476      ENDIF
1477   END SUBROUTINE iom_rp0d
1478
1479   SUBROUTINE iom_rp1d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1480      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1481      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1482      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1483      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1484      REAL(wp)        , INTENT(in), DIMENSION(          :) ::   pvar     ! written field
1485      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1486      LOGICAL, OPTIONAL                                    ::   ldxios   ! xios write flag
1487      LOGICAL :: llx                ! local xios write flag
1488      INTEGER :: ivid   ! variable id
1489
1490      llx = .FALSE.
1491      IF(PRESENT(ldxios)) llx = ldxios
1492      IF( llx ) THEN
1493#ifdef key_iomput
1494      IF( kt == kwrite ) THEN
1495         IF(lwp) write(numout,*) 'RESTART: write (XIOS 1D) ',trim(cdvar)
1496         CALL xios_send_field(trim(cdvar), pvar)
1497      ENDIF
1498#endif
1499      ELSE
1500         IF( kiomid > 0 ) THEN
1501            IF( iom_file(kiomid)%nfid > 0 ) THEN
1502               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1503               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r1d = pvar )
1504            ENDIF
1505         ENDIF
1506      ENDIF
1507   END SUBROUTINE iom_rp1d
1508
1509   SUBROUTINE iom_rp2d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1510      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1511      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1512      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1513      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1514      REAL(wp)        , INTENT(in), DIMENSION(:,    :    ) ::   pvar     ! written field
1515      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1516      LOGICAL, OPTIONAL :: ldxios   ! xios write flag
1517      LOGICAL :: llx
1518      INTEGER :: ivid   ! variable id
1519
1520      llx = .FALSE.
1521      IF(PRESENT(ldxios)) llx = ldxios
1522      IF( llx ) THEN
1523#ifdef key_iomput
1524      IF( kt == kwrite ) THEN
1525         IF(lwp) write(numout,*) 'RESTART: write (XIOS 2D) ',trim(cdvar)
1526         CALL xios_send_field(trim(cdvar), pvar)
1527      ENDIF
1528#endif
1529      ELSE
1530         IF( kiomid > 0 ) THEN
1531            IF( iom_file(kiomid)%nfid > 0 ) THEN
1532               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1533               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r2d = pvar )
1534            ENDIF
1535         ENDIF
1536      ENDIF
1537   END SUBROUTINE iom_rp2d
1538
1539   SUBROUTINE iom_rp3d( kt, kwrite, kiomid, cdvar, pvar, ktype, ldxios )
1540      INTEGER         , INTENT(in)                         ::   kt       ! ocean time-step
1541      INTEGER         , INTENT(in)                         ::   kwrite   ! writing time-step
1542      INTEGER         , INTENT(in)                         ::   kiomid   ! Identifier of the file
1543      CHARACTER(len=*), INTENT(in)                         ::   cdvar    ! time axis name
1544      REAL(wp)        , INTENT(in),       DIMENSION(:,:,:) ::   pvar     ! written field
1545      INTEGER         , INTENT(in), OPTIONAL               ::   ktype    ! variable external type
1546      LOGICAL, OPTIONAL :: ldxios   ! xios write flag
1547      LOGICAL :: llx                 ! local xios write flag
1548      INTEGER :: ivid   ! variable id
1549
1550      llx = .FALSE.
1551      IF(PRESENT(ldxios)) llx = ldxios
1552      IF( llx ) THEN
1553#ifdef key_iomput
1554      IF( kt == kwrite ) THEN
1555         IF(lwp) write(numout,*) 'RESTART: write (XIOS 3D) ',trim(cdvar)
1556         CALL xios_send_field(trim(cdvar), pvar)
1557      ENDIF
1558#endif
1559      ELSE
1560         IF( kiomid > 0 ) THEN
1561            IF( iom_file(kiomid)%nfid > 0 ) THEN
1562               ivid = iom_varid( kiomid, cdvar, ldstop = .FALSE. )
1563               CALL iom_nf90_rstput( kt, kwrite, kiomid, cdvar, ivid, ktype, pv_r3d = pvar )
1564            ENDIF
1565         ENDIF
1566      ENDIF
1567   END SUBROUTINE iom_rp3d
1568
1569
1570  SUBROUTINE iom_delay_rst( cdaction, cdcpnt, kncid )
1571      !!---------------------------------------------------------------------
1572      !!   Routine iom_delay_rst: used read/write restart related to mpp_delay
1573      !!
1574      !!---------------------------------------------------------------------
1575      CHARACTER(len=*), INTENT(in   ) ::   cdaction        !
1576      CHARACTER(len=*), INTENT(in   ) ::   cdcpnt
1577      INTEGER         , INTENT(in   ) ::   kncid
1578      !
1579      INTEGER  :: ji
1580      INTEGER  :: indim
1581      LOGICAL  :: llattexist
1582      REAL(wp), ALLOCATABLE, DIMENSION(:) ::   zreal1d
1583      !!---------------------------------------------------------------------
1584      !
1585      !                                      ===================================
1586      IF( TRIM(cdaction) == 'READ' ) THEN   ! read restart related to mpp_delay !
1587         !                                   ===================================
1588         DO ji = 1, nbdelay
1589            IF ( c_delaycpnt(ji) == cdcpnt ) THEN
1590               CALL iom_chkatt( kncid, 'DELAY_'//c_delaylist(ji), llattexist, indim )
1591               IF( llattexist )  THEN
1592                  ALLOCATE( todelay(ji)%z1d(indim) )
1593                  CALL iom_getatt( kncid, 'DELAY_'//c_delaylist(ji), todelay(ji)%z1d(:) )
1594                  ndelayid(ji) = 0   ! set to 0 to specify that the value was read in the restart
1595               ENDIF
1596           ENDIF
1597         END DO
1598         !                                   ====================================
1599      ELSE                                  ! write restart related to mpp_delay !
1600         !                                   ====================================
1601         DO ji = 1, nbdelay   ! save only ocean delayed global communication variables
1602            IF ( c_delaycpnt(ji) == cdcpnt ) THEN
1603               IF( ASSOCIATED(todelay(ji)%z1d) ) THEN
1604                  CALL mpp_delay_rcv(ji)   ! make sure %z1d is received
1605                  CALL iom_putatt( kncid, 'DELAY_'//c_delaylist(ji), todelay(ji)%z1d(:) )
1606               ENDIF
1607            ENDIF
1608         END DO
1609         !
1610      ENDIF
1611     
1612   END SUBROUTINE iom_delay_rst
1613 
1614   
1615
1616   !!----------------------------------------------------------------------
1617   !!                   INTERFACE iom_put
1618   !!----------------------------------------------------------------------
1619   SUBROUTINE iom_p0d( cdname, pfield0d )
1620      CHARACTER(LEN=*), INTENT(in) ::   cdname
1621      REAL(wp)        , INTENT(in) ::   pfield0d
1622      REAL(wp)        , DIMENSION(jpi,jpj) ::   zz     ! masson
1623#if defined key_iomput
1624      zz(:,:)=pfield0d
1625      CALL xios_send_field(cdname, zz)
1626      !CALL xios_send_field(cdname, (/pfield0d/))
1627#else
1628      IF( .FALSE. )   WRITE(numout,*) cdname, pfield0d   ! useless test to avoid compilation warnings
1629#endif
1630   END SUBROUTINE iom_p0d
1631
1632   SUBROUTINE iom_p1d( cdname, pfield1d )
1633      CHARACTER(LEN=*)          , INTENT(in) ::   cdname
1634      REAL(wp),     DIMENSION(:), INTENT(in) ::   pfield1d
1635#if defined key_iomput
1636      CALL xios_send_field( cdname, RESHAPE( (/pfield1d/), (/1,1,SIZE(pfield1d)/) ) )
1637#else
1638      IF( .FALSE. )   WRITE(numout,*) cdname, pfield1d   ! useless test to avoid compilation warnings
1639#endif
1640   END SUBROUTINE iom_p1d
1641
1642   SUBROUTINE iom_p2d( cdname, pfield2d )
1643      CHARACTER(LEN=*)            , INTENT(in) ::   cdname
1644      REAL(wp),     DIMENSION(:,:), INTENT(in) ::   pfield2d
1645#if defined key_iomput
1646      CALL xios_send_field(cdname, pfield2d)
1647#else
1648      IF( .FALSE. )   WRITE(numout,*) cdname, pfield2d   ! useless test to avoid compilation warnings
1649#endif
1650   END SUBROUTINE iom_p2d
1651
1652   SUBROUTINE iom_p3d( cdname, pfield3d )
1653      CHARACTER(LEN=*)                , INTENT(in) ::   cdname
1654      REAL(wp),       DIMENSION(:,:,:), INTENT(in) ::   pfield3d
1655#if defined key_iomput
1656      CALL xios_send_field( cdname, pfield3d )
1657#else
1658      IF( .FALSE. )   WRITE(numout,*) cdname, pfield3d   ! useless test to avoid compilation warnings
1659#endif
1660   END SUBROUTINE iom_p3d
1661
1662#if defined key_iomput
1663   !!----------------------------------------------------------------------
1664   !!   'key_iomput'                                         XIOS interface
1665   !!----------------------------------------------------------------------
1666
1667   SUBROUTINE iom_set_domain_attr( cdid, ni_glo, nj_glo, ibegin, jbegin, ni, nj,                                               &
1668      &                                    data_dim, data_ibegin, data_ni, data_jbegin, data_nj, lonvalue, latvalue, mask,     &
1669      &                                    nvertex, bounds_lon, bounds_lat, area )
1670      !!----------------------------------------------------------------------
1671      !!----------------------------------------------------------------------
1672      CHARACTER(LEN=*)                  , INTENT(in) ::   cdid
1673      INTEGER                 , OPTIONAL, INTENT(in) ::   ni_glo, nj_glo, ibegin, jbegin, ni, nj
1674      INTEGER                 , OPTIONAL, INTENT(in) ::   data_dim, data_ibegin, data_ni, data_jbegin, data_nj
1675      INTEGER                 , OPTIONAL, INTENT(in) ::   nvertex
1676      REAL(wp), DIMENSION(:)  , OPTIONAL, INTENT(in) ::   lonvalue, latvalue
1677      REAL(wp), DIMENSION(:,:), OPTIONAL, INTENT(in) ::   bounds_lon, bounds_lat, area
1678      LOGICAL , DIMENSION(:)  , OPTIONAL, INTENT(in) ::   mask
1679      !!----------------------------------------------------------------------
1680      !
1681      IF( xios_is_valid_domain     (cdid) ) THEN
1682         CALL xios_set_domain_attr     ( cdid, ni_glo=ni_glo, nj_glo=nj_glo, ibegin=ibegin, jbegin=jbegin, ni=ni, nj=nj,   &
1683            &    data_dim=data_dim, data_ibegin=data_ibegin, data_ni=data_ni, data_jbegin=data_jbegin, data_nj=data_nj ,   &
1684            &    lonvalue_1D=lonvalue, latvalue_1D=latvalue, mask_1D=mask, nvertex=nvertex, bounds_lon_1D=bounds_lon,      &
1685            &    bounds_lat_1D=bounds_lat, area=area, type='curvilinear')
1686      ENDIF
1687      IF( xios_is_valid_domaingroup(cdid) ) THEN
1688         CALL xios_set_domaingroup_attr( cdid, ni_glo=ni_glo, nj_glo=nj_glo, ibegin=ibegin, jbegin=jbegin, ni=ni, nj=nj,   &
1689            &    data_dim=data_dim, data_ibegin=data_ibegin, data_ni=data_ni, data_jbegin=data_jbegin, data_nj=data_nj ,   &
1690            &    lonvalue_1D=lonvalue, latvalue_1D=latvalue, mask_1D=mask, nvertex=nvertex, bounds_lon_1D=bounds_lon,      &
1691            &    bounds_lat_1D=bounds_lat, area=area, type='curvilinear' )
1692      ENDIF
1693      !
1694      CALL xios_solve_inheritance()
1695      !
1696   END SUBROUTINE iom_set_domain_attr
1697
1698
1699   SUBROUTINE iom_set_zoom_domain_attr( cdid, ibegin, jbegin, ni, nj )
1700      !!----------------------------------------------------------------------
1701      !!----------------------------------------------------------------------
1702      CHARACTER(LEN=*), INTENT(in) ::   cdid
1703      INTEGER         , INTENT(in) ::   ibegin, jbegin, ni, nj
1704      !
1705      TYPE(xios_gridgroup) :: gridgroup_hdl
1706      TYPE(xios_grid)      :: grid_hdl
1707      TYPE(xios_domain)    :: domain_hdl 
1708      TYPE(xios_axis)      :: axis_hdl 
1709      CHARACTER(LEN=64)    :: cldomrefid   ! domain_ref name
1710      CHARACTER(len=1)     :: cl1          ! last character of this name
1711      !!----------------------------------------------------------------------
1712      !
1713      IF( xios_is_valid_zoom_domain(cdid) ) THEN
1714         ! define the zoom_domain attributs
1715         CALL xios_set_zoom_domain_attr( cdid, ibegin=ibegin, jbegin=jbegin, ni=ni, nj=nj )
1716         ! define a new 2D grid with this new domain
1717         CALL xios_get_handle("grid_definition", gridgroup_hdl )
1718         CALL xios_add_child(gridgroup_hdl, grid_hdl, TRIM(cdid)//'_2D' )   ! add a new 2D grid to grid_definition
1719         CALL xios_add_child(grid_hdl, domain_hdl, TRIM(cdid) )             ! add its domain
1720         ! define a new 3D grid with this new domain
1721         CALL xios_add_child(gridgroup_hdl, grid_hdl, TRIM(cdid)//'_3D' )   ! add a new 3D grid to grid_definition
1722         CALL xios_add_child(grid_hdl, domain_hdl, TRIM(cdid) )             ! add its domain
1723         ! vertical axis
1724         cl1 = cdid(LEN_TRIM(cdid):)                                        ! last letter of cdid
1725         cl1 = CHAR(ICHAR(cl1)+32)                                          ! from upper to lower case
1726         CALL xios_add_child(grid_hdl, axis_hdl, 'depth'//cl1)              ! add its axis
1727      ENDIF
1728      !     
1729   END SUBROUTINE iom_set_zoom_domain_attr
1730
1731
1732   SUBROUTINE iom_set_axis_attr( cdid, paxis, bounds )
1733      !!----------------------------------------------------------------------
1734      !!----------------------------------------------------------------------
1735      CHARACTER(LEN=*)      , INTENT(in) ::   cdid
1736      REAL(wp), DIMENSION(:)  , OPTIONAL, INTENT(in) ::   paxis
1737      REAL(wp), DIMENSION(:,:), OPTIONAL, INTENT(in) ::   bounds
1738      !!----------------------------------------------------------------------
1739      IF( PRESENT(paxis) ) THEN
1740         IF( xios_is_valid_axis     (cdid) )   CALL xios_set_axis_attr     ( cdid, n_glo=SIZE(paxis), value=paxis )
1741         IF( xios_is_valid_axisgroup(cdid) )   CALL xios_set_axisgroup_attr( cdid, n_glo=SIZE(paxis), value=paxis )
1742      ENDIF
1743      IF( xios_is_valid_axis     (cdid) )   CALL xios_set_axis_attr     ( cdid, bounds=bounds )
1744      IF( xios_is_valid_axisgroup(cdid) )   CALL xios_set_axisgroup_attr( cdid, bounds=bounds )
1745      CALL xios_solve_inheritance()
1746   END SUBROUTINE iom_set_axis_attr
1747
1748
1749   SUBROUTINE iom_set_field_attr( cdid, freq_op, freq_offset )
1750      !!----------------------------------------------------------------------
1751      !!----------------------------------------------------------------------
1752      CHARACTER(LEN=*)             , INTENT(in) ::   cdid
1753      TYPE(xios_duration), OPTIONAL, INTENT(in) ::   freq_op
1754      TYPE(xios_duration), OPTIONAL, INTENT(in) ::   freq_offset
1755      !!----------------------------------------------------------------------
1756      IF( xios_is_valid_field     (cdid) )   CALL xios_set_field_attr     ( cdid, freq_op=freq_op, freq_offset=freq_offset )
1757      IF( xios_is_valid_fieldgroup(cdid) )   CALL xios_set_fieldgroup_attr( cdid, freq_op=freq_op, freq_offset=freq_offset )
1758      CALL xios_solve_inheritance()
1759   END SUBROUTINE iom_set_field_attr
1760
1761
1762   SUBROUTINE iom_set_file_attr( cdid, name, name_suffix )
1763      !!----------------------------------------------------------------------
1764      !!----------------------------------------------------------------------
1765      CHARACTER(LEN=*)          , INTENT(in) ::   cdid
1766      CHARACTER(LEN=*),OPTIONAL , INTENT(in) ::   name, name_suffix
1767      !!----------------------------------------------------------------------
1768      IF( xios_is_valid_file     (cdid) )   CALL xios_set_file_attr     ( cdid, name=name, name_suffix=name_suffix )
1769      IF( xios_is_valid_filegroup(cdid) )   CALL xios_set_filegroup_attr( cdid, name=name, name_suffix=name_suffix )
1770      CALL xios_solve_inheritance()
1771   END SUBROUTINE iom_set_file_attr
1772
1773
1774   SUBROUTINE iom_get_file_attr( cdid, name, name_suffix, output_freq )
1775      !!----------------------------------------------------------------------
1776      !!----------------------------------------------------------------------
1777      CHARACTER(LEN=*)          , INTENT(in ) ::   cdid
1778      CHARACTER(LEN=*),OPTIONAL , INTENT(out) ::   name, name_suffix
1779      TYPE(xios_duration), OPTIONAL , INTENT(out) :: output_freq
1780      LOGICAL                                 ::   llexist1,llexist2,llexist3
1781      !---------------------------------------------------------------------
1782      IF( PRESENT( name        ) )   name = ''          ! default values
1783      IF( PRESENT( name_suffix ) )   name_suffix = ''
1784      IF( PRESENT( output_freq ) )   output_freq = xios_duration(0,0,0,0,0,0)
1785      IF( xios_is_valid_file     (cdid) ) THEN
1786         CALL xios_solve_inheritance()
1787         CALL xios_is_defined_file_attr     ( cdid, name = llexist1, name_suffix = llexist2, output_freq = llexist3)
1788         IF(llexist1)   CALL xios_get_file_attr     ( cdid, name = name )
1789         IF(llexist2)   CALL xios_get_file_attr     ( cdid, name_suffix = name_suffix )
1790         IF(llexist3)   CALL xios_get_file_attr     ( cdid, output_freq = output_freq )
1791      ENDIF
1792      IF( xios_is_valid_filegroup(cdid) ) THEN
1793         CALL xios_solve_inheritance()
1794         CALL xios_is_defined_filegroup_attr( cdid, name = llexist1, name_suffix = llexist2, output_freq = llexist3)
1795         IF(llexist1)   CALL xios_get_filegroup_attr( cdid, name = name )
1796         IF(llexist2)   CALL xios_get_filegroup_attr( cdid, name_suffix = name_suffix )
1797         IF(llexist3)   CALL xios_get_filegroup_attr( cdid, output_freq = output_freq )
1798      ENDIF
1799   END SUBROUTINE iom_get_file_attr
1800
1801
1802   SUBROUTINE iom_set_grid_attr( cdid, mask )
1803      !!----------------------------------------------------------------------
1804      !!----------------------------------------------------------------------
1805      CHARACTER(LEN=*)                   , INTENT(in) ::   cdid
1806      LOGICAL, DIMENSION(:,:,:), OPTIONAL, INTENT(in) ::   mask
1807      !!----------------------------------------------------------------------
1808      IF( xios_is_valid_grid     (cdid) )   CALL xios_set_grid_attr     ( cdid, mask_3D=mask )
1809      IF( xios_is_valid_gridgroup(cdid) )   CALL xios_set_gridgroup_attr( cdid, mask_3D=mask )
1810      CALL xios_solve_inheritance()
1811   END SUBROUTINE iom_set_grid_attr
1812
1813   SUBROUTINE iom_setkt( kt, cdname )
1814      !!----------------------------------------------------------------------
1815      !!----------------------------------------------------------------------
1816      INTEGER         , INTENT(in) ::   kt 
1817      CHARACTER(LEN=*), INTENT(in) ::   cdname
1818      !!----------------------------------------------------------------------
1819      CALL iom_swap( cdname )   ! swap to cdname context
1820      CALL xios_update_calendar(kt)
1821      IF( cdname /= TRIM(cxios_context) )   CALL iom_swap( TRIM(cxios_context) )   ! return back to nemo context
1822   END SUBROUTINE iom_setkt
1823
1824   SUBROUTINE iom_context_finalize( cdname )
1825      !!----------------------------------------------------------------------
1826      !!----------------------------------------------------------------------
1827      CHARACTER(LEN=*), INTENT(in) :: cdname
1828      CHARACTER(LEN=120)           :: clname
1829      !!----------------------------------------------------------------------
1830      clname = cdname
1831      IF( TRIM(Agrif_CFixed()) .NE. '0' ) clname = TRIM(Agrif_CFixed())//"_"//clname 
1832      IF( xios_is_valid_context(clname) ) THEN
1833         CALL iom_swap( cdname )   ! swap to cdname context
1834         CALL xios_context_finalize() ! finalize the context
1835         IF( cdname /= TRIM(cxios_context) ) CALL iom_swap( TRIM(cxios_context) )   ! return back to nemo context
1836      ENDIF
1837      !
1838   END SUBROUTINE iom_context_finalize
1839
1840
1841   SUBROUTINE set_grid( cdgrd, plon, plat, ldxios, ldrxios )
1842      !!----------------------------------------------------------------------
1843      !!                     ***  ROUTINE set_grid  ***
1844      !!
1845      !! ** Purpose :   define horizontal grids
1846      !!----------------------------------------------------------------------
1847      CHARACTER(LEN=1)            , INTENT(in) ::   cdgrd
1848      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   plon
1849      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   plat
1850      !
1851      INTEGER  :: ni, nj
1852      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zmask
1853      LOGICAL, INTENT(IN) :: ldxios, ldrxios
1854      !!----------------------------------------------------------------------
1855      !
1856      ni = nlei-nldi+1
1857      nj = nlej-nldj+1
1858      !
1859      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)
1860      CALL iom_set_domain_attr("grid_"//cdgrd, data_dim=2, data_ibegin = 1-nldi, data_ni = jpi, data_jbegin = 1-nldj, data_nj = jpj)
1861!don't define lon and lat for restart reading context.
1862      IF ( .NOT.ldrxios ) &
1863         CALL iom_set_domain_attr("grid_"//cdgrd, lonvalue = RESHAPE(plon(nldi:nlei, nldj:nlej),(/ ni*nj /)),   &
1864         &                                     latvalue = RESHAPE(plat(nldi:nlei, nldj:nlej),(/ ni*nj /))) 
1865      !
1866      IF ( ln_mskland .AND. (.NOT.ldxios) ) THEN
1867         ! mask land points, keep values on coast line -> specific mask for U, V and W points
1868         SELECT CASE ( cdgrd )
1869         CASE('T')   ;   zmask(:,:,:)       = tmask(:,:,:)
1870         CASE('U')   ;   zmask(2:jpim1,:,:) = tmask(2:jpim1,:,:) + tmask(3:jpi,:,:)   ;   CALL lbc_lnk( 'iom', zmask, 'U', 1. )
1871         CASE('V')   ;   zmask(:,2:jpjm1,:) = tmask(:,2:jpjm1,:) + tmask(:,3:jpj,:)   ;   CALL lbc_lnk( 'iom', zmask, 'V', 1. )
1872         CASE('W')   ;   zmask(:,:,2:jpk  ) = tmask(:,:,1:jpkm1) + tmask(:,:,2:jpk)   ;   zmask(:,:,1) = tmask(:,:,1)
1873         END SELECT
1874         !
1875         CALL iom_set_domain_attr( "grid_"//cdgrd       , mask = RESHAPE(zmask(nldi:nlei,nldj:nlej,1),(/ni*nj    /)) /= 0. )
1876         CALL iom_set_grid_attr  ( "grid_"//cdgrd//"_3D", mask = RESHAPE(zmask(nldi:nlei,nldj:nlej,:),(/ni,nj,jpk/)) /= 0. )
1877      ENDIF
1878      !
1879   END SUBROUTINE set_grid
1880
1881
1882   SUBROUTINE set_grid_bounds( cdgrd, plon_cnr, plat_cnr, plon_pnt, plat_pnt )
1883      !!----------------------------------------------------------------------
1884      !!                   ***  ROUTINE set_grid_bounds  ***
1885      !!
1886      !! ** Purpose :   define horizontal grid corners
1887      !!
1888      !!----------------------------------------------------------------------
1889      CHARACTER(LEN=1)                      , INTENT(in) :: cdgrd
1890      REAL(wp), DIMENSION(jpi,jpj)          , INTENT(in) :: plon_cnr, plat_cnr  ! Lat/lon coord. of a contiguous vertex of cell (i,j)
1891      REAL(wp), DIMENSION(jpi,jpj), OPTIONAL, INTENT(in) :: plon_pnt, plat_pnt  ! Lat/lon coord. of the point of cell (i,j)
1892      !
1893      INTEGER :: ji, jj, jn, ni, nj
1894      INTEGER :: icnr, jcnr                                    ! Offset such that the vertex coordinate (i+icnr,j+jcnr)
1895      !                                                        ! represents the bottom-left corner of cell (i,j)
1896      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: z_bnds      ! Lat/lon coordinates of the vertices of cell (i,j)
1897      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z_fld       ! Working array to determine where to rotate cells
1898      REAL(wp), ALLOCATABLE, DIMENSION(:,:)     :: z_rot       ! Lat/lon working array for rotation of cells
1899      !!----------------------------------------------------------------------
1900      !
1901      ALLOCATE( z_bnds(4,jpi,jpj,2), z_fld(jpi,jpj), z_rot(4,2)  )
1902      !
1903      ! Offset of coordinate representing bottom-left corner
1904      SELECT CASE ( TRIM(cdgrd) )
1905      CASE ('T', 'W')   ;   icnr = -1   ;   jcnr = -1
1906      CASE ('U')        ;   icnr =  0   ;   jcnr = -1
1907      CASE ('V')        ;   icnr = -1   ;   jcnr =  0
1908      END SELECT
1909      !
1910      ni = nlei-nldi+1   ! Dimensions of subdomain interior
1911      nj = nlej-nldj+1
1912      !
1913      z_fld(:,:) = 1._wp
1914      CALL lbc_lnk( 'iom', z_fld, cdgrd, -1. )    ! Working array for location of northfold
1915      !
1916      ! Cell vertices that can be defined
1917      DO jj = 2, jpjm1
1918         DO ji = 2, jpim1
1919            z_bnds(1,ji,jj,1) = plat_cnr(ji+icnr,  jj+jcnr  ) ! Bottom-left
1920            z_bnds(2,ji,jj,1) = plat_cnr(ji+icnr+1,jj+jcnr  ) ! Bottom-right
1921            z_bnds(3,ji,jj,1) = plat_cnr(ji+icnr+1,jj+jcnr+1) ! Top-right
1922            z_bnds(4,ji,jj,1) = plat_cnr(ji+icnr,  jj+jcnr+1) ! Top-left
1923            z_bnds(1,ji,jj,2) = plon_cnr(ji+icnr,  jj+jcnr  ) ! Bottom-left
1924            z_bnds(2,ji,jj,2) = plon_cnr(ji+icnr+1,jj+jcnr  ) ! Bottom-right
1925            z_bnds(3,ji,jj,2) = plon_cnr(ji+icnr+1,jj+jcnr+1) ! Top-right
1926            z_bnds(4,ji,jj,2) = plon_cnr(ji+icnr,  jj+jcnr+1) ! Top-left
1927         END DO
1928      END DO
1929      !
1930      ! Cell vertices on boundries
1931      DO jn = 1, 4
1932         CALL lbc_lnk( 'iom', z_bnds(jn,:,:,1), cdgrd, 1., pval=999._wp )
1933         CALL lbc_lnk( 'iom', z_bnds(jn,:,:,2), cdgrd, 1., pval=999._wp )
1934      END DO
1935      !
1936      ! Zero-size cells at closed boundaries if cell points provided,
1937      ! otherwise they are closed cells with unrealistic bounds
1938      IF( PRESENT(plon_pnt) .AND. PRESENT(plat_pnt) ) THEN
1939         IF( (nbondi == -1 .OR. nbondi == 2) .AND. .NOT. (jperio == 1 .OR. jperio == 4 .OR. jperio == 6) ) THEN
1940            DO jn = 1, 4        ! (West or jpni = 1), closed E-W
1941               z_bnds(jn,1,:,1) = plat_pnt(1,:)  ;  z_bnds(jn,1,:,2) = plon_pnt(1,:)
1942            END DO
1943         ENDIF
1944         IF( (nbondi == 1 .OR. nbondi == 2) .AND. .NOT. (jperio == 1 .OR. jperio == 4 .OR. jperio == 6) ) THEN
1945            DO jn = 1, 4        ! (East or jpni = 1), closed E-W
1946               z_bnds(jn,nlci,:,1) = plat_pnt(nlci,:)  ;  z_bnds(jn,nlci,:,2) = plon_pnt(nlci,:)
1947            END DO
1948         ENDIF
1949         IF( nbondj == -1 .OR. (nbondj == 2 .AND. jperio /= 2) ) THEN
1950            DO jn = 1, 4        ! South or (jpnj = 1, not symmetric)
1951               z_bnds(jn,:,1,1) = plat_pnt(:,1)  ;  z_bnds(jn,:,1,2) = plon_pnt(:,1)
1952            END DO
1953         ENDIF
1954         IF( (nbondj == 1 .OR. nbondj == 2) .AND. jperio  < 3 ) THEN
1955            DO jn = 1, 4        ! (North or jpnj = 1), no north fold
1956               z_bnds(jn,:,nlcj,1) = plat_pnt(:,nlcj)  ;  z_bnds(jn,:,nlcj,2) = plon_pnt(:,nlcj)
1957            END DO
1958         ENDIF
1959      ENDIF
1960      !
1961      IF( (nbondj == 1 .OR. nbondj == 2) .AND. jperio >= 3 ) THEN    ! Rotate cells at the north fold
1962         DO jj = 1, jpj
1963            DO ji = 1, jpi
1964               IF( z_fld(ji,jj) == -1. ) THEN
1965                  z_rot(1,:) = z_bnds(3,ji,jj,:) ; z_rot(2,:) = z_bnds(4,ji,jj,:)
1966                  z_rot(3,:) = z_bnds(1,ji,jj,:) ; z_rot(4,:) = z_bnds(2,ji,jj,:)
1967                  z_bnds(:,ji,jj,:) = z_rot(:,:)
1968               ENDIF
1969            END DO
1970         END DO
1971      ELSE IF( nbondj == 2 .AND. jperio == 2 ) THEN                  ! Invert cells at the symmetric equator
1972         DO ji = 1, jpi
1973            z_rot(1:2,:) = z_bnds(3:4,ji,1,:)
1974            z_rot(3:4,:) = z_bnds(1:2,ji,1,:)
1975            z_bnds(:,ji,1,:) = z_rot(:,:)
1976         END DO
1977      ENDIF
1978      !
1979      CALL iom_set_domain_attr("grid_"//cdgrd, bounds_lat = RESHAPE(z_bnds(:,nldi:nlei,nldj:nlej,1),(/ 4,ni*nj /)),           &
1980          &                                    bounds_lon = RESHAPE(z_bnds(:,nldi:nlei,nldj:nlej,2),(/ 4,ni*nj /)), nvertex=4 )
1981      !
1982      DEALLOCATE( z_bnds, z_fld, z_rot ) 
1983      !
1984   END SUBROUTINE set_grid_bounds
1985
1986
1987   SUBROUTINE set_grid_znl( plat )
1988      !!----------------------------------------------------------------------
1989      !!                     ***  ROUTINE set_grid_znl  ***
1990      !!
1991      !! ** Purpose :   define grids for zonal mean
1992      !!
1993      !!----------------------------------------------------------------------
1994      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   plat
1995      !
1996      INTEGER  :: ni, nj, ix, iy
1997      REAL(wp), DIMENSION(:), ALLOCATABLE  ::   zlon
1998      !!----------------------------------------------------------------------
1999      !
2000      ni=nlei-nldi+1       ! define zonal mean domain (jpj*jpk)
2001      nj=nlej-nldj+1
2002      ALLOCATE( zlon(ni*nj) )       ;       zlon(:) = 0._wp
2003      !
2004      CALL dom_ngb( -168.53, 65.03, ix, iy, 'T' ) !  i-line that passes through Bering Strait: Reference latitude (used in plots)
2005!      CALL dom_ngb( 180., 90., ix, iy, 'T' ) !  i-line that passes near the North Pole : Reference latitude (used in plots)
2006      CALL iom_set_domain_attr("gznl", ni_glo=jpiglo, nj_glo=jpjglo, ibegin=nimpp+nldi-2, jbegin=njmpp+nldj-2, ni=ni, nj=nj)
2007      CALL iom_set_domain_attr("gznl", data_dim=2, data_ibegin = 1-nldi, data_ni = jpi, data_jbegin = 1-nldj, data_nj = jpj)
2008      CALL iom_set_domain_attr("gznl", lonvalue = zlon,   &
2009         &                             latvalue = RESHAPE(plat(nldi:nlei, nldj:nlej),(/ ni*nj /))) 
2010      CALL iom_set_zoom_domain_attr("znl_T", ibegin=ix-1, jbegin=0, ni=1, nj=jpjglo)
2011      CALL iom_set_zoom_domain_attr("znl_W", ibegin=ix-1, jbegin=0, ni=1, nj=jpjglo)
2012      !
2013      CALL iom_update_file_name('ptr')
2014      !
2015   END SUBROUTINE set_grid_znl
2016
2017
2018   SUBROUTINE set_scalar
2019      !!----------------------------------------------------------------------
2020      !!                     ***  ROUTINE set_scalar  ***
2021      !!
2022      !! ** Purpose :   define fake grids for scalar point
2023      !!
2024      !!----------------------------------------------------------------------
2025      REAL(wp), DIMENSION(1)   ::   zz = 1.
2026      !!----------------------------------------------------------------------
2027      !
2028      CALL iom_set_domain_attr('scalarpoint', ni_glo=jpnij, nj_glo=1, ibegin=narea-1, jbegin=0, ni=1, nj=1)
2029      CALL iom_set_domain_attr('scalarpoint', data_dim=2, data_ibegin = 1, data_ni = 1, data_jbegin = 1, data_nj = 1)
2030      !
2031      zz = REAL( narea, wp )
2032      CALL iom_set_domain_attr('scalarpoint', lonvalue=zz, latvalue=zz)
2033      !
2034   END SUBROUTINE set_scalar
2035
2036
2037   SUBROUTINE set_xmlatt
2038      !!----------------------------------------------------------------------
2039      !!                     ***  ROUTINE set_xmlatt  ***
2040      !!
2041      !! ** Purpose :   automatic definitions of some of the xml attributs...
2042      !!
2043      !!----------------------------------------------------------------------
2044      CHARACTER(len=1),DIMENSION( 3) ::   clgrd                    ! suffix name
2045      CHARACTER(len=256)             ::   clsuff                   ! suffix name
2046      CHARACTER(len=1)               ::   cl1                      ! 1 character
2047      CHARACTER(len=2)               ::   cl2                      ! 2 characters
2048      CHARACTER(len=3)               ::   cl3                      ! 3 characters
2049      INTEGER                        ::   ji, jg                   ! loop counters
2050      INTEGER                        ::   ix, iy                   ! i-,j- index
2051      REAL(wp)        ,DIMENSION(11) ::   zlontao                  ! longitudes of tao    moorings
2052      REAL(wp)        ,DIMENSION( 7) ::   zlattao                  ! latitudes  of tao    moorings
2053      REAL(wp)        ,DIMENSION( 4) ::   zlonrama                 ! longitudes of rama   moorings
2054      REAL(wp)        ,DIMENSION(11) ::   zlatrama                 ! latitudes  of rama   moorings
2055      REAL(wp)        ,DIMENSION( 3) ::   zlonpira                 ! longitudes of pirata moorings
2056      REAL(wp)        ,DIMENSION( 9) ::   zlatpira                 ! latitudes  of pirata moorings
2057      TYPE(xios_duration)            ::   f_op, f_of
2058      !!----------------------------------------------------------------------
2059      !
2060      ! frequency of the call of iom_put (attribut: freq_op)
2061      f_op%timestep = 1        ;  f_of%timestep =  0  ; CALL iom_set_field_attr('field_definition', freq_op=f_op, freq_offset=f_of)
2062      f_op%timestep = 2        ;  f_of%timestep =  0  ; CALL iom_set_field_attr('trendT_even'     , freq_op=f_op, freq_offset=f_of)
2063      f_op%timestep = 2        ;  f_of%timestep = -1  ; CALL iom_set_field_attr('trendT_odd'      , freq_op=f_op, freq_offset=f_of)
2064      f_op%timestep = nn_fsbc  ;  f_of%timestep =  0  ; CALL iom_set_field_attr('SBC'             , freq_op=f_op, freq_offset=f_of)
2065      f_op%timestep = nn_fsbc  ;  f_of%timestep =  0  ; CALL iom_set_field_attr('SBC_scalar'      , freq_op=f_op, freq_offset=f_of)
2066      f_op%timestep = nn_dttrc ;  f_of%timestep =  0  ; CALL iom_set_field_attr('ptrc_T'          , freq_op=f_op, freq_offset=f_of)
2067      f_op%timestep = nn_dttrc ;  f_of%timestep =  0  ; CALL iom_set_field_attr('diad_T'          , freq_op=f_op, freq_offset=f_of)
2068
2069      ! output file names (attribut: name)
2070      DO ji = 1, 9
2071         WRITE(cl1,'(i1)') ji 
2072         CALL iom_update_file_name('file'//cl1)
2073      END DO
2074      DO ji = 1, 99
2075         WRITE(cl2,'(i2.2)') ji 
2076         CALL iom_update_file_name('file'//cl2)
2077      END DO
2078      DO ji = 1, 999
2079         WRITE(cl3,'(i3.3)') ji 
2080         CALL iom_update_file_name('file'//cl3)
2081      END DO
2082
2083      ! Zooms...
2084      clgrd = (/ 'T', 'U', 'W' /) 
2085      DO jg = 1, SIZE(clgrd)                                                                   ! grid type
2086         cl1 = clgrd(jg)
2087         ! Equatorial section (attributs: jbegin, ni, name_suffix)
2088         CALL dom_ngb( 0., 0., ix, iy, cl1 )
2089         CALL iom_set_zoom_domain_attr('Eq'//cl1, ibegin=0, jbegin=iy-1, ni=jpiglo, nj=1 )
2090         CALL iom_get_file_attr   ('Eq'//cl1, name_suffix = clsuff             )
2091         CALL iom_set_file_attr   ('Eq'//cl1, name_suffix = TRIM(clsuff)//'_Eq')
2092         CALL iom_update_file_name('Eq'//cl1)
2093      END DO
2094      ! TAO moorings (attributs: ibegin, jbegin, name_suffix)
2095      zlontao = (/ 137.0, 147.0, 156.0, 165.0, -180.0, -170.0, -155.0, -140.0, -125.0, -110.0, -95.0 /)
2096      zlattao = (/  -8.0,  -5.0,  -2.0,   0.0,    2.0,    5.0,    8.0 /)
2097      CALL set_mooring( zlontao, zlattao )
2098      ! RAMA moorings (attributs: ibegin, jbegin, name_suffix)
2099      zlonrama = (/  55.0,  67.0, 80.5, 90.0 /)
2100      zlatrama = (/ -16.0, -12.0, -8.0, -4.0, -1.5, 0.0, 1.5, 4.0, 8.0, 12.0, 15.0 /)
2101      CALL set_mooring( zlonrama, zlatrama )
2102      ! PIRATA moorings (attributs: ibegin, jbegin, name_suffix)
2103      zlonpira = (/ -38.0, -23.0, -10.0 /)
2104      zlatpira = (/ -19.0, -14.0,  -8.0, 0.0, 4.0, 8.0, 12.0, 15.0, 20.0 /)
2105      CALL set_mooring( zlonpira, zlatpira )
2106      !
2107   END SUBROUTINE set_xmlatt
2108
2109
2110   SUBROUTINE set_mooring( plon, plat )
2111      !!----------------------------------------------------------------------
2112      !!                     ***  ROUTINE set_mooring  ***
2113      !!
2114      !! ** Purpose :   automatic definitions of moorings xml attributs...
2115      !!
2116      !!----------------------------------------------------------------------
2117      REAL(wp), DIMENSION(:), INTENT(in) ::   plon, plat   ! longitudes/latitudes oft the mooring
2118      !
2119!!$      CHARACTER(len=1),DIMENSION(4) ::   clgrd = (/ 'T', 'U', 'V', 'W' /)   ! suffix name
2120      CHARACTER(len=1),DIMENSION(1) ::   clgrd = (/ 'T' /)        ! suffix name
2121      CHARACTER(len=256)            ::   clname                   ! file name
2122      CHARACTER(len=256)            ::   clsuff                   ! suffix name
2123      CHARACTER(len=1)              ::   cl1                      ! 1 character
2124      CHARACTER(len=6)              ::   clon,clat                ! name of longitude, latitude
2125      INTEGER                       ::   ji, jj, jg               ! loop counters
2126      INTEGER                       ::   ix, iy                   ! i-,j- index
2127      REAL(wp)                      ::   zlon, zlat
2128      !!----------------------------------------------------------------------
2129      DO jg = 1, SIZE(clgrd)
2130         cl1 = clgrd(jg)
2131         DO ji = 1, SIZE(plon)
2132            DO jj = 1, SIZE(plat)
2133               zlon = plon(ji)
2134               zlat = plat(jj)
2135               ! modifications for RAMA moorings
2136               IF( zlon ==  67. .AND. zlat ==  15. )   zlon =  65.
2137               IF( zlon ==  90. .AND. zlat <=  -4. )   zlon =  95.
2138               IF( zlon ==  95. .AND. zlat ==  -4. )   zlat =  -5.
2139               ! modifications for PIRATA moorings
2140               IF( zlon == -38. .AND. zlat == -19. )   zlon = -34.
2141               IF( zlon == -38. .AND. zlat == -14. )   zlon = -32.
2142               IF( zlon == -38. .AND. zlat ==  -8. )   zlon = -30.
2143               IF( zlon == -38. .AND. zlat ==   0. )   zlon = -35.
2144               IF( zlon == -23. .AND. zlat ==  20. )   zlat =  21.
2145               IF( zlon == -10. .AND. zlat == -14. )   zlat = -10.
2146               IF( zlon == -10. .AND. zlat ==  -8. )   zlat =  -6.
2147               IF( zlon == -10. .AND. zlat ==   4. ) THEN   ;   zlon = 0.   ;   zlat = 0.   ;   ENDIF
2148               CALL dom_ngb( zlon, zlat, ix, iy, cl1 )
2149               IF( zlon >= 0. ) THEN 
2150                  IF( zlon == REAL(NINT(zlon), wp) ) THEN   ;   WRITE(clon, '(i3,  a)') NINT( zlon), 'e'
2151                  ELSE                                      ;   WRITE(clon, '(f5.1,a)')       zlon , 'e'
2152                  ENDIF
2153               ELSE             
2154                  IF( zlon == REAL(NINT(zlon), wp) ) THEN   ;   WRITE(clon, '(i3,  a)') NINT(-zlon), 'w'
2155                  ELSE                                      ;   WRITE(clon, '(f5.1,a)')      -zlon , 'w'
2156                  ENDIF
2157               ENDIF
2158               IF( zlat >= 0. ) THEN 
2159                  IF( zlat == REAL(NINT(zlat), wp) ) THEN   ;   WRITE(clat, '(i2,  a)') NINT( zlat), 'n'
2160                  ELSE                                      ;   WRITE(clat, '(f4.1,a)')       zlat , 'n'
2161                  ENDIF
2162               ELSE             
2163                  IF( zlat == REAL(NINT(zlat), wp) ) THEN   ;   WRITE(clat, '(i2,  a)') NINT(-zlat), 's'
2164                  ELSE                                      ;   WRITE(clat, '(f4.1,a)')      -zlat , 's'
2165                  ENDIF
2166               ENDIF
2167               clname = TRIM(ADJUSTL(clat))//TRIM(ADJUSTL(clon))
2168               CALL iom_set_zoom_domain_attr(TRIM(clname)//cl1, ibegin= ix-1, jbegin= iy-1, ni=1, nj=1)
2169
2170               CALL iom_get_file_attr   (TRIM(clname)//cl1, name_suffix = clsuff                         )
2171               CALL iom_set_file_attr   (TRIM(clname)//cl1, name_suffix = TRIM(clsuff)//'_'//TRIM(clname))
2172               CALL iom_update_file_name(TRIM(clname)//cl1)
2173            END DO
2174         END DO
2175      END DO
2176     
2177   END SUBROUTINE set_mooring
2178
2179   
2180   SUBROUTINE iom_update_file_name( cdid )
2181      !!----------------------------------------------------------------------
2182      !!                     ***  ROUTINE iom_update_file_name  ***
2183      !!
2184      !! ** Purpose :   
2185      !!
2186      !!----------------------------------------------------------------------
2187      CHARACTER(LEN=*)          , INTENT(in) ::   cdid
2188      !
2189      CHARACTER(LEN=256) ::   clname
2190      CHARACTER(LEN=20)  ::   clfreq
2191      CHARACTER(LEN=20)  ::   cldate
2192      INTEGER            ::   idx
2193      INTEGER            ::   jn
2194      INTEGER            ::   itrlen
2195      INTEGER            ::   iyear, imonth, iday, isec
2196      REAL(wp)           ::   zsec
2197      LOGICAL            ::   llexist
2198      TYPE(xios_duration)   ::   output_freq 
2199      !!----------------------------------------------------------------------
2200      !
2201      DO jn = 1, 2
2202         !
2203         output_freq = xios_duration(0,0,0,0,0,0)
2204         IF( jn == 1 )   CALL iom_get_file_attr( cdid, name        = clname, output_freq = output_freq )
2205         IF( jn == 2 )   CALL iom_get_file_attr( cdid, name_suffix = clname )
2206         !
2207         IF ( TRIM(clname) /= '' ) THEN 
2208            !
2209            idx = INDEX(clname,'@expname@') + INDEX(clname,'@EXPNAME@')
2210            DO WHILE ( idx /= 0 ) 
2211               clname = clname(1:idx-1)//TRIM(cexper)//clname(idx+9:LEN_TRIM(clname))
2212               idx = INDEX(clname,'@expname@') + INDEX(clname,'@EXPNAME@')
2213            END DO
2214            !
2215            idx = INDEX(clname,'@freq@') + INDEX(clname,'@FREQ@')
2216            DO WHILE ( idx /= 0 ) 
2217              IF ( output_freq%timestep /= 0) THEN
2218                  WRITE(clfreq,'(I18,A2)')INT(output_freq%timestep),'ts' 
2219                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2220              ELSE IF ( output_freq%second /= 0 ) THEN
2221                  WRITE(clfreq,'(I19,A1)')INT(output_freq%second),'s' 
2222                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2223              ELSE IF ( output_freq%minute /= 0 ) THEN
2224                  WRITE(clfreq,'(I18,A2)')INT(output_freq%minute),'mi' 
2225                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2226              ELSE IF ( output_freq%hour /= 0 ) THEN
2227                  WRITE(clfreq,'(I19,A1)')INT(output_freq%hour),'h' 
2228                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2229              ELSE IF ( output_freq%day /= 0 ) THEN
2230                  WRITE(clfreq,'(I19,A1)')INT(output_freq%day),'d' 
2231                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2232              ELSE IF ( output_freq%month /= 0 ) THEN   
2233                  WRITE(clfreq,'(I19,A1)')INT(output_freq%month),'m' 
2234                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2235              ELSE IF ( output_freq%year /= 0 ) THEN   
2236                  WRITE(clfreq,'(I19,A1)')INT(output_freq%year),'y' 
2237                  itrlen = LEN_TRIM(ADJUSTL(clfreq))
2238              ELSE
2239                  CALL ctl_stop('error in the name of file id '//TRIM(cdid),   &
2240                     & ' attribute output_freq is undefined -> cannot replace @freq@ in '//TRIM(clname) )
2241              ENDIF
2242              clname = clname(1:idx-1)//TRIM(ADJUSTL(clfreq))//clname(idx+6:LEN_TRIM(clname))
2243              idx = INDEX(clname,'@freq@') + INDEX(clname,'@FREQ@')
2244            END DO
2245            !
2246            idx = INDEX(clname,'@startdate@') + INDEX(clname,'@STARTDATE@')
2247            DO WHILE ( idx /= 0 ) 
2248               cldate = iom_sdate( fjulday - rdt / rday )
2249               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+11:LEN_TRIM(clname))
2250               idx = INDEX(clname,'@startdate@') + INDEX(clname,'@STARTDATE@')
2251            END DO
2252            !
2253            idx = INDEX(clname,'@startdatefull@') + INDEX(clname,'@STARTDATEFULL@')
2254            DO WHILE ( idx /= 0 ) 
2255               cldate = iom_sdate( fjulday - rdt / rday, ldfull = .TRUE. )
2256               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+15:LEN_TRIM(clname))
2257               idx = INDEX(clname,'@startdatefull@') + INDEX(clname,'@STARTDATEFULL@')
2258            END DO
2259            !
2260            idx = INDEX(clname,'@enddate@') + INDEX(clname,'@ENDDATE@')
2261            DO WHILE ( idx /= 0 ) 
2262               cldate = iom_sdate( fjulday + rdt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE. )
2263               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+9:LEN_TRIM(clname))
2264               idx = INDEX(clname,'@enddate@') + INDEX(clname,'@ENDDATE@')
2265            END DO
2266            !
2267            idx = INDEX(clname,'@enddatefull@') + INDEX(clname,'@ENDDATEFULL@')
2268            DO WHILE ( idx /= 0 ) 
2269               cldate = iom_sdate( fjulday + rdt / rday * REAL( nitend - nit000, wp ), ld24 = .TRUE., ldfull = .TRUE. )
2270               clname = clname(1:idx-1)//TRIM(cldate)//clname(idx+13:LEN_TRIM(clname))
2271               idx = INDEX(clname,'@enddatefull@') + INDEX(clname,'@ENDDATEFULL@')
2272            END DO
2273            !
2274            IF( jn == 1 .AND. TRIM(Agrif_CFixed()) /= '0' )   clname = TRIM(Agrif_CFixed())//"_"//TRIM(clname)
2275            IF( jn == 1 )   CALL iom_set_file_attr( cdid, name        = clname )
2276            IF( jn == 2 )   CALL iom_set_file_attr( cdid, name_suffix = clname )
2277            !
2278         ENDIF
2279         !
2280      END DO
2281      !
2282   END SUBROUTINE iom_update_file_name
2283
2284
2285   FUNCTION iom_sdate( pjday, ld24, ldfull )
2286      !!----------------------------------------------------------------------
2287      !!                     ***  ROUTINE iom_sdate  ***
2288      !!
2289      !! ** Purpose :   send back the date corresponding to the given julian day
2290      !!----------------------------------------------------------------------
2291      REAL(wp), INTENT(in   )           ::   pjday    ! julian day
2292      LOGICAL , INTENT(in   ), OPTIONAL ::   ld24     ! true to force 24:00 instead of 00:00
2293      LOGICAL , INTENT(in   ), OPTIONAL ::   ldfull   ! true to get the compleate date: yyyymmdd_hh:mm:ss
2294      !
2295      CHARACTER(LEN=20) ::   iom_sdate
2296      CHARACTER(LEN=50) ::   clfmt                         !  format used to write the date
2297      INTEGER           ::   iyear, imonth, iday, ihour, iminute, isec
2298      REAL(wp)          ::   zsec
2299      LOGICAL           ::   ll24, llfull
2300      !!----------------------------------------------------------------------
2301      !
2302      IF( PRESENT(ld24) ) THEN   ;   ll24 = ld24
2303      ELSE                       ;   ll24 = .FALSE.
2304      ENDIF
2305      !
2306      IF( PRESENT(ldfull) ) THEN   ;   llfull = ldfull
2307      ELSE                         ;   llfull = .FALSE.
2308      ENDIF
2309      !
2310      CALL ju2ymds( pjday, iyear, imonth, iday, zsec )
2311      isec = NINT(zsec)
2312      !
2313      IF ( ll24 .AND. isec == 0 ) THEN   ! 00:00 of the next day -> move to 24:00 of the current day
2314         CALL ju2ymds( pjday - 1., iyear, imonth, iday, zsec )
2315         isec = 86400
2316      ENDIF
2317      !
2318      IF( iyear < 10000 ) THEN   ;   clfmt = "i4.4,2i2.2"                ! format used to write the date
2319      ELSE                       ;   WRITE(clfmt, "('i',i1,',2i2.2')") INT(LOG10(REAL(iyear,wp))) + 1
2320      ENDIF
2321      !
2322!$AGRIF_DO_NOT_TREAT     
2323      ! needed in the conv
2324      IF( llfull ) THEN
2325         clfmt = TRIM(clfmt)//",'_',i2.2,':',i2.2,':',i2.2"
2326         ihour   = isec / 3600
2327         isec    = MOD(isec, 3600)
2328         iminute = isec / 60
2329         isec    = MOD(isec, 60)
2330         WRITE(iom_sdate, '('//TRIM(clfmt)//')') iyear, imonth, iday, ihour, iminute, isec    ! date of the end of run
2331      ELSE
2332         WRITE(iom_sdate, '('//TRIM(clfmt)//')') iyear, imonth, iday                          ! date of the end of run
2333      ENDIF
2334!$AGRIF_END_DO_NOT_TREAT     
2335      !
2336   END FUNCTION iom_sdate
2337
2338#else
2339   !!----------------------------------------------------------------------
2340   !!   NOT 'key_iomput'                               a few dummy routines
2341   !!----------------------------------------------------------------------
2342
2343   SUBROUTINE iom_setkt( kt, cdname )
2344      INTEGER         , INTENT(in)::   kt 
2345      CHARACTER(LEN=*), INTENT(in) ::   cdname
2346      IF( .FALSE. )   WRITE(numout,*) kt, cdname   ! useless test to avoid compilation warnings
2347   END SUBROUTINE iom_setkt
2348
2349   SUBROUTINE iom_context_finalize( cdname )
2350      CHARACTER(LEN=*), INTENT(in) ::   cdname
2351      IF( .FALSE. )   WRITE(numout,*)  cdname   ! useless test to avoid compilation warnings
2352   END SUBROUTINE iom_context_finalize
2353
2354#endif
2355
2356   LOGICAL FUNCTION iom_use( cdname )
2357      !!----------------------------------------------------------------------
2358      !!----------------------------------------------------------------------
2359      CHARACTER(LEN=*), INTENT(in) ::   cdname
2360      !!----------------------------------------------------------------------
2361#if defined key_iomput
2362      iom_use = xios_field_is_active( cdname )
2363#else
2364      iom_use = .FALSE.
2365#endif
2366   END FUNCTION iom_use
2367   
2368   !!======================================================================
2369END MODULE iom
Note: See TracBrowser for help on using the repository browser.