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_cdf.f90 in utils/tools/SIREN/src – NEMO

source: utils/tools/SIREN/src/iom_cdf.f90 @ 12080

Last change on this file since 12080 was 12080, checked in by jpaul, 4 years ago

update nemo trunk

File size: 99.1 KB
Line 
1!----------------------------------------------------------------------
2! NEMO system team, System and Interface for oceanic RElocable Nesting
3!----------------------------------------------------------------------
4!
5! DESCRIPTION:
6!> @brief NETCDF Input/Output manager :  Library to read Netcdf input files
7!>
8!> @details
9!>    to open netcdf file:<br/>
10!> @code
11!>    CALL iom_cdf_open(td_file)
12!> @endcode
13!>       - td_file is file structure (see @ref file)
14!>
15!>    to write in netcdf file:<br/>
16!> @code
17!>    CALL  iom_cdf_write_file(td_file)
18!> @endcode
19!>
20!>    to close netcdf file:<br/>
21!> @code
22!>    CALL iom_cdf_close(tl_file)
23!> @endcode
24!>
25!>    to read one dimension in netcdf file:<br/>
26!> @code
27!>    tl_dim = iom_cdf_read_dim(tl_file, id_dimid)
28!> @endcode
29!>    or
30!> @code
31!>    tl_dim = iom_cdf_read_dim(tl_file, cd_name)
32!> @endcode
33!>       - id_dimid is dimension id<br/>
34!>       - cd_name is dimension name
35!>
36!>    to read one attribute in netcdf file:<br/>
37!> @code
38!>    tl_att = iom_cdf_read_att(tl_file, id_varid, id_attid)
39!> @endcode
40!>    or
41!> @code
42!>    tl_att = iom_cdf_read_att(tl_file, id_varid, cd_name)
43!> @endcode
44!>       - id_varid is variable id
45!>       - id_attid is attribute id<br/>
46!>       - cd_name is attribute name
47!>   
48!>    to read one variable in netcdf file:<br/>
49!> @code
50!>    tl_var = iom_cdf_read_var(td_file, id_varid, [id_start, id_count])
51!> @endcode
52!>    or
53!> @code
54!>    tl_var = iom_cdf_read_var(td_file, cd_name, [id_start, [id_count,]])
55!> @endcode
56!>       - id_varid is variabale id
57!>       - cd_name is variabale name
58!>       - id_start is a integer(4) 1D array of index from which the data
59!>          values will be read [optional]
60!>       - id_count is a integer(4) 1D array of the number of indices selected
61!>          along each dimension [optional]
62!>
63!> @author
64!> J.Paul
65!>
66!> @date November, 2013 - Initial Version
67!> @date August, 2017
68!> - permit to write header and variable independantly
69!>
70!> @note Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
71!----------------------------------------------------------------------
72MODULE iom_cdf
73
74   USE netcdf                          ! nf90 library
75   USE global                          ! global parameter
76   USE kind                            ! F90 kind parameter
77   USE fct                             ! basic useful function
78   USE logger                          ! log file manager
79   USE att                             ! attribute manage
80   USE dim                             ! dimension manager
81   USE var                             ! variable manager
82   USE file                            ! file manager
83
84   IMPLICIT NONE
85   ! NOTE_avoid_public_variables_if_possible
86
87   ! function and subroutine
88   PUBLIC :: iom_cdf_open           !< open or create netcdf file, return file structure
89   PUBLIC :: iom_cdf_close          !< close netcdf file
90   PUBLIC :: iom_cdf_read_dim       !< read one dimension in an opened netcdf file, return dimension structure
91   PUBLIC :: iom_cdf_read_att       !< read one attribute in an opened netcdf file, return attribute structure
92   PUBLIC :: iom_cdf_read_var       !< read one variable  in an opened netcdf file, return variable  structure
93   PUBLIC :: iom_cdf_fill_var       !< fill variable value in an opened netcdf file
94   PUBLIC :: iom_cdf_write_header   !< write header in an opened netcdf file
95   PUBLIC :: iom_cdf_write_var      !< write variables in an opened netcdf file
96
97   PRIVATE :: iom_cdf__check           ! provides a simple interface to netcdf error message
98   PRIVATE :: iom_cdf__get_info        ! get global information in an opened netcdf file
99   PRIVATE :: iom_cdf__get_file_dim    ! read dimension on an opened netcdf file, and reorder it
100   PRIVATE :: iom_cdf__get_file_att    ! read global attribute on an opened netcdf file
101   PRIVATE :: iom_cdf__get_file_var    ! read information about variable on an opened netcdf file
102   PRIVATE :: iom_cdf__read_dim_id     ! read one dimension in an opened netcdf file, given dimension id.
103   PRIVATE :: iom_cdf__read_dim_name   ! read one dimension in an opened netcdf file, given dimension name.
104   PRIVATE :: iom_cdf__read_att_name   ! read variable or global attribute in an opened netcdf file, given attribute name.
105   PRIVATE :: iom_cdf__read_att_id     ! read variable or global attribute in an opened netcdf file, given attribute id.
106   PRIVATE :: iom_cdf__read_var_id     ! read variable value in an opened netcdf file, given variable id.
107   PRIVATE :: iom_cdf__read_var_name   ! read variable value in an opened netcdf file, given variable name or standard name.
108   PRIVATE :: iom_cdf__read_var_meta   ! read metadata of a variable in an opened netcdf file.
109   PRIVATE :: iom_cdf__read_var_dim    ! read variable dimension in an opened netcdf file.
110   PRIVATE :: iom_cdf__read_var_att    ! read variable attributes in an opened netcdf file.
111   PRIVATE :: iom_cdf__read_var_value  ! read variable value in an opened netcdf file.
112   PRIVATE :: iom_cdf__write_dim_def   ! write dimension definition in an opened netcdf file.
113   PRIVATE :: iom_cdf__write_att_def   ! write attribute definition in an opened netcdf file.
114   PRIVATE :: iom_cdf__write_var_def   ! write variable definition in an opened netcdf file.
115   PRIVATE :: iom_cdf__write_var       ! write a variable in an opened netcdf file.
116   PRIVATE :: iom_cdf__write_var_value ! put variable value in an opened netcdf file.
117   PRIVATE :: iom_cdf__fill_var_id     ! fill variable value in an opened netcdf file, given variable id
118   PRIVATE :: iom_cdf__fill_var_name   ! fill variable value in an opened netcdf file, given variable name
119   PRIVATE :: iom_cdf__fill_var_all    ! fill all variable value in an opened netcdf file
120   PRIVATE :: iom_cdf__del_coord_var   ! remove coordinate variable from an opened netcdf file
121
122   INTERFACE iom_cdf_read_var
123      MODULE PROCEDURE iom_cdf__read_var_id
124      MODULE PROCEDURE iom_cdf__read_var_name
125   END INTERFACE iom_cdf_read_var
126
127   INTERFACE iom_cdf_fill_var
128      MODULE PROCEDURE iom_cdf__fill_var_id
129      MODULE PROCEDURE iom_cdf__fill_var_name
130      MODULE PROCEDURE iom_cdf__fill_var_all
131   END INTERFACE iom_cdf_fill_var
132
133   INTERFACE iom_cdf_read_dim
134      MODULE PROCEDURE iom_cdf__read_dim_id
135      MODULE PROCEDURE iom_cdf__read_dim_name
136   END INTERFACE iom_cdf_read_dim
137
138   INTERFACE iom_cdf_read_att
139      MODULE PROCEDURE iom_cdf__read_att_id
140      MODULE PROCEDURE iom_cdf__read_att_name
141   END INTERFACE iom_cdf_read_att
142
143CONTAINS
144   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
145   SUBROUTINE iom_cdf__check(id_status, cd_msg)
146   !-------------------------------------------------------------------
147   !> @brief This subroutine provides a simple interface to
148   !> netcdf error message
149   !>
150   !> @author J.Paul
151   !> @date November, 2013 - Initial Version
152   !> @date May, 2015
153   !> - add optional message to netcdf error message
154   !>
155   !> @param[in] id_status error status
156   !> @param[in] cd_msg    message
157   !-------------------------------------------------------------------
158
159      IMPLICIT NONE
160
161      ! Argument     
162      INTEGER(i4)     , INTENT(IN)           :: id_status
163      CHARACTER(LEN=*), INTENT(IN), OPTIONAL :: cd_msg
164      ! local variable
165      CHARACTER(LEN=lc) :: cl_msg
166      !----------------------------------------------------------------
167
168      cl_msg=""
169      IF( PRESENT(cd_msg) ) cl_msg=cd_msg
170
171      IF( id_status /= NF90_NOERR )THEN
172         CALL logger_error(TRIM(cl_msg)//TRIM(NF90_STRERROR(id_status)))
173      ENDIF
174
175   END SUBROUTINE iom_cdf__check
176   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
177   SUBROUTINE iom_cdf_open(td_file)
178   !-------------------------------------------------------------------
179   !> @brief This subroutine open a netcdf file in read or write mode.
180   !> @details
181   !> if try to open a file in write mode that did not exist, create it.<br/>
182   !> if file already exist, get information about0:<br/>
183   !> - the number of variables
184   !> - the number of dimensions
185   !> - the number of global attributes
186   !> - the ID of the unlimited dimension
187   !> - the file format
188   !> Finally it read dimensions, and 'longitude' variable to compute East-West
189   !> overlap.
190   !>
191   !> @author J.Paul
192   !> @date November, 2013 - Initial Version
193   !> @date November, 2018
194   !> - write netcdf file as netcdf4
195   !>
196   !> @param[inout] td_file   file structure
197   !-------------------------------------------------------------------
198
199      IMPLICIT NONE
200
201      ! Argument     
202      TYPE(TFILE), INTENT(INOUT)  :: td_file
203
204      ! local variable
205      LOGICAL     :: ll_exist
206      LOGICAL     :: ll_open
207
208      INTEGER(i4) :: il_status
209      INTEGER(i4) :: il_oldmode
210      !----------------------------------------------------------------
211
212      ! check file existence
213      INQUIRE(FILE=TRIM(td_file%c_name), EXIST=ll_exist, OPENED=ll_open)
214      ! ll_open do not work for netcdf file, always return FALSE
215      IF( .NOT. ll_exist .OR. TRIM(td_file%c_type) /= 'cdf' )THEN
216
217         IF( .NOT. td_file%l_wrt )THEN
218
219            CALL logger_fatal( " IOM CDF OPEN: can not open file "//&
220            &               TRIM(td_file%c_name) )
221 
222         ELSE
223
224            CALL logger_info( " IOM CDF CREATE: file "//TRIM(td_file%c_name) )
225
226            il_status = NF90_CREATE(TRIM(td_file%c_name),&
227               &                    cmode=NF90_NETCDF4,  &
228               &                    ncid=td_file%i_id)
229            CALL iom_cdf__check(il_status," IOM CDF CREATE: ")
230            il_status = NF90_SET_FILL(td_file%i_id,   &
231               &                      NF90_NOFILL,    &
232               &                      il_oldmode)
233            CALL iom_cdf__check(il_status," IOM CDF SET FILL: ")
234
235            td_file%l_def=.TRUE.
236            CALL logger_debug( " IOM CDF CREATE: td_file%l_def"//fct_str(td_file%l_def))
237
238         ENDIF
239
240      ELSE
241
242         IF( td_file%i_id /= 0 )THEN
243
244            CALL logger_error( " IOM CDF OPEN: file "//&
245            &               TRIM(td_file%c_name)//" already opened")
246
247         ELSE
248 
249            IF( .NOT. td_file%l_wrt )THEN
250
251               CALL logger_info( " IOM CDF OPEN: file "//&
252               &              TRIM(td_file%c_name)//" in read only mode" )
253
254               il_status = NF90_OPEN( TRIM(td_file%c_name), &
255               &                      NF90_NOWRITE,         &
256               &                      td_file%i_id)
257               CALL iom_cdf__check(il_status," IOM CDF OPEN: ")
258
259            ELSE
260
261               CALL logger_info( "IOM CDF OPEN: file "//&
262               &              TRIM(td_file%c_name)//" in write mode" )
263
264               il_status = NF90_OPEN( TRIM(td_file%c_name), &
265               &                      NF90_WRITE,           &
266               &                      td_file%i_id)
267               CALL iom_cdf__check(il_status,"IOM CDF OPEN: ")
268
269            ENDIF
270
271            ! get general information about file
272            CALL iom_cdf__get_info(td_file)
273
274            ! read dimension in file
275            CALL iom_cdf__get_file_dim(td_file) 
276
277            ! read global attribute in file
278            CALL iom_cdf__get_file_att(td_file)
279
280            ! get information about variables in file
281            CALL iom_cdf__get_file_var(td_file)
282
283            ! remove dimension variable from list of variable
284            CALL iom_cdf__del_coord_var(td_file)
285
286         ENDIF
287
288      ENDIF
289
290   END SUBROUTINE iom_cdf_open
291   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
292   SUBROUTINE iom_cdf_close(td_file)
293   !-------------------------------------------------------------------
294   !> @brief This subroutine close netcdf file.
295   !>
296   !> @author J.Paul
297   !> @date November, 2013 - Initial Version
298   !>
299   !> @param[inout] td_file   file structure
300   !-------------------------------------------------------------------
301
302      IMPLICIT NONE
303
304      ! Argument     
305      TYPE(TFILE), INTENT(INOUT) :: td_file
306
307      ! local variable
308      INTEGER(i4) :: il_status
309      !----------------------------------------------------------------
310
311      ! check if file opened
312      IF( td_file%i_id == 0 )THEN
313
314         CALL logger_error( &
315         &  " IOM CDF CLOSE: no id associated to file "//TRIM(td_file%c_name))
316
317      ELSE
318         CALL logger_info( &
319         &  " IOM CDF CLOSE: file "//TRIM(td_file%c_name))
320
321         il_status = NF90_CLOSE(td_file%i_id)
322         CALL iom_cdf__check(il_status,"IOM CDF CLOSE: ")
323
324         td_file%i_id = 0
325
326      ENDIF
327
328   END SUBROUTINE iom_cdf_close
329   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
330   SUBROUTINE iom_cdf__get_info(td_file)
331   !-------------------------------------------------------------------
332   !> @brief This subroutine get global information in an opened netcdf
333   !> file.
334   !> @details
335   !> It gets the number of variables, the number of dimensions,
336   !> the number of global attributes, the ID of the unlimited dimension
337   !> and finally the format version and filled file strucuture with it.
338   !>
339   !> @author J.Paul
340   !> @date November, 2013 - Initial Version
341   !> @date October, 2016
342   !> - define cdf4 as cdf.
343   !>
344   !> @param[inout] td_file   file structure
345   !-------------------------------------------------------------------
346
347      IMPLICIT NONE
348
349      ! Argument     
350      TYPE(TFILE), INTENT(INOUT) :: td_file
351
352      ! local variable
353      INTEGER(i4) :: il_fmt   ! format version
354      INTEGER(i4) :: il_status
355      !----------------------------------------------------------------
356
357      CALL logger_trace( &
358      &  " IOM CDF GET INFO: about netcdf file "//TRIM(td_file%c_name))
359
360      il_status=NF90_INQUIRE(td_file%i_id, td_file%i_ndim, &
361      &     td_file%i_nvar, td_file%i_natt, td_file%i_uldid, il_fmt)
362      CALL iom_cdf__check(il_status,"IOM CDF GET INFO: ")
363
364      SELECT CASE(il_fmt)
365         CASE(nf90_format_classic, nf90_format_64bit)
366            td_file%c_type='cdf'
367         CASE(nf90_format_netcdf4,nf90_format_netcdf4_classic)
368            td_file%c_type='cdf'
369      END SELECT
370      CALL logger_debug("IOM CDF GET INFO: type "//TRIM(td_file%c_type)) 
371
372      ! record header infos
373      td_file%i_rhd=1
374
375   END SUBROUTINE iom_cdf__get_info
376   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
377   SUBROUTINE iom_cdf__get_file_dim(td_file)
378   !-------------------------------------------------------------------
379   !> @brief This subroutine read dimension on an opened netcdf file, and
380   !> reorder dimension to ('x', 'y', 'z', 't').
381   !> The dimension structure inside file structure is then completed.
382   !>
383   !> @author J.Paul
384   !> @date November, 2013 - Initial Version
385   !> @date October, 2016
386   !> - check unknown dimension
387   !> @date January, 2019
388   !> - clean dimension structure
389   !>
390   !> @param[inout] td_file   file structure
391   !-------------------------------------------------------------------
392
393      IMPLICIT NONE
394
395      ! Argument     
396      TYPE(TFILE), INTENT(INOUT) :: td_file
397
398      ! local variable
399      TYPE(TDIM) :: tl_dim
400
401      ! loop indices
402      INTEGER(i4) :: ji
403      INTEGER(i4) :: ii
404      !----------------------------------------------------------------
405
406      ! clean dimension
407      DO ji=1,ip_maxdim
408         CALL dim_clean(td_file%t_dim(ji))
409      ENDDO
410
411      IF( td_file%i_ndim > 0 )THEN
412
413         ii=1
414         DO ji = 1, td_file%i_ndim
415            ! read dimension information
416            tl_dim=iom_cdf_read_dim( td_file, ji)
417            ! sname == 'u' if dimension is unknown (not to be used)
418            IF( TRIM(tl_dim%c_sname) /= 'u' )THEN
419               IF( ii > ip_maxdim )THEN
420                  CALL logger_fatal("IOM CDF OPEN: too much dimension "//&
421                  & "to be read. you could choose dimension to be used. see "//&
422                  & " configuration file")
423               ENDIF
424               td_file%t_dim(ii)=dim_copy(tl_dim)
425               ii=ii+1
426            ENDIF
427            ! clean
428            CALL dim_clean(tl_dim)
429         ENDDO
430
431         ! inform unlimited dimension
432         IF( td_file%i_uldid == -1 )THEN
433            CALL logger_warn( &
434            &  " IOM CDF GET FILE DIM: there is no unlimited dimension in file "//&
435            &  TRIM(td_file%c_name))
436         !ELSE
437         !   td_file%t_dim( td_file%i_uldid )%l_uld=.TRUE.
438         ENDIF
439
440      ELSE
441
442         CALL logger_warn( &
443         &  " IOM CDF GET FILE DIM: there is no dimension in file "//&
444         &  TRIM(td_file%c_name))
445
446      ENDIF
447
448      ! reorder dimension to ('x','y','z','t')
449      CALL dim_reorder(td_file%t_dim(:))
450
451   END SUBROUTINE iom_cdf__get_file_dim
452   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
453   SUBROUTINE iom_cdf__get_file_att(td_file)
454   !-------------------------------------------------------------------
455   !> @brief This subroutine read global attribute on an opened netcdf
456   !> file.
457   !> The attribute structure inside file structure is then completed.
458   !>
459   !> @author J.Paul
460   !> @date November, 2013 - Initial Version
461   !> @date September, 2014
462   !> - use attribute periodicity read from the file if present.
463   !> @date January, 2019
464   !> - clean attribute structure
465   !>
466   !> @param[inout] td_file   file structure
467   !-------------------------------------------------------------------
468
469      IMPLICIT NONE
470
471      ! Argument     
472      TYPE(TFILE), INTENT(INOUT) :: td_file
473
474      ! local variable
475      TYPE(TATT) :: tl_att
476
477      ! loop indices
478      INTEGER(i4) :: ji
479      INTEGER(i4) :: ii
480      !----------------------------------------------------------------
481      CALL logger_trace("IOM CDF GET FILE ATT : get attr(s) in &
482         &     file "//TRIM(td_file%c_name))
483
484      IF( td_file%i_natt > 0 )THEN
485         IF(ASSOCIATED(td_file%t_att))THEN
486            CALL att_clean(td_file%t_att(:))
487            DEALLOCATE(td_file%t_att)
488         ENDIF
489         ALLOCATE(td_file%t_att(td_file%i_natt))
490
491         ii=1
492         DO ji = 1, td_file%i_natt
493            ! read global attribute
494            tl_att=iom_cdf_read_att( td_file, NF90_GLOBAL, ji)
495            IF( .NOT. att_is_dummy(tl_att) )THEN
496               td_file%t_att(ii)=att_copy(tl_att)
497               ii=ii+1
498            ENDIF
499            ! clean
500            CALL att_clean(tl_att)
501         ENDDO
502
503      ELSE
504         CALL logger_debug( &
505         &  " IOM CDF GET FILE ATT: there is no global attribute in file "//&
506         &  TRIM(td_file%c_name))
507      ENDIF
508
509   END SUBROUTINE iom_cdf__get_file_att
510   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
511   SUBROUTINE iom_cdf__get_file_var(td_file)
512   !-------------------------------------------------------------------
513   !> @brief This subroutine read information about variable of an
514   !> opened netcdf file.
515   !> The variable structure inside file structure is then completed.
516   !> @note variable value are not read !
517   !>
518   !> @author J.Paul
519   !> @date November, 2013 - Initial Version
520   !> @date September, 2015
521   !> - manage useless (dummy) variable
522   !> @date January, 2016
523   !> - increment n3d for 4D variable
524   !> @date October, 2016
525   !> - check variable to be used (variable's dimension allowed and variable
526   !> not "dummy")
527   !>
528   !> @param[inout] td_file   file structure
529   !-------------------------------------------------------------------
530
531      IMPLICIT NONE
532
533      ! Argument     
534      TYPE(TFILE), INTENT(INOUT) :: td_file
535
536      ! local variable
537      INTEGER(i4) :: il_attid
538      INTEGER(i4) :: il_nvar
539
540      TYPE(TVAR), DIMENSION(:), ALLOCATABLE  :: tl_var
541
542      ! loop indices
543      INTEGER(i4) :: ji
544      INTEGER(i4) :: ii
545      !----------------------------------------------------------------
546
547      IF( td_file%i_nvar > 0 )THEN
548
549         IF(ASSOCIATED(td_file%t_var))THEN
550            CALL var_clean(td_file%t_var(:))
551            DEALLOCATE(td_file%t_var)
552         ENDIF
553
554         il_nvar=td_file%i_nvar
555         ALLOCATE(tl_var(il_nvar))
556         DO ji = 1, il_nvar
557           ! read variable information
558           tl_var(ji)=iom_cdf__read_var_meta( td_file, ji) 
559         ENDDO
560
561         ! update number of variable used
562         td_file%i_nvar=COUNT(tl_var(:)%l_use)
563
564         ALLOCATE(td_file%t_var(td_file%i_nvar))
565
566         ii=0
567         DO ji = 1, il_nvar
568            IF( tl_var(ji)%l_use )THEN
569               ii=ii+1
570               td_file%t_var(ii)=var_copy(tl_var(ji))
571               SELECT CASE(td_file%t_var(ii)%i_ndim)
572                  CASE(0)
573                     td_file%i_n0d=td_file%i_n0d+1
574                  CASE(1)
575                     td_file%i_n1d=td_file%i_n1d+1
576                     td_file%i_rhd=td_file%i_rhd+1
577                  CASE(2)
578                     td_file%i_n2d=td_file%i_n2d+1
579                     td_file%i_rhd=td_file%i_rhd+1
580                  CASE(3,4)
581                     td_file%i_n3d=td_file%i_n3d+1
582                     td_file%i_rhd=td_file%i_rhd+td_file%t_dim(3)%i_len
583               END SELECT
584
585               ! look for depth id
586               IF( INDEX(TRIM(fct_lower(td_file%t_var(ii)%c_name)),'depth')/=0 )THEN
587                  IF( td_file%i_depthid == 0 )THEN
588                     td_file%i_depthid=ji
589                  ELSE
590                     IF( td_file%i_depthid /= ji )THEN
591                        CALL logger_error("IOM CDF GET FILE VAR: find more"//&
592                           &  " than one depth variable in file "//&
593                           &  TRIM(td_file%c_name) )
594                     ENDIF
595                  ENDIF
596               ENDIF
597
598               ! look for time id
599               IF( INDEX(TRIM(fct_lower(td_file%t_var(ii)%c_name)),'time')/=0 )THEN
600                  IF( td_file%i_timeid == 0 )THEN
601                     td_file%i_timeid=ji
602                  ELSE
603                     IF( td_file%i_timeid /= ji )THEN
604                        CALL logger_warn("IOM CDF GET FILE VAR: find more "//&
605                        &                 "than one time variable in file "//&
606                        &                 TRIM(td_file%c_name)//". see "//&
607                        &                 "dummy.cfg configuration file to"//&
608                        &                 " not used dummy variables.")
609                     ENDIF
610                     il_attid=0
611                     IF( ASSOCIATED(td_file%t_var(ii)%t_att) )THEN
612                        il_attid=att_get_id(td_file%t_var(ii)%t_att(:),'calendar')
613                     ENDIF
614                     IF( il_attid /= 0 )THEN
615                        td_file%i_timeid=ji
616                     ENDIF
617                  ENDIF
618               ENDIF
619
620            ENDIF
621         ENDDO
622
623         CALL var_clean(tl_var(:))
624         DEALLOCATE(tl_var)
625
626      ELSE
627         CALL logger_debug( &
628         &  " IOM CDF GET FILE VAR: there is no variable in file "//&
629         &  TRIM(td_file%c_name))
630      ENDIF
631
632   END SUBROUTINE iom_cdf__get_file_var
633   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
634   SUBROUTINE iom_cdf__del_coord_var(td_file)
635   !-------------------------------------------------------------------
636   !> @brief This subroutine delete coordinate variable from an
637   !> opened netcdf file if present.
638   !>
639   !> @author J.Paul
640   !> @date November, 2013 - Initial Version
641   !>
642   !> @param[inout] td_file   file structure
643   !-------------------------------------------------------------------
644
645      IMPLICIT NONE
646
647      ! Argument     
648      TYPE(TFILE), INTENT(INOUT) :: td_file
649
650      ! local variable
651      CHARACTER(LEN=lc) :: cl_name
652      CHARACTER(LEN=lc) :: cl_sname
653
654      ! loop indices
655      INTEGER(i4) :: ji
656      INTEGER(i4) :: jj
657      !----------------------------------------------------------------
658
659      IF( td_file%i_nvar > 0 )THEN
660         DO ji=td_file%i_nvar,1,-1
661            cl_name=TRIM(td_file%t_var(ji)%c_name)
662            DO jj=1,ip_maxdim
663               IF( td_file%t_dim(jj)%l_use )THEN
664                  cl_sname=fct_upper(td_file%t_dim(jj)%c_sname)
665                  IF( TRIM(cl_name) == TRIM(cl_sname) )THEN
666                     CALL file_del_var(td_file,TRIM(cl_name))
667                     EXIT
668                  ENDIF
669               ENDIF
670            ENDDO
671         ENDDO
672      ELSE
673         CALL logger_debug( &
674         &  " IOM CDF DEL VAR DIM: there is no variable in file "//&
675         &  TRIM(td_file%c_name))
676      ENDIF
677
678   END SUBROUTINE iom_cdf__del_coord_var
679   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
680   FUNCTION iom_cdf__read_dim_id(td_file, id_dimid) &
681         & RESULT (tf_dim)
682   !-------------------------------------------------------------------
683   !> @brief This function read one dimension in an opened netcdf file,
684   !> given dimension id.
685   !>
686   !> @author J.Paul
687   !> @date November, 2013 - Initial Version
688   !> @date February, 2015 - create unused dimension, when reading dimension
689   !> of length less or equal to zero
690   !>
691   !> @param[in] td_file   file structure
692   !> @param[in] id_dimid  dimension id
693   !> @return  dimension structure
694   !-------------------------------------------------------------------
695
696      IMPLICIT NONE
697
698      ! Argument     
699      TYPE(TFILE), INTENT(IN) :: td_file
700      INTEGER(i4), INTENT(IN) :: id_dimid
701
702      ! function
703      TYPE(TDIM)              :: tf_dim
704
705      ! local variable
706      INTEGER(i4)       :: il_status
707      INTEGER(i4)       :: il_len
708      CHARACTER(LEN=lc) :: cl_name
709      LOGICAL           :: ll_use
710      !----------------------------------------------------------------
711
712      ! check if file opened
713      IF( td_file%i_id == 0 )THEN
714
715         CALL logger_error( &
716         &  " IOM CDF READ DIM: no id associated to file "//TRIM(td_file%c_name))
717
718      ELSE
719
720         CALL logger_trace( &
721         &  " IOM CDF READ DIM: dimension "//TRIM(fct_str(id_dimid))//&
722         &  " in file "//TRIM(td_file%c_name))
723
724         il_status=NF90_INQUIRE_DIMENSION(td_file%i_id, id_dimid, &
725         &                                cl_name, il_len )
726         CALL iom_cdf__check(il_status,"IOM CDF READ DIM: ")
727
728         ll_use=.TRUE.
729         IF( il_len <= 0 )THEN
730            CALL logger_warn( &
731         &  " IOM CDF READ DIM: dimension "//TRIM(fct_str(id_dimid))//&
732         &  " in file "//TRIM(td_file%c_name)//" is less or equel to zero")
733            il_len=1
734            ll_use=.FALSE.
735         ENDIF
736         tf_dim=dim_init(cl_name, il_len, ld_use=ll_use)
737
738      ENDIF
739
740      tf_dim%i_id=id_dimid
741
742   END FUNCTION iom_cdf__read_dim_id
743   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
744   FUNCTION iom_cdf__read_dim_name(td_file, cd_name) &
745         & RESULT (tf_dim)
746   !-------------------------------------------------------------------
747   !> @brief This function read one dimension in an opened netcdf file,
748   !> given dimension name.
749   !>
750   !> @author J.Paul
751   !> @date November, 2013 - Initial Version
752   !>
753   !> @param[in] td_file   file structure
754   !> @param[in] cd_name   dimension name
755   !> @return  dimension structure
756   !-------------------------------------------------------------------
757
758      IMPLICIT NONE
759
760      ! Argument     
761      TYPE(TFILE),      INTENT(IN) :: td_file
762      CHARACTER(LEN=*), INTENT(IN) :: cd_name
763
764      ! function
765      TYPE(TDIM)                   :: tf_dim 
766
767      ! local variable
768      INTEGER(i4) :: il_status
769      INTEGER(i4) :: il_dimid
770      !----------------------------------------------------------------
771
772      ! check if file opened
773      IF( td_file%i_id == 0 )THEN
774
775         CALL logger_error( &
776         &  " IOM CDF READ DIM: no id associated to file "//&
777         &  TRIM(td_file%c_name))
778
779      ELSE     
780
781         il_status=NF90_INQ_DIMID( td_file%i_id, TRIM(ADJUSTL(cd_name)), &
782         &                         il_dimid)
783         CALL iom_cdf__check(il_status,"IOM CDF READ DIM: ")
784
785         tf_dim=iom_cdf_read_dim(td_file, il_dimid)
786
787      ENDIF
788
789   END FUNCTION iom_cdf__read_dim_name
790   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
791   FUNCTION iom_cdf__read_att_name(td_file, id_varid, cd_name) &
792         & RESULT (tf_att)
793   !-------------------------------------------------------------------
794   !> @brief This function read variable or global attribute in an opened
795   !> netcdf file, given attribute name.
796   !>
797   !> @author J.Paul
798   !> @date November, 2013 - Initial Version
799   !> @date November 2017
800   !> - check if cl_value is not bug
801   !>
802   !> @param[in] td_file   file structure
803   !> @param[in] id_varid  variable id. use NF90_GLOBAL to read global
804   !> attribute in a file
805   !> @param[in] cd_name   attribute name
806   !> @return  attribute structure
807   !-------------------------------------------------------------------
808
809      IMPLICIT NONE
810
811      ! Argument     
812      TYPE(TFILE),      INTENT(IN) :: td_file
813      INTEGER(i4),      INTENT(IN) :: id_varid
814      CHARACTER(LEN=*), INTENT(IN) :: cd_name
815   
816      ! function
817      TYPE(TATT)                   :: tf_att
818
819      ! local variable
820      CHARACTER(LEN=lc) :: cl_name
821
822      INTEGER(i4) :: il_status
823      INTEGER(i4) :: il_attid
824      INTEGER(i4) :: il_type
825      INTEGER(i4) :: il_len
826
827      CHARACTER(LEN=lc) :: cl_value
828     
829      INTEGER(i1), DIMENSION(:), ALLOCATABLE :: bl_value
830      INTEGER(i2), DIMENSION(:), ALLOCATABLE :: sl_value
831      INTEGER(i4), DIMENSION(:), ALLOCATABLE :: il_value
832      REAL(sp)   , DIMENSION(:), ALLOCATABLE :: rl_value
833      REAL(dp)   , DIMENSION(:), ALLOCATABLE :: dl_value
834      !----------------------------------------------------------------
835      ! check if file opened
836      IF( td_file%i_id == 0 )THEN
837
838         CALL logger_error( &
839            &  " IOM CDF READ ATT: no id associated to file "//TRIM(td_file%c_name))
840
841      ELSE     
842
843         cl_name=TRIM(ADJUSTL(cd_name))
844
845         ! inquire attribute
846         IF( id_varid == NF90_GLOBAL )THEN
847
848            CALL logger_trace( &
849               &  " IOM CDF READ ATT: inquire global attribute "//&
850               &  " in file "//TRIM(td_file%c_name))
851
852         ELSE
853
854            CALL logger_trace( &
855               &  " IOM CDF READ ATT: inquire attribute "//&
856               &  " of variable "//TRIM(fct_str(id_varid))//&
857               &  " in file "//TRIM(td_file%c_name))
858
859         ENDIF
860
861         il_status=NF90_INQUIRE_ATTRIBUTE(td_file%i_id, id_varid,  &
862            &                             cl_name,&
863            &                             il_type,&
864            &                             il_len, &
865            &                             il_attid )
866         CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")
867
868         !! get attribute value
869         CALL logger_debug( " IOM CDF READ ATT: get attribute "//&
870            &            TRIM(cl_name)//" in file "//TRIM(td_file%c_name))
871
872         SELECT CASE( il_type )
873
874            CASE(NF90_CHAR)
875               CALL logger_debug( " IOM CDF READ ATT: get NF90_CHAR ")
876
877               ! check string lengths
878               IF( LEN(cl_value) < il_len )THEN
879
880                  CALL logger_warn( &
881                     &  " IOM CDF READ ATT: not enough space to put "//&
882                     &  "attribute "//TRIM(cl_name) )
883
884               ELSE
885
886                  ! Read the attributes
887                  il_status=NF90_GET_ATT(td_file%i_id, id_varid, &
888                     &                   cl_name, &
889                     &                   cl_value )
890                  CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")
891
892                  ! check cl_value
893                  CALL logger_trace("IOM CDF READ ATT: value = "//TRIM(cl_value))
894                  IF( LLT(cl_value,'') ) cl_value = ''
895                  tf_att=att_init(cl_name, cl_value)
896
897               ENDIF
898         
899            CASE(NF90_BYTE)
900               CALL logger_debug( " IOM CDF READ ATT: get NF90_BYTE ")
901
902               ALLOCATE( bl_value(il_len), &
903                  &      stat=il_status)
904               IF(il_status /= 0 )THEN
905
906                  CALL logger_error( "IOM CDF READ ATT: "//&
907                     &  "not enough space to put attribute "//TRIM(cl_name) )
908
909               ELSE
910
911                  ! Read the attributes
912                  il_status=NF90_GET_ATT(td_file%i_id, id_varid, &
913                     &                   cl_name, &
914                     &                   bl_value(:))
915                  CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")   
916
917                  tf_att=att_init(cl_name, bl_value(:))
918
919               ENDIF
920
921               DEALLOCATE(bl_value)
922
923            CASE(NF90_SHORT)
924               CALL logger_debug( " IOM CDF READ ATT: get NF90_SHORT ")
925
926               ALLOCATE( sl_value(il_len), &
927                  &      stat=il_status)
928               IF(il_status /= 0 )THEN
929
930                  CALL logger_error( &
931                     &  " IOM CDF READ ATT: not enough space to put "//&
932                     &  "attribute "//TRIM(cl_name) )
933
934               ELSE
935
936                  ! Read the attributes
937                  il_status=NF90_GET_ATT(td_file%i_id, id_varid, &
938                     &                   cl_name, &
939                     &                   sl_value(:))
940                  CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")   
941
942                  tf_att=att_init(cl_name, sl_value(:))
943
944               ENDIF
945
946               DEALLOCATE(sl_value)
947
948            CASE(NF90_INT)
949               CALL logger_debug( " IOM CDF READ ATT: get NF90_INT ")
950
951               ALLOCATE( il_value(il_len), &
952                  &      stat=il_status)
953               IF(il_status /= 0 )THEN
954
955                  CALL logger_error( &
956                     &  " IOM CDF READ ATT: not enough space to put "//&
957                     &  "attribute "//TRIM(cl_name) )
958
959               ELSE
960
961                  ! Read the attributes
962                  il_status=NF90_GET_ATT(td_file%i_id, id_varid, &
963                     &                   cl_name, &
964                     &                   il_value(:))
965                  CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")   
966
967                  tf_att=att_init(cl_name, il_value(:))
968               ENDIF
969
970               DEALLOCATE(il_value)
971
972            CASE(NF90_FLOAT)
973               CALL logger_debug( " IOM CDF READ ATT: get NF90_FLOAT ")
974
975               ALLOCATE( rl_value(il_len), &
976                  &      stat=il_status)
977               IF(il_status /= 0 )THEN
978
979                  CALL logger_error( &
980                  &  " IOM CDF READ ATT: not enough space to put "//&
981                  &  "attribute "//TRIM(cl_name) )
982
983               ELSE
984
985                  ! Read the attributes
986                  il_status=NF90_GET_ATT(td_file%i_id, id_varid, &
987                     &                   cl_name, &
988                     &                   rl_value(:))
989                  CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")   
990
991                  tf_att=att_init(cl_name, rl_value(:))
992
993               ENDIF
994
995               DEALLOCATE(rl_value)
996
997            CASE(NF90_DOUBLE)
998               CALL logger_debug( " IOM CDF READ ATT: get NF90_DOUBLE ")
999
1000               ALLOCATE( dl_value(il_len), &
1001                  &      stat=il_status)
1002               IF(il_status /= 0 )THEN
1003
1004                  CALL logger_error( &
1005                     &  " IOM CDF READ ATT: not enough space to put "//&
1006                     &  "attribute "//TRIM(cl_name) )
1007
1008               ELSE
1009
1010                  ! Read the attributes
1011                  il_status=NF90_GET_ATT(td_file%i_id, id_varid, &
1012                     &                   cl_name, &
1013                     &                   dl_value(:))
1014                  CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")   
1015
1016                  tf_att=att_init(cl_name, dl_value(:))
1017
1018               ENDIF
1019
1020               DEALLOCATE(dl_value)
1021
1022         END SELECT
1023
1024         tf_att%i_id=il_attid
1025
1026      ENDIF
1027
1028   END FUNCTION iom_cdf__read_att_name
1029   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1030   FUNCTION iom_cdf__read_att_id(td_file, id_varid, id_attid) &
1031         & RESULT (tf_att)
1032   !-------------------------------------------------------------------
1033   !> @brief This function read variable or global attribute in an opened
1034   !> netcdf file, given attribute id.
1035   !>
1036   !> @author J.Paul
1037   !> @date November, 2013 - Initial Version
1038   !>
1039   !> @param[in] td_file   file structure
1040   !> @param[in] id_varid  variable id. use NF90_GLOBAL to read global
1041   !> attribute in a file
1042   !> @param[in] id_attid  attribute id
1043   !> @return  attribute structure
1044   !-------------------------------------------------------------------
1045
1046      IMPLICIT NONE
1047
1048      ! Argument     
1049      TYPE(TFILE), INTENT(IN) :: td_file
1050      INTEGER(i4), INTENT(IN) :: id_varid
1051      INTEGER(i4), INTENT(IN) :: id_attid
1052
1053      ! function
1054      TYPE(TATT)              :: tf_att
1055
1056      ! local variable
1057      INTEGER(i4)       :: il_status
1058      CHARACTER(LEN=lc) :: cl_name
1059      !----------------------------------------------------------------
1060      ! check if file opened
1061      IF( td_file%i_id == 0 )THEN
1062
1063         CALL logger_error( &
1064            &  "IOM CDF READ ATT: no id associated to file "//TRIM(td_file%c_name))
1065
1066      ELSE
1067
1068         CALL logger_trace( " IOM CDF READ ATT ID: get attribute "//&
1069            &               TRIM(fct_str(id_attid))//" of var "//&
1070            &               TRIM(fct_str(id_varid))//" in file "//&
1071            &               TRIM(td_file%c_name) &
1072         )
1073
1074         ! get attribute name
1075         il_status=NF90_INQ_ATTNAME(td_file%i_id, id_varid, id_attid, cl_name)
1076         CALL iom_cdf__check(il_status,"IOM CDF READ ATT: ")
1077
1078         ! read attribute
1079         tf_att=iom_cdf__read_att_name(td_file, id_varid, cl_name)
1080
1081      ENDIF
1082
1083   END FUNCTION iom_cdf__read_att_id
1084   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1085   FUNCTION iom_cdf__read_var_id(td_file, id_varid, id_start, id_count) &
1086         & RESULT (tf_var)
1087   !-------------------------------------------------------------------
1088   !> @brief This function read variable value in an opened
1089   !> netcdf file, given variable id.
1090   !> @details
1091   !> Optionaly, start indices and number of indices selected along each dimension
1092   !> could be specify in a 4 dimension array (/'x','y','z','t'/)
1093   !>
1094   !> @author J.Paul
1095   !> @date November, 2013 - Initial Version
1096   !>
1097   !> @param[in] td_file   file structure
1098   !> @param[in] id_varid  variable id
1099   !> @param[in] id_start  index in the variable from which the data values
1100   !> will be read
1101   !> @param[in] id_count  number of indices selected along each dimension
1102   !> @return  variable structure
1103   !-------------------------------------------------------------------
1104
1105      IMPLICIT NONE
1106
1107      ! Argument     
1108      TYPE(TFILE),               INTENT(IN) :: td_file
1109      INTEGER(i4),               INTENT(IN) :: id_varid
1110      INTEGER(i4), DIMENSION(:), INTENT(IN), OPTIONAL :: id_start
1111      INTEGER(i4), DIMENSION(:), INTENT(IN), OPTIONAL :: id_count
1112
1113      ! function
1114      TYPE(TVAR)                            :: tf_var
1115
1116      ! local variable
1117      INTEGER(i4), DIMENSION(1) :: il_ind
1118      !----------------------------------------------------------------
1119      ! check if file opened
1120      IF( td_file%i_id == 0 )THEN
1121
1122         CALL logger_error( &
1123         &  " IOM CDF READ VAR: no id associated to file "//TRIM(td_file%c_name))
1124
1125      ELSE
1126
1127         ! look for variable index
1128         il_ind(:)=MINLOC(td_file%t_var(:)%i_id,mask=(td_file%t_var(:)%i_id==id_varid))
1129         IF( il_ind(1) /= 0 )THEN
1130
1131            tf_var=var_copy(td_file%t_var(il_ind(1)))
1132
1133            !!! read variable value
1134            CALL iom_cdf__read_var_value(td_file, tf_var, id_start, id_count)
1135
1136         ELSE
1137            CALL logger_error( &
1138            &  " IOM CDF READ VAR: there is no variable with id "//&
1139            &  TRIM(fct_str(id_varid))//" in file "//TRIM(td_file%c_name))
1140         ENDIF
1141
1142      ENDIF
1143   END FUNCTION iom_cdf__read_var_id
1144   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1145   FUNCTION iom_cdf__read_var_name(td_file, cd_name, id_start, id_count) &
1146         & RESULT (tf_var)
1147   !-------------------------------------------------------------------
1148   !> @brief This function read variable value in an opened
1149   !> netcdf file, given variable name or standard name.
1150   !> @details
1151   !> Optionaly, start indices and number of indices selected along each dimension
1152   !> could be specify in a 4 dimension array (/'x','y','z','t'/)
1153   !>
1154   !> look first for variable name. If it doesn't
1155   !> exist in file, look for variable standard name.<br/>
1156   !>
1157   !> @author J.Paul
1158   !> @date November, 2013 - Initial Version
1159   !>
1160   !> @param[in] td_file   file structure
1161   !> @param[in] cd_name   variable name or standard name.
1162   !> @param[in] id_start  index in the variable from which the data values will be read
1163   !> @param[in] id_count  number of indices selected along each dimension
1164   !> @return  variable structure
1165   !-------------------------------------------------------------------
1166
1167      IMPLICIT NONE
1168
1169      ! Argument     
1170      TYPE(TFILE)     ,                INTENT(IN) :: td_file
1171      CHARACTER(LEN=*),                INTENT(IN), OPTIONAL :: cd_name
1172      INTEGER(i4)     , DIMENSION(:),  INTENT(IN), OPTIONAL :: id_start
1173      INTEGER(i4)     , DIMENSION(:),  INTENT(IN), OPTIONAL :: id_count
1174
1175      ! function
1176      TYPE(TVAR)                                  :: tf_var
1177
1178      ! local variable
1179      INTEGER(i4)       :: il_varid
1180      !----------------------------------------------------------------
1181      ! check if file opened
1182      IF( td_file%i_id == 0 )THEN
1183
1184         CALL logger_error( &
1185         &  " IOM CDF READ VAR: no id associated to file "//TRIM(td_file%c_name))
1186
1187      ELSE
1188
1189         IF( .NOT. PRESENT(cd_name) )THEN
1190
1191            CALL logger_error( &
1192            &  " IOM CDF READ VAR: you must specify a variable to read "//&
1193            &  " in file "//TRIM(td_file%c_name))
1194
1195         ELSE
1196
1197            il_varid=var_get_index(td_file%t_var(:), cd_name)
1198            IF( il_varid /= 0 )THEN
1199
1200               tf_var=var_copy(td_file%t_var(il_varid))
1201
1202               !!! read variable value
1203               CALL iom_cdf__read_var_value( td_file, tf_var, id_start, id_count)
1204
1205            ELSE
1206
1207               CALL logger_error( &
1208               &  " IOM CDF READ VAR: there is no variable with "//&
1209               &  " name or standard name "//TRIM(cd_name)//&
1210               &  " in file "//TRIM(td_file%c_name) )
1211            ENDIF
1212
1213         ENDIF
1214
1215      ENDIF
1216 
1217   END FUNCTION iom_cdf__read_var_name
1218   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1219   SUBROUTINE iom_cdf__fill_var_all(td_file, id_start, id_count)
1220   !-------------------------------------------------------------------
1221   !> @brief This subroutine fill all variable value from an opened
1222   !> netcdf file.
1223   !> @details
1224   !> Optionaly, start indices and number of indices selected along each dimension
1225   !> could be specify in a 4 dimension array (/'x','y','z','t'/)
1226   !>
1227   !> @author J.Paul
1228   !> @date November, 2013 - Initial Version
1229   !>
1230   !> @param[inout] td_file   file structure
1231   !> @param[in] id_start     index in the variable from which the data values
1232   !> will be read
1233   !> @param[in] id_count     number of indices selected along each dimension
1234   !-------------------------------------------------------------------
1235
1236      IMPLICIT NONE
1237
1238      ! Argument     
1239      TYPE(TFILE),               INTENT(INOUT) :: td_file
1240      INTEGER(i4), DIMENSION(:), INTENT(IN   ),  OPTIONAL :: id_start
1241      INTEGER(i4), DIMENSION(:), INTENT(IN   ),  OPTIONAL :: id_count
1242
1243      ! local variable
1244
1245      ! loop indices
1246      INTEGER(i4) :: ji
1247      !----------------------------------------------------------------
1248      ! check if file opened
1249      IF( td_file%i_id == 0 )THEN
1250
1251         CALL logger_error( &
1252         &  " IOM CDF FILL VAR: no id associated to file "//TRIM(td_file%c_name))
1253
1254      ELSE
1255
1256         DO ji=1,td_file%i_nvar
1257            CALL iom_cdf_fill_var(td_file, td_file%t_var(ji)%i_id, &
1258            &                     id_start, id_count)
1259         ENDDO
1260
1261      ENDIF
1262
1263   END SUBROUTINE iom_cdf__fill_var_all
1264   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1265   SUBROUTINE iom_cdf__fill_var_id(td_file, id_varid, id_start, id_count)
1266   !-------------------------------------------------------------------
1267   !> @brief This subroutine fill variable value in an opened
1268   !> netcdf file, given variable id.
1269   !> @details
1270   !> Optionaly, start indices and number of indices selected along each dimension
1271   !> could be specify in a 4 dimension array (/'x','y','z','t'/)
1272   !
1273   !> @author J.Paul
1274   !> @date November, 2013 - Initial Version
1275   !
1276   !> @param[inout] td_file   file structure
1277   !> @param[in] id_varid     variable id
1278   !> @param[in] id_start     index in the variable from which the data values
1279   !> will be read
1280   !> @param[in] id_count     number of indices selected along each dimension
1281   !-------------------------------------------------------------------
1282
1283      IMPLICIT NONE
1284
1285      ! Argument     
1286      TYPE(TFILE),               INTENT(INOUT) :: td_file
1287      INTEGER(i4),               INTENT(IN)    :: id_varid
1288      INTEGER(i4), DIMENSION(:), INTENT(IN),  OPTIONAL :: id_start
1289      INTEGER(i4), DIMENSION(:), INTENT(IN),  OPTIONAL :: id_count
1290
1291      ! local variable
1292      INTEGER(i4), DIMENSION(1) :: il_varid
1293
1294      ! loop indices
1295      INTEGER(i4) :: ji
1296      !----------------------------------------------------------------
1297      ! check if file opened
1298      IF( td_file%i_id == 0 )THEN
1299
1300         CALL logger_error( &
1301         &  "IOM CDF FILL VAR: no id associated to file "//TRIM(td_file%c_name))
1302
1303      ELSE
1304
1305         ! look for variable id
1306         il_varid(:)=MINLOC( td_file%t_var(:)%i_id, &
1307         &                 mask=(td_file%t_var(:)%i_id==id_varid))
1308         IF( il_varid(1) /= 0 )THEN
1309
1310            !!! read variable value
1311            CALL iom_cdf__read_var_value(td_file, td_file%t_var(il_varid(1)), &
1312            &                            id_start, id_count)
1313
1314            DO ji=1,td_file%i_nvar
1315               CALL logger_debug(" IOM CDF FILL VAR: var id "//&
1316               &     TRIM(td_file%t_var(ji)%c_name)//" "//&
1317               &     TRIM(fct_str(td_file%t_var(ji)%i_id)) )
1318            ENDDO
1319         ELSE
1320            CALL logger_error( &
1321            &  " IOM CDF FILL VAR: there is no variable with id "//&
1322            &  TRIM(fct_str(id_varid))//" in file "//TRIM(td_file%c_name))
1323         ENDIF
1324
1325      ENDIF
1326
1327   END SUBROUTINE iom_cdf__fill_var_id
1328   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1329   SUBROUTINE iom_cdf__fill_var_name(td_file, cd_name, id_start, id_count)
1330   !-------------------------------------------------------------------
1331   !> @brief This subroutine fill variable value in an opened
1332   !> netcdf file, given variable name or standard name.
1333   !> @details
1334   !> Optionaly, start indices and number of indices selected along each dimension
1335   !> could be specify in a 4 dimension array (/'x','y','z','t'/)
1336   !>
1337   !> look first for variable name. If it doesn't
1338   !> exist in file, look for variable standard name.<br/>
1339   !>
1340   !> @author J.Paul
1341   !> @date November, 2013 - Initial Version
1342   !>
1343   !> @param[inout] td_file   file structure
1344   !> @param[in] cd_name      variable name or standard name
1345   !> @param[in] id_start     index in the variable from which the data values will be read
1346   !> @param[in] id_count     number of indices selected along each dimension
1347   !-------------------------------------------------------------------
1348
1349      IMPLICIT NONE
1350
1351      ! Argument     
1352      TYPE(TFILE),                   INTENT(INOUT) :: td_file
1353      CHARACTER(LEN=*),              INTENT(IN)    :: cd_name
1354      INTEGER(i4),     DIMENSION(:), INTENT(IN),  OPTIONAL :: id_start
1355      INTEGER(i4),     DIMENSION(:), INTENT(IN),  OPTIONAL :: id_count
1356
1357      ! local variable
1358      INTEGER(i4)       :: il_varid
1359      !----------------------------------------------------------------
1360      ! check if file opened
1361      IF( td_file%i_id == 0 )THEN
1362
1363         CALL logger_error( &
1364         &  "IOM CDF FILL VAR: no id associated to file "//TRIM(td_file%c_name))
1365
1366      ELSE
1367
1368            il_varid=var_get_index(td_file%t_var(:), cd_name)
1369            IF( il_varid /= 0 )THEN
1370
1371               !!! read variable value
1372               CALL iom_cdf__read_var_value(td_file, td_file%t_var(il_varid), &
1373               &                            id_start, id_count)
1374
1375            ELSE
1376
1377               CALL logger_error( &
1378               &  "IOM CDF FILL VAR: there is no variable with "//&
1379               &  "name or standard name"//TRIM(cd_name)//&
1380               &  " in file "//TRIM(td_file%c_name))
1381            ENDIF
1382
1383      ENDIF
1384     
1385   END SUBROUTINE iom_cdf__fill_var_name
1386   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1387   FUNCTION iom_cdf__read_var_meta(td_file, id_varid) &
1388         & RESULT (tf_var)
1389   !-------------------------------------------------------------------
1390   !> @brief This function read metadata of a variable in an opened
1391   !> netcdf file.
1392   !>
1393   !> @note variable value are not read
1394   !>
1395   !> @author J.Paul
1396   !> @date November, 2013 - Initial Version
1397   !> @date September, 2014
1398   !> - force to use FillValue=1.e20 if no FillValue for coordinate variable.
1399   !> @date September, 2015
1400   !> - manage useless (dummy) attribute
1401   !>
1402   !> @param[in] td_file   file structure
1403   !> @param[in] id_varid  variable id
1404   !> @return  variable structure
1405   !-------------------------------------------------------------------
1406
1407      IMPLICIT NONE
1408
1409      ! Argument     
1410      TYPE(TFILE), INTENT(IN) :: td_file
1411      INTEGER(i4), INTENT(IN) :: id_varid
1412
1413      ! function
1414      TYPE(TVAR)              :: tf_var
1415
1416      ! local variable
1417      CHARACTER(LEN=lc)                                       :: cl_name
1418
1419      INTEGER(i4)                                             :: il_status
1420      INTEGER(i4)                                             :: il_type
1421      INTEGER(i4)                                             :: il_ndim
1422      INTEGER(i4)                                             :: il_natt
1423      INTEGER(i4)                                             :: il_attid
1424      INTEGER(i4), DIMENSION(NF90_MAX_VAR_DIMS)               :: il_dimid
1425
1426      TYPE(TDIM) , DIMENSION(ip_maxdim)                       :: tl_dim
1427      TYPE(TATT)                                              :: tl_fill
1428      TYPE(TATT) , DIMENSION(:)                 , ALLOCATABLE :: tl_att
1429      TYPE(TATT) , DIMENSION(:)                 , ALLOCATABLE :: tl_tmp
1430
1431      ! loop indices
1432      INTEGER(i4) :: ji
1433      !----------------------------------------------------------------
1434      ! check if file opened
1435      IF( td_file%i_id == 0 )THEN
1436
1437         CALL logger_error( &
1438            &  " IOM CDF READ VAR META: no id associated to file "//&
1439            &   TRIM(td_file%c_name))
1440
1441      ELSE
1442
1443         ! inquire variable
1444         CALL logger_debug( &
1445            &  " IOM CDF READ VAR META: inquire variable "//&
1446            &  TRIM(fct_str(id_varid))//&
1447            &  " in file "//TRIM(td_file%c_name))
1448         
1449         il_dimid(:)=0
1450
1451         il_status=NF90_INQUIRE_VARIABLE( td_file%i_id, id_varid,        &
1452            &                             cl_name,    &
1453            &                             il_type,    &
1454            &                             il_ndim,    &
1455            &                             il_dimid(:),&
1456            &                             il_natt )
1457         CALL iom_cdf__check(il_status,"IOM CDF READ VAR META: ")
1458
1459         !!! fill variable dimension structure
1460         tl_dim(:)=iom_cdf__read_var_dim( td_file, il_ndim, cl_name, il_dimid(:) )
1461
1462         IF( il_natt /= 0 )THEN
1463            ALLOCATE( tl_att(il_natt) )
1464            !!! fill variable attribute structure
1465            tl_att(:)=iom_cdf__read_var_att(td_file, id_varid, il_natt)
1466
1467            !! look for _FillValue. if none add one
1468            il_attid=att_get_id(tl_att(:),'_FillValue')
1469            IF( il_attid == 0 )THEN
1470               CALL logger_info("IOM CDF READ VAR META: no _FillValue for variable "//&
1471                  &  TRIM(cl_name)//" in file "//TRIM(td_file%c_name) )
1472
1473               il_attid=att_get_id(tl_att(:),'missing_value')
1474               IF( il_attid /= 0 )THEN
1475                  ! create attribute _FillValue
1476                  CALL logger_info("IOM CDF READ VAR META: assume _FillValue is equal to "//&
1477                     &             "missing_value for variable "//TRIM(cl_name) )
1478                  tl_fill=att_init('_FillValue',tl_att(il_attid)%d_value(:), &
1479                     &             id_type=tl_att(il_attid)%i_type)
1480               ELSE
1481                  ! create attribute _FillValue
1482                  SELECT CASE(TRIM(fct_lower(cl_name)))
1483                     CASE DEFAULT
1484                        CALL logger_info("IOM CDF READ VAR META: assume _FillValue is equal to "//&
1485                           &             "zero for variable "//TRIM(cl_name) )
1486                        tl_fill=att_init('_FillValue',0.)
1487                     CASE('nav_lon','nav_lat', 'nav_lev', &
1488                        &  'glamt','glamu','glamv','glamf', &
1489                        &  'gphit','gphiu','gphiv','gphif')
1490                        CALL logger_info("IOM CDF READ VAR META: assume _FillValue is equal to "//&
1491                           &             "dummy fillValue (1.e20) for variable "//TRIM(cl_name) )
1492                        tl_fill=att_init('_FillValue',1.e20)
1493                  END SELECT
1494               ENDIF
1495
1496               ALLOCATE( tl_tmp(il_natt) )
1497               ! save read attribut
1498               tl_tmp(:)=att_copy(tl_att(:))
1499               ! change number of attribute in array
1500               CALL att_clean(tl_att(:))
1501               DEALLOCATE( tl_att )
1502               ALLOCATE( tl_att(il_natt+1) )
1503               ! copy read attribut
1504               tl_att(1:il_natt)=att_copy(tl_tmp(:))
1505               ! clean
1506               CALL att_clean(tl_tmp(:))
1507               DEALLOCATE( tl_tmp )
1508
1509               ! create attribute _FillValue
1510               tl_att(il_natt+1)=att_copy(tl_fill)
1511
1512            ENDIF
1513
1514         ELSE
1515            ALLOCATE(tl_att(il_natt+1) )
1516            ! create attribute _FillValue
1517            SELECT CASE(TRIM(fct_lower(cl_name)))
1518               CASE DEFAULT
1519                  CALL logger_info("IOM CDF READ VAR META: assume _FillValue is equal to "//&
1520                     &             "zero for variable "//TRIM(cl_name) )
1521                  tl_fill=att_init('_FillValue',0.)
1522               CASE('nav_lon','nav_lat', &
1523                  &  'glamt','glamu','glamv','glamf', &
1524                  &  'gphit','gphiu','gphiv','gphif')
1525                  CALL logger_info("IOM CDF READ VAR META: assume _FillValue is equal to "//&
1526                     &             "dummy fillValue (1.e20) for variable "//TRIM(cl_name) )
1527                  tl_fill=att_init('_FillValue',1.e20)
1528            END SELECT           
1529            ! create attribute _FillValue
1530            tl_att(il_natt+1)=att_copy(tl_fill)
1531         ENDIF
1532
1533         !! initialize variable
1534         tf_var=var_init( cl_name, il_type, tl_dim(:), tl_att(:), id_id=id_varid )
1535
1536         !! look for dummy attribute
1537         DO ji=il_natt,1,-1
1538            IF( att_is_dummy(tl_att(ji)) )THEN
1539               CALL var_del_att(tf_var, tl_att(ji))
1540            ENDIF
1541         ENDDO
1542
1543         !! check if variable is dummy
1544         IF( var_is_dummy(tf_var) )THEN
1545            tf_var%l_use=.FALSE.
1546         ENDIF
1547
1548         !! check if all dimensions are allowed
1549         DO ji=1,il_ndim
1550            IF( ALL(td_file%t_dim(:)%i_id /= il_dimid(ji)) )THEN
1551               tf_var%l_use=.FALSE.
1552            ENDIF
1553         ENDDO
1554
1555         ! clean
1556         CALL dim_clean(tl_dim(:))
1557         CALL att_clean(tl_fill)
1558         CALL att_clean(tl_att(:))
1559         DEALLOCATE( tl_att )
1560
1561      ENDIF
1562
1563   END FUNCTION iom_cdf__read_var_meta
1564   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1565   FUNCTION iom_cdf__read_var_dim(td_file, id_ndim, cd_name, id_dimid) &
1566         & RESULT (tf_dim)
1567   !-------------------------------------------------------------------
1568   !> @brief This subroutine read variable dimension
1569   !> in an opened netcdf file.
1570   !>
1571   !> @details
1572   !> the number of dimension can't exceed 4,
1573   !> and should be 'x', 'y', 'z', 't' (whatever their order).<br/>
1574   !> If the number of dimension read is less than 4, the array of dimension
1575   !> strucure is filled with unused dimension.<br/>
1576   !> So the array of dimension structure of a variable is always compose of 4
1577   !> dimension (use or not).
1578   !>
1579   !> @warn dummy dimension are not used.
1580   !>
1581   !> @author J.Paul
1582   !> @date November, 2013 - Initial Version
1583   !> @date July, 2015
1584   !> - Bug fix: use order to disorder table (see dim_init)
1585   !> @date September, 2015
1586   !> - check dummy dimension
1587   !> @date April, 2018
1588   !> - Bug fix : use dimid to look for the index of the dimension, and not as
1589   !> dimension index
1590   !>
1591   !> @param[in] td_file   file structure
1592   !> @param[in] id_ndim   number of dimension
1593   !> @param[in] cd_name   variable name
1594   !> @param[in] id_dimid  array of dimension id
1595   !> @return array dimension structure
1596   !-------------------------------------------------------------------
1597
1598      IMPLICIT NONE
1599
1600      ! Argument     
1601      TYPE(TFILE),               INTENT(IN) :: td_file
1602      INTEGER(i4),               INTENT(IN) :: id_ndim
1603      CHARACTER(LEN=*)         , INTENT(IN) :: cd_name
1604      INTEGER(i4), DIMENSION(:), INTENT(IN) :: id_dimid
1605
1606      ! function
1607      TYPE(TDIM), DIMENSION(ip_maxdim)      :: tf_dim
1608
1609      ! local variable
1610      INTEGER(i4), DIMENSION(ip_maxdim) :: il_xyzt2
1611      INTEGER(i4)                       :: il_idx
1612
1613      TYPE(TDIM) , DIMENSION(ip_maxdim) :: tl_dim
1614
1615      ! loop indices
1616      INTEGER(i4) :: ji
1617      INTEGER(i4) :: jj
1618      INTEGER(i4) :: ii
1619      !----------------------------------------------------------------
1620
1621      IF( id_ndim == 0 )THEN
1622
1623         tl_dim(:)%l_use=.FALSE.
1624
1625         ! reorder dimension to ('x','y','z','t')
1626         CALL dim_reorder(tl_dim(:))
1627
1628         tf_dim(:)=dim_copy(tl_dim(:))
1629
1630         ! clean
1631         CALL dim_clean(tl_dim(:))
1632
1633      ELSE IF( id_ndim > 0 )THEN
1634
1635         ii=1
1636         DO ji = 1, id_ndim
1637
1638            ! check if dimension to be used, is allowed
1639            IF( ANY(td_file%t_dim(:)%i_id == id_dimid(ji)) )THEN
1640               IF( ii > ip_maxdim )THEN
1641                  CALL logger_error(" IOM CDF READ VAR DIM: "//&
1642                     &  "too much dimensions for variable "//&
1643                     &  TRIM(cd_name)//". check dummy configuration file.")
1644               ENDIF
1645
1646               CALL logger_debug( " IOM CDF READ VAR DIM: get variable "//&
1647                  &  "dimension "//TRIM(fct_str(ji)) )
1648
1649               ! look for dimension index
1650               DO jj=1,ip_maxdim
1651                  IF( td_file%t_dim(jj)%i_id == id_dimid(ji) )THEN
1652                     il_idx=jj
1653                     EXIT
1654                  ENDIF
1655               ENDDO
1656               !il_xyzt2(ii)=td_file%t_dim(id_dimid(ji))%i_xyzt2
1657               il_xyzt2(ii)=td_file%t_dim(il_idx)%i_xyzt2
1658
1659               ! read dimension information
1660               tl_dim(ii) = dim_init( td_file%t_dim(il_xyzt2(ii))%c_name, &
1661                  &                   td_file%t_dim(il_xyzt2(ii))%i_len )
1662
1663               ii=ii+1
1664            ELSE
1665               CALL logger_debug(" IOM CDF READ VAR DIM: dummy variable "//&
1666                  &              "dimension "//TRIM(fct_str(ji))//" not used.")
1667            ENDIF
1668         ENDDO
1669
1670         ! reorder dimension to ('x','y','z','t')
1671         CALL dim_reorder(tl_dim(:))
1672 
1673         tf_dim(:)=dim_copy(tl_dim(:))
1674
1675         ! clean
1676         CALL dim_clean(tl_dim(:))
1677
1678      ENDIF
1679
1680   END FUNCTION iom_cdf__read_var_dim
1681   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1682   FUNCTION iom_cdf__read_var_att(td_file, id_varid, id_natt) &
1683         & RESULT (tf_att)
1684   !-------------------------------------------------------------------
1685   !> @brief This subroutine read variable attributes
1686   !> in an opened netcdf file.
1687   !>
1688   !> @author J.Paul
1689   !> @date November, 2013 - Initial Version
1690   !>
1691   !> @param[in] td_file   file structure
1692   !> @param[in] id_varid  variable id
1693   !> @param[in] id_natt   number of attributes
1694   !> @return array of attribute structure
1695   !-------------------------------------------------------------------
1696
1697      IMPLICIT NONE
1698
1699      ! Argument     
1700      TYPE(TFILE), INTENT(IN) :: td_file
1701      INTEGER(i4), INTENT(IN) :: id_varid
1702      INTEGER(i4), INTENT(IN) :: id_natt     
1703
1704      ! function
1705      TYPE(TATT), DIMENSION(id_natt) :: tf_att
1706
1707      ! local variable
1708
1709      ! loop indices
1710      INTEGER(i4) :: ji
1711      !----------------------------------------------------------------
1712
1713      IF( id_natt > 0 )THEN
1714     
1715         ! read attributes
1716         DO ji = 1, id_natt
1717            CALL logger_trace( " IOM CDF READ VAR ATT: get attribute "//&
1718               &               TRIM(fct_str(ji))//" in file "//&
1719               &               TRIM(td_file%c_name) &
1720            )
1721
1722            tf_att(ji)=iom_cdf_read_att(td_file, id_varid, ji)
1723
1724         ENDDO
1725
1726      ELSE
1727
1728         CALL logger_debug( " IOM CDF READ VAR ATT: no attribute for variable " )
1729
1730      ENDIF
1731
1732   END FUNCTION iom_cdf__read_var_att
1733   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1734   SUBROUTINE iom_cdf__read_var_value(td_file, td_var, id_start, id_count)
1735   !-------------------------------------------------------------------
1736   !> @brief This subroutine read variable value
1737   !> in an opened netcdf file.
1738   !> @details
1739   !> Optionaly, start indices and number of indices selected along each dimension
1740   !> could be specify in a 4 dimension array (/'x','y','z','t'/)
1741   !>
1742   !> @author J.Paul
1743   !> @date November, 2013 - Initial Version
1744   !> @date June, 2015
1745   !> - use scale factor and offset, as soon as read variable value
1746   !> @date January, 2019
1747   !> - read array in netcdf file, level by level, and time step by time step
1748   !> - apply scale factor and offset, level by level, and time step by time step
1749   !>
1750   !> @param[in] td_file   file structure
1751   !> @param[inout] td_var variable structure
1752   !> @param[in] id_start  index in the variable from which the data values will be read
1753   !> @param[in] id_count  number of indices selected along each dimension
1754   !> @return variable structure completed
1755   !-------------------------------------------------------------------
1756
1757      IMPLICIT NONE
1758
1759      ! Argument     
1760      TYPE(TFILE),               INTENT(IN)    :: td_file
1761      TYPE(TVAR) ,               INTENT(INOUT) :: td_var
1762      INTEGER(i4), DIMENSION(:), INTENT(IN),   OPTIONAL :: id_start
1763      INTEGER(i4), DIMENSION(:), INTENT(IN),   OPTIONAL :: id_count
1764
1765      ! local variable
1766      INTEGER(i4)                                    :: il_status
1767      INTEGER(i4)                                    :: il_tmp1
1768      INTEGER(i4)                                    :: il_tmp2
1769      INTEGER(i4)                                    :: il_varid
1770      INTEGER(i4), DIMENSION(ip_maxdim)              :: il_start
1771      INTEGER(i4), DIMENSION(ip_maxdim)              :: il_count
1772      INTEGER(i4), DIMENSION(ip_maxdim)              :: il_start_ord
1773      INTEGER(i4), DIMENSION(ip_maxdim)              :: il_count_ord
1774      INTEGER(i4), DIMENSION(ip_maxdim)              :: il_start_tmp
1775      INTEGER(i4), DIMENSION(ip_maxdim)              :: il_count_tmp
1776
1777      REAL(dp)   , DIMENSION(:,:,:,:)  , ALLOCATABLE :: dl_value
1778      REAL(dp)   , DIMENSION(:,:,:,:)  , ALLOCATABLE :: dl_tmp
1779
1780      ! loop indices
1781      INTEGER(i4) :: ji
1782      INTEGER(i4) :: jj
1783      INTEGER(i4) :: jk
1784      INTEGER(i4) :: jl
1785      !----------------------------------------------------------------
1786
1787      ! check if variable in file structure
1788      il_varid=var_get_id(td_file%t_var(:),TRIM(td_var%c_name))
1789      IF( il_varid /= 0 )THEN
1790
1791         ! check id_count and id_start optionals parameters...
1792         IF( (       PRESENT(id_start)  .AND. (.NOT. PRESENT(id_count))) .OR. &
1793             ((.NOT. PRESENT(id_start)) .AND.        PRESENT(id_count) ) )THEN
1794            CALL logger_warn( "IOM CDF READ VAR VALUE: id_start and id_count"//&
1795               & " should be both specify")
1796         ENDIF
1797         IF( PRESENT(id_start).AND.PRESENT(id_count) )THEN
1798
1799            IF( SIZE(id_start(:)) /= ip_maxdim .OR. &
1800            &   SIZE(id_count(:)) /= ip_maxdim )THEN
1801               CALL logger_error("IOM CDF READ VAR: dimension of array start"//&
1802                  &  " or count are invalid to read variable "//&
1803                  &  TRIM(td_var%c_name)//" in file "//TRIM(td_file%c_name) )
1804            ENDIF
1805
1806            ! change dimension order from ('x','y','z','t')
1807            il_start(:)=dim_reorder_xyzt2(td_var%t_dim, id_start(:))
1808            il_count(:)=dim_reorder_xyzt2(td_var%t_dim, id_count(:))
1809
1810            ! keep ordered array ('x','y','z','t')
1811            il_start_ord(:)=id_start(:)
1812            il_count_ord(:)=id_count(:)
1813
1814         ELSE
1815
1816            ! change dimension order from ('x','y','z','t')
1817            il_start(:)=(/1,1,1,1/)
1818            il_count(:)=dim_reorder_xyzt2(td_var%t_dim(:),td_var%t_dim(:)%i_len)
1819
1820            ! keep ordered array ('x','y','z','t')
1821            il_start_ord(:)=(/1,1,1,1/)
1822            il_count_ord(:)=td_var%t_dim(:)%i_len
1823
1824         ENDIF
1825
1826         ! check dimension
1827         IF( .NOT. ALL(il_start_ord(:)>=(/1,1,1,1/)) )THEN
1828
1829            CALL logger_error( "IOM CDF READ VAR VALUE: start indices should"//&
1830            &  " be greater than or equal to 1")
1831
1832         ENDIF
1833
1834         IF(.NOT.ALL(il_start_ord(:)+il_count_ord(:)-1 <= &
1835            &  (/td_var%t_dim( 1 )%i_len,&
1836            &    td_var%t_dim( 2 )%i_len,&
1837            &    td_var%t_dim( 3 )%i_len,&
1838            &    td_var%t_dim( 4 )%i_len &
1839            &                                            /)) )THEN
1840
1841            DO ji = 1, ip_maxdim
1842               il_tmp1=il_start_ord(ji)+il_count_ord(ji)-1
1843               il_tmp2=td_var%t_dim(ji)%i_len
1844               CALL logger_debug( "IOM CDF READ VAR VALUE: start + count -1:"//&
1845               &  TRIM(fct_str(il_tmp1))//" variable dimension"//&
1846               &  TRIM(fct_str(il_tmp2)))
1847            ENDDO
1848            CALL logger_error( "IOM CDF READ VAR VALUE: start + count exceed "//&
1849            &  "variable dimension for "//TRIM(td_var%c_name) )
1850
1851         ELSE
1852
1853            ! Allocate space to hold variable value (disorder)
1854            ALLOCATE(dl_value( il_count(1), &
1855               &               il_count(2), &
1856               &               il_count(3), &
1857               &               il_count(4)),&
1858               &               stat=il_status)
1859            IF( il_status /= 0 )THEN
1860
1861              CALL logger_fatal( &
1862               &  "IOM CDF READ VAR VALUE: not enough space to put variable "//&
1863               &  TRIM(td_var%c_name))
1864
1865            ENDIF
1866
1867            ! read values
1868            CALL logger_debug( &
1869            &  "IOM CDF READ VAR VALUE: get variable "//TRIM(td_var%c_name)//&
1870            &  " in file "//TRIM(td_file%c_name))
1871
1872            il_start_tmp(:)=il_start(:)
1873            il_count_tmp(:)=il_count(:)
1874            DO jl=il_start(jp_L),il_start(jp_L)+il_count(jp_L)-1
1875               il_start_tmp(jp_L)=jl
1876               il_count_tmp(jp_L) = 1
1877               DO jk=il_start(jp_K),il_start(jp_K)+il_count(jp_K)-1
1878                  il_start_tmp(jp_K)=jk
1879                  il_count_tmp(jp_K)=1 
1880                  il_status = NF90_GET_VAR( td_file%i_id, il_varid,           &
1881                  &                                       dl_value(:,:,jk,jl),&
1882                  &                                       start = il_start_tmp(:),&
1883                  &                                       count = il_count_tmp(:) )
1884            !   il_status = NF90_GET_VAR( td_file%i_id, il_varid,           &
1885            !   &                                       dl_value(:,:,:,:),&
1886               CALL iom_cdf__check(il_status,"IOM CDF READ VAR VALUE: ")
1887               ENDDO
1888            ENDDO
1889
1890            ! Allocate space to hold variable value in structure
1891            IF( ASSOCIATED(td_var%d_value) )THEN
1892               DEALLOCATE(td_var%d_value)   
1893            ENDIF
1894 
1895            ! new dimension length
1896            td_var%t_dim(:)%i_len=il_count_ord(:)
1897
1898!>   dummy patch for pgf95
1899            ALLOCATE( dl_tmp( td_var%t_dim(1)%i_len, &
1900            &                 td_var%t_dim(2)%i_len, &
1901            &                 td_var%t_dim(3)%i_len, &
1902            &                 td_var%t_dim(4)%i_len),&
1903            &        stat=il_status)
1904            IF(il_status /= 0 )THEN
1905
1906               CALL logger_fatal( &
1907               &  "IOM CDF READ VAR VALUE: not enough space to put variable "//&
1908               &  TRIM(td_var%c_name)//&
1909               &  " in variable structure")
1910            ENDIF
1911            dl_tmp(:,:,:,:)=td_var%d_fill
1912
1913            ! reshape values to be ordered as ('x','y','z','t')
1914            dl_tmp(:,:,:,:)=dim_reshape_2xyzt(td_var%t_dim(:), &
1915            &                                 dl_value(:,:,:,:))
1916
1917            DEALLOCATE(dl_value)
1918
1919            ALLOCATE(td_var%d_value( td_var%t_dim(1)%i_len, &
1920            &                        td_var%t_dim(2)%i_len, &
1921            &                        td_var%t_dim(3)%i_len, &
1922            &                        td_var%t_dim(4)%i_len),&
1923            &        stat=il_status)
1924            IF(il_status /= 0 )THEN
1925
1926               CALL logger_fatal( &
1927               &  "IOM CDF READ VAR VALUE: not enough space to put variable "//&
1928               &  TRIM(td_var%c_name)//&
1929               &  " in variable structure")
1930
1931            ENDIF
1932!            ! FillValue by default
1933!            td_var%d_value(:,:,:,:)=td_var%d_fill
1934!
1935!            ! reshape values to be ordered as ('x','y','z','t')
1936!            td_var%d_value(:,:,:,:)=dim_reshape_2xyzt(td_var%t_dim(:), &
1937!            &                                         dl_value(:,:,:,:))
1938!
1939!            DEALLOCATE(dl_value)
1940
1941            DO jl=1,td_var%t_dim(jp_L)%i_len
1942               DO jk=1,td_var%t_dim(jp_K)%i_len
1943                  DO jj=1,td_var%t_dim(jp_J)%i_len
1944                     DO ji=1,td_var%t_dim(jp_I)%i_len
1945                        td_var%d_value(ji,jj,jk,jl)=dl_tmp(ji,jj,jk,jl)
1946                     ENDDO
1947                  ENDDO
1948               ENDDO
1949            ENDDO
1950            DEALLOCATE(dl_tmp)
1951!<   dummy patch for pgf95
1952
1953            ! force to change _FillValue to avoid mistake
1954            ! with dummy zero _FillValue
1955            IF( td_var%d_fill == 0._dp )THEN
1956               CALL var_chg_FillValue(td_var)
1957            ENDIF
1958
1959            ! use scale factor and offset
1960            DO jl=1,td_var%t_dim(jp_L)%i_len
1961               DO jk=1,td_var%t_dim(jp_K)%i_len
1962                  WHERE( td_var%d_value(:,:,jk,jl) /= td_var%d_fill )
1963                     td_var%d_value(:,:,jk,jl) = &
1964                     &  td_var%d_value(:,:,jk,jl)*td_var%d_scf + td_var%d_ofs
1965                  END WHERE
1966               ENDDO
1967            ENDDO
1968
1969         ENDIF
1970      ELSE
1971         CALL logger_error( &
1972         &  "IOM CDF READ VAR VALUE: no variable "//TRIM(td_var%c_name)//&
1973         &  " in file structure "//TRIM(td_file%c_name))
1974      ENDIF
1975
1976   END SUBROUTINE iom_cdf__read_var_value
1977   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1978   SUBROUTINE iom_cdf_write_header(td_file, cd_dimorder, td_dim)
1979   !-------------------------------------------------------------------
1980   !> @brief This subroutine write file header in an opened netcdf file.
1981   !>
1982   !> @details
1983   !> optionally, you could specify dimension order (default 'xyzt'),
1984   !> and/or dimension structure to be used.
1985   !>
1986   !> @author J.Paul
1987   !> @date November, 2013 - Initial Version
1988   !> @date July, 2015
1989   !> - add dimension order option
1990   !> @date August, 2017
1991   !> - split write_file into write_header and write_var
1992   !> - add dimension structure as optional argument
1993   !> @date September, 2017
1994   !> - do not check variable dimension if dimension forced
1995   !>
1996   !> @param[inout] td_file   file structure
1997   !-------------------------------------------------------------------
1998
1999      IMPLICIT NONE
2000
2001      ! Argument     
2002      TYPE(TFILE)                           , INTENT(INOUT) :: td_file
2003      CHARACTER(LEN=*)                      , INTENT(IN   ), OPTIONAL :: cd_dimorder
2004      TYPE(TDIM )     , DIMENSION(ip_maxdim), INTENT(IN   ), OPTIONAL :: td_dim
2005
2006      ! local variable
2007      INTEGER(i4)                            :: il_status
2008      INTEGER(i4), DIMENSION(:), ALLOCATABLE :: il_value
2009
2010      CHARACTER(LEN=lc)                      :: cl_dimorder
2011
2012      LOGICAL                                :: ll_chkdim
2013
2014      TYPE(TVAR), DIMENSION(ip_maxdim)       :: tl_var
2015
2016      TYPE(TDIM), DIMENSION(ip_maxdim)       :: tl_dim
2017
2018      ! loop indices
2019      INTEGER(i4) :: ji
2020      INTEGER(i4) :: jj
2021      INTEGER(i4) :: jvar
2022      !----------------------------------------------------------------
2023
2024      cl_dimorder='xyzt'
2025      IF( PRESENT(cd_dimorder) ) cl_dimorder=TRIM(cd_dimorder)
2026
2027      ! check if file opened
2028      IF( td_file%i_id == 0 )THEN
2029
2030         CALL logger_error( &
2031         &  " IOM CDF WRITE HEADER: no id associated to file "//&
2032         &  TRIM(td_file%c_name)//". Check if file is opened.")
2033
2034      ELSE
2035         IF( td_file%l_wrt )THEN
2036
2037            ! remove dummy variable
2038            CALL file_del_var(td_file,'no0d')
2039            CALL file_del_var(td_file,'no1d')
2040            CALL file_del_var(td_file,'no2d')
2041            CALL file_del_var(td_file,'no3d')
2042
2043            DO ji = 1, td_file%i_nvar
2044               CALL var_check_dim( td_file%t_var(ji) )
2045            ENDDO
2046
2047            IF( PRESENT(td_dim) )THEN
2048               ! special case to rebuild mpp layout
2049               DO ji=1,ip_maxdim
2050                  IF( td_dim(ji)%l_use ) CALL file_move_dim(td_file, td_dim(ji))
2051               ENDDO
2052            ELSE
2053               ! save useful dimension
2054               IF( ASSOCIATED(td_file%t_var) )THEN
2055                  tl_dim(:) = var_max_dim( td_file%t_var(:) )
2056
2057                  DO ji=1,ip_maxdim
2058                     IF( tl_dim(ji)%l_use ) CALL file_move_dim(td_file, tl_dim(ji))
2059                  ENDDO
2060                  ! clean
2061                  CALL dim_clean(tl_dim(:))
2062               ENDIF
2063            ENDIF
2064
2065            ! Enter in define mode
2066            IF( .NOT. td_file%l_def )THEN
2067               CALL logger_debug( &
2068               &  " IOM CDF WRITE HEADER: Enter define mode, file "//&
2069               &  TRIM(td_file%c_name))
2070
2071               ! Enter define mode
2072               il_status=NF90_REDEF(td_file%i_id)
2073               CALL iom_cdf__check(il_status,"IOM CDF WRITE HEADER: ")
2074
2075               td_file%l_def=.TRUE.
2076            ENDIF           
2077
2078            ! write dimension definition in header of file
2079            IF( TRIM(cl_dimorder) /= 'xyzt' )THEN
2080               CALL dim_reorder(td_file%t_dim(:),TRIM(cl_dimorder))
2081               DO jvar=1,td_file%i_nvar
2082                  CALL logger_debug("VAR REORDER: "//TRIM(td_file%t_var(jvar)%c_name))
2083                  CALL var_reorder(td_file%t_var(jvar),TRIM(cl_dimorder))
2084               ENDDO
2085            ENDIF
2086
2087            ! write dimension in file
2088            DO ji = 1, ip_maxdim
2089               IF( td_file%t_dim(ji)%l_use )THEN
2090                  CALL iom_cdf__write_dim_def(td_file, td_file%t_dim(ji))
2091
2092                  ! write dimension variable
2093                  ALLOCATE(il_value(td_file%t_dim(ji)%i_len))
2094                  il_value(:)=(/(jj,jj=1,td_file%t_dim(ji)%i_len)/)
2095
2096                  tl_var(ji)=var_init( fct_upper(td_file%t_dim(ji)%c_sname), &
2097                  &                    il_value(:),                          &
2098                  &                    td_dim=td_file%t_dim(ji) )
2099
2100                  DEALLOCATE(il_value)
2101
2102                  ! do not use FillValue for dimension variable
2103                  CALL var_del_att(tl_var(ji), "_FillValue")
2104                   
2105                  ! write dimension variable definition in header of file
2106                  CALL iom_cdf__write_var_def(td_file,tl_var(ji))
2107                 
2108               ENDIF
2109            ENDDO
2110
2111            ! write global attibute in file
2112            DO ji = 1, td_file%i_natt
2113               CALL iom_cdf__write_att_def(td_file, NF90_GLOBAL, td_file%t_att(ji))
2114            ENDDO
2115 
2116            ! write variable definition in header of file
2117            ll_chkdim=.TRUE.
2118            IF( PRESENT(td_dim) )THEN
2119               ! special case to rebuild mpp layout
2120               ! do not check dimension length of variable
2121               ll_chkdim=.FALSE.
2122            ENDIF
2123            DO ji=1,td_file%i_nvar
2124               CALL iom_cdf__write_var_def(td_file, td_file%t_var(ji),ll_chkdim)
2125            ENDDO
2126
2127            ! Leave define mode
2128            IF( td_file%l_def )THEN
2129               CALL logger_debug( &
2130               &  " IOM CDF WRITE HEADER: Leave define mode, file "//&
2131               &  TRIM(td_file%c_name))
2132
2133               ! Leave define mode
2134               il_status=NF90_ENDDEF(td_file%i_id)
2135               CALL iom_cdf__check(il_status,"IOM CDF WRITE HEADER: ")
2136
2137               td_file%l_def=.FALSE.
2138            ENDIF
2139
2140            ! write dimension variable in file
2141            DO ji = 1, ip_maxdim
2142               IF( td_file%t_dim(ji)%l_use )THEN
2143                  ! do not use FillValue for dimension variable
2144                  CALL var_del_att(tl_var(ji), "_FillValue")
2145                  CALL iom_cdf__write_var(td_file,tl_var(ji))
2146               ENDIF
2147            ENDDO
2148            ! clean
2149            CALL var_clean(tl_var(:))
2150
2151         ELSE
2152
2153            CALL logger_error( &
2154            &  "IOM CDF WRITE HEADER: try to write in file "//TRIM(td_file%c_name)//&
2155            &  ", not opened in write mode")
2156
2157         ENDIF
2158      ENDIF
2159
2160   END SUBROUTINE iom_cdf_write_header
2161   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2162   SUBROUTINE iom_cdf_write_var(td_file, cd_dimorder, id_start, id_count)!, ld_first)
2163   !-------------------------------------------------------------------
2164   !> @brief This subroutine write variable(s) in an opened netcdf file.
2165   !>
2166   !> @details
2167   !> optionally, you could specify dimension order (default 'xyzt')
2168   !>
2169   !> @author J.Paul
2170   !> @date November, 2013 - Initial Version
2171   !> @date July, 2015
2172   !> - add dimension order option
2173   !> @date August, 2017
2174   !> - add start and count array as optional argument
2175   !>
2176   !> @param[inout] td_file   file structure
2177   !> @param[in] td_var       array of variable structure
2178   !> @param[in] cd_dimorder  dimension order
2179   !> @param[in] id_start     index in the variable from which the data values
2180   !> will be read
2181   !> @param[in] id_count     number of indices selected along each dimension
2182   !-------------------------------------------------------------------
2183
2184      IMPLICIT NONE
2185
2186      ! Argument     
2187      TYPE(TFILE)                   , INTENT(INOUT) :: td_file
2188      CHARACTER(LEN=*)              , INTENT(IN   ), OPTIONAL :: cd_dimorder
2189      INTEGER(i4)     , DIMENSION(:), INTENT(IN   ), OPTIONAL :: id_start
2190      INTEGER(i4)     , DIMENSION(:), INTENT(IN   ), OPTIONAL :: id_count
2191      ! local variable
2192      CHARACTER(LEN=lc)                      :: cl_dimorder
2193
2194      ! loop indices
2195      INTEGER(i4) :: ji
2196      !----------------------------------------------------------------
2197
2198      cl_dimorder='xyzt'
2199      IF( PRESENT(cd_dimorder) ) cl_dimorder=TRIM(cd_dimorder)
2200
2201      ! check if file opened
2202      IF( td_file%i_id == 0 )THEN
2203
2204         CALL logger_error( &
2205         &  " IOM CDF WRITE VAR: no id associated to file "//TRIM(td_file%c_name))
2206
2207      ELSE
2208         IF( td_file%l_wrt )THEN
2209
2210            ! write variable in file
2211            DO ji = 1, td_file%i_nvar
2212
2213               ! change dimension order
2214               IF( TRIM(cl_dimorder) /= 'xyzt' )THEN
2215                  CALL logger_debug("VAR REORDER: "//TRIM(td_file%t_var(ji)%c_name))
2216                  CALL var_reorder(td_file%t_var(ji),TRIM(cl_dimorder))
2217               ENDIF
2218
2219               IF( ASSOCIATED(td_file%t_var(ji)%d_value) )THEN
2220                  ! write variable value in file
2221                  CALL iom_cdf__write_var_value( td_file, td_file%t_var(ji), &
2222                     &                           id_start, id_count)
2223               ENDIF
2224            ENDDO
2225
2226         ELSE
2227
2228            CALL logger_error( &
2229            &  "IOM CDF WRITE VAR: try to write in file "//TRIM(td_file%c_name)//&
2230            &  ", not opened in write mode")
2231
2232         ENDIF
2233      ENDIF
2234
2235   END SUBROUTINE iom_cdf_write_var
2236   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2237   SUBROUTINE iom_cdf__write_dim_def(td_file, td_dim)
2238   !-------------------------------------------------------------------
2239   !> @brief This subroutine define a dimension in the header of a netcdf
2240   !> file.
2241   !>
2242   !> @author J.Paul
2243   !> @date November, 2013 - Initial Version
2244   !> @date August, 2017
2245   !> - rename in write_dim_def
2246   !> - do not check define mode here anymore
2247   !>
2248   !> @param[inout] td_file   file structure
2249   !> @param[inout] td_dim    dimension structure
2250   !-------------------------------------------------------------------
2251
2252      IMPLICIT NONE
2253
2254      ! Argument     
2255      TYPE(TFILE), INTENT(INOUT) :: td_file
2256      TYPE(TDIM),  INTENT(INOUT) :: td_dim
2257
2258      ! local variable
2259      INTEGER(i4) :: il_status
2260      !----------------------------------------------------------------
2261
2262      IF( td_dim%l_use )THEN
2263         IF( td_dim%l_uld )THEN
2264            ! write unlimited dimension
2265            CALL logger_trace( &
2266            &  "IOM CDF WRITE FILE DIM: write unlimited dimension "//&
2267            &  TRIM(td_dim%c_name)//" in file "//TRIM(td_file%c_name))
2268
2269            il_status=NF90_DEF_DIM(td_file%i_id, fct_upper(td_dim%c_sname), &
2270            &                      NF90_UNLIMITED, td_dim%i_id)
2271            CALL iom_cdf__check(il_status,"IOM CDF WRITE FILE DIM: ")
2272
2273         ELSE
2274            ! write not unlimited dimension
2275            CALL logger_debug( &
2276            &  "IOM CDF WRITE FILE DIM: write dimension "//TRIM(td_dim%c_name)//&
2277            &  " in file "//TRIM(td_file%c_name))
2278           
2279            CALL logger_debug("IOM CDF WRITE FILE DIM: id "//TRIM(fct_str(td_file%i_id))//&
2280               & " sname "//TRIM(td_dim%c_sname))
2281            il_status=NF90_DEF_DIM(td_file%i_id, fct_upper(td_dim%c_sname), &
2282            &                      td_dim%i_len, td_dim%i_id)
2283            CALL iom_cdf__check(il_status,"IOM CDF WRITE FILE DIM: ")
2284
2285         ENDIF
2286      ENDIF
2287
2288   END SUBROUTINE iom_cdf__write_dim_def
2289   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2290   SUBROUTINE iom_cdf__write_att_def(td_file, id_varid, td_att)
2291   !-------------------------------------------------------------------
2292   !> @brief This subroutine write a variable attribute in
2293   !> an opened netcdf file.
2294   !>
2295   !> @author J.Paul
2296   !> @date November, 2013 - Initial Version
2297   !> @date August, 2017
2298   !> - rename in write_att_def
2299   !> - do not check define mode here anymore
2300   !>
2301   !> @param[inout] td_file   file structure
2302   !> @param[in] id_varid     variable id. use NF90_GLOBAL to write
2303   !> global attribute in a file
2304   !> @param[in] td_att       attribute structure
2305   !-------------------------------------------------------------------
2306
2307      IMPLICIT NONE
2308
2309      ! Argument     
2310      TYPE(TFILE), INTENT(INOUT) :: td_file
2311      INTEGER(i4), INTENT(IN)    :: id_varid
2312      TYPE(TATT),  INTENT(IN)    :: td_att
2313
2314      ! local variable
2315      INTEGER(i4) :: il_status
2316      !----------------------------------------------------------------
2317
2318      !! put attribute value
2319      CALL logger_trace( &
2320      &  "IOM CDF WRITE FILE ATT: write attribute "//TRIM(td_att%c_name)//&
2321      &  " of variable "//TRIM(fct_str(id_varid))//&
2322      &  " in file "//TRIM(td_file%c_name))
2323      SELECT CASE( td_att%i_type )
2324
2325         CASE(NF90_CHAR)
2326            ! put the attribute
2327            il_status = NF90_PUT_ATT(td_file%i_id, id_varid, &
2328            &  td_att%c_name, td_att%c_value )
2329            CALL iom_cdf__check(il_status,"IOM CDF WRITE FILE ATT: ")
2330
2331         CASE(NF90_BYTE, NF90_SHORT, NF90_INT, NF90_FLOAT, NF90_DOUBLE)
2332            ! put the attribute
2333            il_status = NF90_PUT_ATT(td_file%i_id, id_varid, &
2334            &  td_att%c_name, td_att%d_value )
2335            CALL iom_cdf__check(il_status,"IOM CDF WRITE FILE ATT: ")
2336
2337      END SELECT
2338
2339   END SUBROUTINE iom_cdf__write_att_def
2340   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2341   SUBROUTINE iom_cdf__write_var(td_file, td_var, id_start, id_count)
2342   !-------------------------------------------------------------------
2343   !> @brief This subroutine write a variable in an opened netcdf file.
2344   !>
2345   !> @author J.Paul
2346   !> @date November, 2013 - Initial Version
2347   !> @date September, 2015
2348   !> - do not force to use zero as FillValue for any meshmask variable
2349   !> @date August, 2017
2350   !> - add start and count array as optional argument
2351   !> - variable definition now done in write_var_def
2352   !>
2353   !> @param[inout] td_file   file structure
2354   !> @param[inout] td_var    variable structure
2355   !> @param[in] id_start     index in the variable from which the data
2356   !> values will be read
2357   !> @param[in] id_count     number of indices selected along each dimension
2358   !-------------------------------------------------------------------
2359
2360      IMPLICIT NONE
2361
2362      ! Argument     
2363      TYPE(TFILE)              , INTENT(INOUT) :: td_file
2364      TYPE(TVAR)               , INTENT(INOUT) :: td_var
2365      INTEGER(i4), DIMENSION(:), INTENT(IN   ), OPTIONAL :: id_start
2366      INTEGER(i4), DIMENSION(:), INTENT(IN   ), OPTIONAL :: id_count
2367
2368      ! local variable
2369      !----------------------------------------------------------------
2370
2371      IF( ASSOCIATED(td_var%d_value) )THEN
2372         
2373         ! write variable value in file
2374         CALL iom_cdf__write_var_value(td_file, td_var, id_start, id_count)
2375      ENDIF
2376
2377   END SUBROUTINE iom_cdf__write_var
2378   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2379   SUBROUTINE iom_cdf__write_var_def(td_file, td_var,ld_chkdim)
2380   !-------------------------------------------------------------------
2381   !> @brief This subroutine define variable in the header of a netcdf
2382   !> file.
2383   !>
2384   !> @author J.Paul
2385   !> @date November, 2013 - Initial Version
2386   !> @date September, 2017
2387   !> - add option to not check dimension length
2388   !> @date August, 2017
2389   !> - extract from write_var
2390   !> @date November, 2019
2391   !> - compress 2D,3D, and 4D variables in netcdf4 files
2392   !> - set no_fill mode on every variables
2393   !>
2394   !> @param[in] td_file   file structure
2395   !> @param[in] td_var    variable structure
2396   !> @param[in] ld_chkdim check dimension length
2397   !> @return  variable id
2398   !-------------------------------------------------------------------
2399
2400      IMPLICIT NONE
2401
2402      ! Argument     
2403      TYPE(TFILE), INTENT(IN   ) :: td_file
2404      TYPE(TVAR),  INTENT(INOUT) :: td_var
2405      LOGICAL    , INTENT(IN   ), OPTIONAL :: ld_chkdim
2406
2407      ! local variable
2408      INTEGER(i4)                       :: il_status
2409      INTEGER(i4)                       :: il_ind
2410      INTEGER(i4)                       :: il_deflate
2411      INTEGER(i4), DIMENSION(ip_maxdim) :: il_dimid
2412
2413      LOGICAL                           :: ll_chg
2414      LOGICAL                           :: ll_defdim
2415
2416      ! loop indices
2417      INTEGER(i4) :: ji
2418      INTEGER(i4) :: jj
2419      !----------------------------------------------------------------
2420
2421      ll_defdim=file_check_var_dim(td_file, td_var, ld_chkdim )
2422      IF( ll_defdim )THEN
2423         ll_chg=.TRUE.
2424         DO ji=1,ip_maxdim
2425            IF( TRIM(fct_lower(cp_dimorder(ji:ji))) == &
2426            &   TRIM(fct_lower(td_var%c_name)) )THEN
2427               ll_chg=.FALSE.
2428               CALL logger_trace(TRIM(fct_lower(td_var%c_name))//' is var dimension')
2429               EXIT
2430            ENDIF
2431         ENDDO
2432         ! ugly patch until NEMO do not force to use 0. as FillValue
2433         IF( ll_chg )THEN
2434            ! not a dimension variable
2435            ! change FillValue
2436            SELECT CASE( TRIM(fct_lower(td_var%c_name)) )
2437               CASE DEFAULT
2438                  CALL var_chg_FillValue(td_var,0._dp)
2439               CASE('nav_lon','nav_lat', &
2440                  & 'glamt','glamu','glamv','glamf', &
2441                  & 'gphit','gphiu','gphiv','gphif', &
2442                  & 'e1t','e1u','e1v','e1f',         &
2443                  & 'e2t','e2u','e2v','e2f','ff',    &
2444                  & 'gcost','gcosu','gcosv','gcosf', &
2445                  & 'gsint','gsinu','gsinv','gsinf', &
2446                  & 'mbathy','misf','isf_draft','isfdraft',     &
2447                  & 'hbatt','hbatu','hbatv','hbatf', &
2448                  & 'gsigt','gsigu','gsigv','gsigf', &
2449                  & 'e3t_0','e3u_0','e3v_0','e3w_0', &
2450                  & 'e3f_0','gdepw_1d','gdept_1d',   &
2451                  & 'e3tp','e3wp','gdepw_0','rx1',   &
2452                  & 'gdept_0','gdepu','gdepv',       &
2453                  & 'hdept','hdepw','e3w_1d','e3t_1d',&
2454                  & 'tmask','umask','vmask','fmask'  )
2455                  ! do not change for coordinates and meshmask variables
2456               !CASE('no3','o2','po4','si', &
2457               !   & 'solub','ndepo','dust','fr_par', &
2458               !   & 'alk','dic','doc','fer' )
2459               !   ! do not change for BGC variables
2460            END SELECT
2461         ENDIF         
2462
2463         ! forced to use float type
2464         IF( td_var%d_unf /= 1. .AND. td_var%i_type==NF90_SHORT )THEN
2465            td_var%i_type=NF90_FLOAT
2466         ENDIF
2467
2468         IF( ALL( .NOT. td_var%t_dim(:)%l_use ) )THEN
2469            CALL logger_debug( &
2470            &  "IOM CDF WRITE VAR DEF scalar: define variable "//&
2471            &  TRIM(td_var%c_name)//" in file "//TRIM(td_file%c_name))
2472            ! scalar value
2473            il_status = NF90_DEF_VAR(td_file%i_id,          &
2474               &                     TRIM(td_var%c_name),   &
2475               &                     td_var%i_type,         &
2476               &                     varid=td_var%i_id) 
2477            CALL iom_cdf__check(il_status,"IOM CDF WRITE VAR DEF: ")
2478         ELSE
2479
2480            ! check which dimension use
2481            jj=0
2482            il_dimid(:)=0
2483            ! reorder dimension, so unused dimension won't be written
2484            DO ji = 1,  ip_maxdim
2485               IF( td_var%t_dim(ji)%l_use )THEN
2486                  jj=jj+1
2487                  il_dimid(jj)=dim_get_id(td_file%t_dim(:),td_var%t_dim(ji)%c_name)
2488               ENDIF
2489            ENDDO
2490
2491            CALL logger_debug( &
2492            &  "IOM CDF WRITE VAR DEF: define dimension to be used for variable "//&
2493            &  TRIM(td_var%c_name)//" in file "//TRIM(td_file%c_name))
2494
2495            DO ji=1,jj
2496               CALL logger_debug("IOM CDF WRITE VAR DEF: dimname : "//TRIM(td_var%t_dim(ji)%c_name))
2497               CALL logger_debug("IOM CDF WRITE VAR DEF: dimid "//TRIM(fct_str(il_dimid(ji))) )
2498            ENDDO
2499
2500            il_deflate=0
2501            ! compress 2D,3D, and 4D variables
2502            if( jj > 1 ) il_deflate=1
2503            CALL logger_debug("IOM CDF WRITE VAR DEF: deflate = "//TRIM(fct_str(il_deflate)))
2504            il_status = NF90_DEF_VAR(td_file%i_id,       &
2505               &                     TRIM(td_var%c_name),&
2506               &                     td_var%i_type,      &
2507               &                     il_dimid(1:jj),     &
2508               &                     varid=td_var%i_id,  &
2509               &                     deflate_level=il_deflate )
2510            CALL iom_cdf__check(il_status,"IOM CDF WRITE VAR DEF: ")
2511         ENDIF
2512         CALL logger_debug("IOM CDF WRITE VAR DEF: type = "//TRIM(fct_str(td_var%i_type)))
2513
2514         ! remove useless attribute
2515         il_ind=att_get_index( td_var%t_att(:), "ew_overlap" )
2516         IF( il_ind /= 0 )THEN
2517            IF( td_var%t_att(il_ind)%d_value(1) == -1 )THEN
2518               CALL var_del_att(td_var, td_var%t_att(il_ind))
2519            ENDIF
2520         ENDIF
2521
2522         DO ji = 1, td_var%i_natt
2523            CALL logger_debug( &
2524            &  " IOM CDF WRITE VAR DEF: put attribute "//TRIM(td_var%t_att(ji)%c_name)//&
2525            &  " for variable "//TRIM(td_var%c_name)//&
2526            &  " in file "//TRIM(td_file%c_name) )
2527
2528            ! forced FillValue to have same type than variable
2529            IF( TRIM(td_var%t_att(ji)%c_name) == '_FillValue' )THEN
2530
2531               td_var%t_att(ji)%i_type=td_var%i_type
2532
2533               SELECT CASE(td_var%t_att(ji)%i_type)
2534                  CASE(NF90_BYTE)
2535                     il_status = NF90_DEF_VAR_FILL(td_file%i_id,  &
2536                        &                          td_var%i_id,   &
2537                        &                          1_i4,          &
2538                        &                          INT(td_var%d_fill,i1))
2539                  CASE(NF90_SHORT)
2540                     il_status = NF90_DEF_VAR_FILL(td_file%i_id,  &
2541                        &                          td_var%i_id,   &
2542                        &                          1_i4,          &
2543                        &                          INT(td_var%d_fill,i2))
2544                  CASE(NF90_INT)
2545                     il_status = NF90_DEF_VAR_FILL(td_file%i_id,  &
2546                        &                          td_var%i_id,   &
2547                        &                          1_i4,          &
2548                        &                          INT(td_var%d_fill,i4))
2549                  CASE(NF90_FLOAT)
2550                     il_status = NF90_DEF_VAR_FILL(td_file%i_id,  &
2551                        &                          td_var%i_id,   &
2552                        &                          1_i4,          &
2553                        &                          REAL(td_var%d_fill,sp))
2554                  CASE(NF90_DOUBLE)
2555                     il_status = NF90_DEF_VAR_FILL(td_file%i_id,  &
2556                        &                          td_var%i_id,   &
2557                        &                          1_i4,          &
2558                        &                          REAL(td_var%d_fill,dp))
2559               END SELECT
2560               CALL iom_cdf__check(il_status,"IOM CDF WRITE VAR DEF FILL: ")
2561
2562            ELSE
2563
2564               SELECT CASE(td_var%t_att(ji)%i_type)
2565                  CASE(NF90_CHAR)
2566                     IF( TRIM(td_var%t_att(ji)%c_value) /= '' )THEN
2567                        il_status = NF90_PUT_ATT(td_file%i_id, td_var%i_id,      &
2568                        &                        TRIM(td_var%t_att(ji)%c_name),  &
2569                        &                        TRIM(td_var%t_att(ji)%c_value)  )
2570                     ENDIF
2571                  CASE(NF90_BYTE)
2572                     il_status = NF90_PUT_ATT(td_file%i_id,                   &
2573                     &                        td_var%i_id,                    &
2574                     &                        TRIM(td_var%t_att(ji)%c_name),  &
2575                     &                        INT(td_var%t_att(ji)%d_value(:),i1))
2576                  CASE(NF90_SHORT)
2577                     il_status = NF90_PUT_ATT(td_file%i_id,                   &
2578                     &                        td_var%i_id,                    &
2579                     &                        TRIM(td_var%t_att(ji)%c_name),  &
2580                     &                        INT(td_var%t_att(ji)%d_value(:),i2))
2581                  CASE(NF90_INT)
2582                     il_status = NF90_PUT_ATT(td_file%i_id,                   &
2583                     &                        td_var%i_id,                    &
2584                     &                        TRIM(td_var%t_att(ji)%c_name),  &
2585                     &                        INT(td_var%t_att(ji)%d_value(:),i4))
2586                  CASE(NF90_FLOAT)
2587                     il_status = NF90_PUT_ATT(td_file%i_id,                   &
2588                     &                        td_var%i_id,                    &
2589                     &                        TRIM(td_var%t_att(ji)%c_name),  &
2590                     &                        REAL(td_var%t_att(ji)%d_value(:),sp))
2591                  CASE(NF90_DOUBLE)
2592                     il_status = NF90_PUT_ATT(td_file%i_id,                   &
2593                     &                        td_var%i_id,                    &
2594                     &                        TRIM(td_var%t_att(ji)%c_name),  &
2595                     &                        REAL(td_var%t_att(ji)%d_value(:),dp))
2596               END SELECT
2597               CALL iom_cdf__check(il_status,"IOM CDF WRITE VAR DEF: ")
2598
2599            ENDIF
2600
2601         ENDDO
2602
2603      ENDIF
2604      CALL logger_debug("IOM CDF WRITE VAR DEF: type = "//TRIM(fct_str(td_var%i_type)))
2605
2606   END SUBROUTINE iom_cdf__write_var_def
2607   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2608   SUBROUTINE iom_cdf__write_var_value(td_file, td_var, id_start, id_count)
2609   !-------------------------------------------------------------------
2610   !> @brief This subroutine put variable value in an opened netcdf file.
2611   !>
2612   !> @details
2613   !> The variable is written in the type define in variable structure.
2614   !> Only dimension used are printed, and fillValue in array are
2615   !> replaced by default fill values defined in module netcdf for each type.
2616   !>
2617   !> @author J.Paul
2618   !> @date November, 2013 - Initial Version
2619   !> @date June, 2015
2620   !> - reuse scale factor and offset, before writing variable
2621   !> @August, 2017
2622   !> - use start and count array to write variable value
2623   !>
2624   !> @param[in] td_file   file structure
2625   !> @param[in] td_var    variable structure
2626   !> @param[in] id_start  index in the variable from which the data
2627   !> values will be read
2628   !> @param[in] id_count  number of indices selected along each dimension
2629   !-------------------------------------------------------------------
2630
2631      IMPLICIT NONE
2632
2633      ! Argument     
2634      TYPE(TFILE)                   , INTENT(IN   ) :: td_file
2635      TYPE(TVAR)                    , INTENT(IN   ) :: td_var
2636      INTEGER(i4)     , DIMENSION(:), INTENT(IN   ), OPTIONAL :: id_start
2637      INTEGER(i4)     , DIMENSION(:), INTENT(IN   ), OPTIONAL :: id_count
2638
2639      ! local variable
2640      INTEGER(i4)                       :: il_status
2641      INTEGER(i4), DIMENSION(ip_maxdim) :: il_order
2642      INTEGER(i4), DIMENSION(ip_maxdim) :: il_shape
2643      REAL(dp),    DIMENSION(:,:,:,:), ALLOCATABLE :: dl_value
2644
2645      INTEGER(i4), DIMENSION(ip_maxdim) :: il_start
2646      INTEGER(i4), DIMENSION(ip_maxdim) :: il_count     
2647
2648      INTEGER(i4), DIMENSION(ip_maxdim) :: il_start_ord
2649      INTEGER(i4), DIMENSION(ip_maxdim) :: il_count_ord     
2650      ! loop indices
2651      INTEGER(i4) :: ji, jj
2652      !----------------------------------------------------------------
2653
2654      ! check which dimension use
2655      CALL logger_debug( &
2656      &  "IOM CDF WRITE VAR VALUE: get dimension to be used for variable "//&
2657      &  TRIM(td_var%c_name)//" in file "//TRIM(td_file%c_name))
2658 
2659      il_start(:)=1
2660      IF( PRESENT(id_start) ) il_start(:)=id_start(:)
2661
2662      il_count(:)=td_var%t_dim(:)%i_len
2663      IF( PRESENT(id_count) ) il_count(:)=id_count(:)
2664 
2665      ! use scale factor and offset
2666      WHERE( td_var%d_value(:,:,:,:) /= td_var%d_fill )
2667         td_var%d_value(:,:,:,:) = &
2668         &  (td_var%d_value(:,:,:,:)-td_var%d_ofs)/td_var%d_scf
2669      END WHERE
2670
2671      jj=0
2672      DO ji = 1, ip_maxdim
2673         IF( td_var%t_dim(ji)%l_use )THEN
2674            jj=jj+1
2675            il_order(jj)=ji
2676            il_shape(jj)=td_var%t_dim(ji)%i_len
2677         ENDIF
2678      ENDDO
2679      ! dimension not use
2680      DO ji = 1, ip_maxdim
2681         IF( .NOT. td_var%t_dim(ji)%l_use )THEN
2682            jj=jj+1
2683            il_order(jj)=ji
2684            il_shape(jj)=td_var%t_dim(ji)%i_len
2685         ENDIF
2686      ENDDO
2687
2688      ALLOCATE( dl_value( il_shape(1),il_shape(2),il_shape(3),il_shape(4) ) )
2689
2690      ! reshape array, so useless dimension won't be written
2691      dl_value(:,:,:,:)=RESHAPE(source=td_var%d_value(:,:,:,:),&
2692      &                         SHAPE = il_shape(:), &
2693      &                         ORDER = il_order(:))
2694
2695      DO ji = 1, ip_maxdim
2696         il_start_ord(il_order(ji))=il_start(ji)
2697         il_count_ord(il_order(ji))=il_count(ji)
2698      ENDDO
2699
2700      ! put value
2701      CALL logger_debug( &
2702      &  "IOM CDF WRITE VAR VALUE: put "//TRIM(td_var%c_name)//" value "//&
2703      &  "in file "//TRIM(td_file%c_name))
2704      il_status = NF90_PUT_VAR( td_file%i_id, td_var%i_id, dl_value(:,:,:,:), &
2705      &                      start=il_start_ord(:), &
2706      &                      count=il_count_ord(:) )
2707      CALL iom_cdf__check(il_status,"IOM CDF WRITE VAR VALUE ("//&
2708      &  TRIM(td_var%c_name)//") :" )
2709
2710      DEALLOCATE( dl_value )
2711
2712   END SUBROUTINE iom_cdf__write_var_value
2713   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2714END MODULE iom_cdf
Note: See TracBrowser for help on using the repository browser.