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

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

source: utils/tools_AGRIF_CMEMS_2020/DOMAINcfg/src/iom.F90 @ 10727

Last change on this file since 10727 was 10727, checked in by rblod, 5 years ago

new nesting tools (attempt) and brutal cleaning of DOMAINcfg, see ticket #2129

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