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 @ 14712

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

update: cf changelog inside documentation

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