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

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

source: NEMO/trunk/src/OCE/IOM/iom.F90 @ 13214

Last change on this file since 13214 was 13214, checked in by smasson, 4 years ago

trunk: Mid-year merge, merge back dev_r12563_ASINTER-06_ABL_improvement

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