SURFEX v8.1
General documentation of Surfex
mode_read_cdf.F90
Go to the documentation of this file.
1 !SFX_LIC Copyright 1994-2014 CNRS, Meteo-France and Universite Paul Sabatier
2 !SFX_LIC This is part of the SURFEX software governed by the CeCILL-C licence
3 !SFX_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
4 !SFX_LIC for details. version 1.
6 !===================================================================
7 !
8 !
9 USE modi_abor1_sfx
10 !
11 USE yomhook ,ONLY : lhook, dr_hook
12 USE parkind1 ,ONLY : jprb
13 !
14 CONTAINS
15 !-------------------------------------------------------------------
16 !-------------------------------------------------------------------
17 ! ####################
18  SUBROUTINE handle_err_cdf(status,line)
19 ! ####################
20 USE netcdf
21 !
22 IMPLICIT NONE
23 INTEGER, INTENT(IN) :: status
24  CHARACTER(*), INTENT(IN) :: line
25 REAL(KIND=JPRB) :: ZHOOK_HANDLE
26 !
27 !
28 IF (lhook) CALL dr_hook('MODE_READ_CDF:HANDLE_ERR_CDF',0,zhook_handle)
29 IF (status /= nf90_noerr) THEN
30  CALL abor1_sfx('MODE_READ_NETCDF: HANDLE_ERR_CDF:'//trim(line))
31 END IF
32 IF (lhook) CALL dr_hook('MODE_READ_CDF:HANDLE_ERR_CDF',1,zhook_handle)
33 END SUBROUTINE handle_err_cdf
34 !-------------------------------------------------------------------
35 !-------------------------------------------------------------------
36 ! ####################
37  SUBROUTINE get1dcdf(KCDF_ID,IDVAR,PMISSVALUE,PVALU1D)
38 ! ####################
39 !
40 USE netcdf
41 !
42 IMPLICIT NONE
43 !
44 INTEGER,INTENT(IN) :: KCDF_ID !netcdf file identifiant
45 INTEGER,INTENT(IN) :: IDVAR !variable to read identifiant
46 REAL, INTENT(OUT) :: PMISSVALUE !undefined value
47 REAL,DIMENSION(:),INTENT(OUT) :: PVALU1D !value array
48 !
49 integer :: status
50 character(len=80) :: HACTION
51 integer,save :: NDIMS=1
52 integer :: KVARTYPE
53 integer,DIMENSION(:),ALLOCATABLE :: NVARDIMID,NVARDIMLEN
54 character(len=80),DIMENSION(:),ALLOCATABLE :: NVARDIMNAM
55 integer :: JLOOP
56 integer :: NGATTS
57 character(len=80),DIMENSION(:),ALLOCATABLE :: HNAME
58 REAL,DIMENSION(:),ALLOCATABLE :: ZVALU1D !value array
59 REAL(KIND=JPRB) :: ZHOOK_HANDLE
60 !
61 !
62 IF (lhook) CALL dr_hook('MODE_READ_CDF:GET1DCDF',0,zhook_handle)
63 pmissvalue=-9999.9
64 ALLOCATE(nvardimid(ndims))
65 ALLOCATE(nvardimlen(ndims))
66 ALLOCATE(nvardimnam(ndims))
67 nvardimid(:)=0
68 nvardimlen(:)=0
69 nvardimnam(:)=' '
70 !
71 haction='get variable type'
72 status=nf90_inquire_variable(kcdf_id,idvar,xtype=kvartype)
73 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
74 !write(0,*) 'variable type = ',KVARTYPE
75 !
76 haction='get variable dimensions identifiant'
77 status=nf90_inquire_variable(kcdf_id,idvar,dimids=nvardimid)
78 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
79 !write(0,*) 'variable dimension ',NDIMS,' identifiant ',NVARDIMID(NDIMS)
80 !
81 haction='get variable dimensions name'
82 status=nf90_inquire_dimension(kcdf_id,nvardimid(ndims),name=nvardimnam(ndims))
83 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
84 !
85 haction='get variable dimensions length'
86 status=nf90_inquire_dimension(kcdf_id,nvardimid(ndims),len=nvardimlen(ndims))
87 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
88 !write(0,*) 'variable dimension ',NDIMS,' named ',NVARDIMNAM(NDIMS),&
89 ! &'has a length of',NVARDIMLEN(NDIMS)
90 !
91 haction='get attributs'
92 !status=nf90_inq_natts(KCDF_ID,NGATTS)
93 status=nf90_inquire_variable(kcdf_id,idvar,natts=ngatts)
94 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
95 !write(0,*) 'number of attributes = ',NGATTS
96 allocate(hname(1:ngatts))
97 !
98 DO jloop=1,ngatts
99  status=nf90_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
100  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
101  !write(0,*) 'attributes names = ', hname(JLOOP)
102  if (trim(hname(jloop))=='missing_value') then
103  !write(0,*) 'missing value search '
104  haction='get missing value'
105  status=nf90_get_att(kcdf_id,idvar,"missing_value",pmissvalue)
106  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
107  !write(0,*) 'missing value = ',PMISSVALUE
108  endif
109 ENDDO
110 !
111 ALLOCATE(zvalu1d(1:nvardimlen(ndims)))
112 zvalu1d=0.
113 !
114 IF (kvartype>=5) then
115  haction='get variable values (1D)'
116  status=nf90_get_var(kcdf_id,idvar,zvalu1d(:))
117  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
118 ENDIF
119 pvalu1d(:)=zvalu1d(:)
120 IF (ALLOCATED(zvalu1d )) DEALLOCATE(zvalu1d)
121 IF (lhook) CALL dr_hook('MODE_READ_CDF:GET1DCDF',1,zhook_handle)
122 !
123 END SUBROUTINE get1dcdf
124 !-------------------------------------------------------------------
125 !-------------------------------------------------------------------
126 ! ####################
127  SUBROUTINE get2dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
128  PMISSVALUE,PVALU2D)
129 ! ####################
130 !
131 USE netcdf
132 !
133 IMPLICIT NONE
134 !
135 INTEGER,INTENT(IN) :: KCDF_ID !netcdf file identifiant
136 INTEGER,INTENT(IN) :: IDVAR !variable to read identifiant
137 REAL,DIMENSION(:),INTENT(OUT) :: PDIM1,PDIM2 !dimensions for PVALU2D array
138  CHARACTER(len=80),INTENT(OUT) :: HDIM1NAME,HDIM2NAME !dimensions names
139 REAL, INTENT(OUT) :: PMISSVALUE
140 REAL,DIMENSION(:,:),INTENT(OUT) :: PVALU2D !value array
141 !
142 integer :: status
143 character(len=80) :: HACTION
144 integer,save :: NDIMS=2
145 integer :: KVARTYPE
146 integer,DIMENSION(:),ALLOCATABLE :: NVARDIMID,NVARDIMLEN
147 character(len=80),DIMENSION(:),ALLOCATABLE :: NVARDIMNAM
148 integer :: JLOOP
149 integer :: NGATTS
150 character(len=80),DIMENSION(:),ALLOCATABLE :: HNAME
151 real :: ZMISS1,ZMISS2
152 REAL,DIMENSION(:,:),ALLOCATABLE :: ZVALU2D !value array
153 REAL(KIND=JPRB) :: ZHOOK_HANDLE
154 !
155 !
156 IF (lhook) CALL dr_hook('MODE_READ_CDF:GET2DCDF',0,zhook_handle)
157 pmissvalue=-9999.9
158 ALLOCATE(nvardimid(ndims))
159 ALLOCATE(nvardimlen(ndims))
160 ALLOCATE(nvardimnam(ndims))
161 nvardimid(:)=0
162 nvardimlen(:)=0
163 nvardimnam(:)=' '
164 !
165 haction='get variable type'
166 status=nf90_inquire_variable(kcdf_id,idvar,xtype=kvartype)
167 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
168 !write(0,*) 'variable type = ',KVARTYPE
169 !
170 haction='get variable dimensions identifiant'
171 status=nf90_inquire_variable(kcdf_id,idvar,dimids=nvardimid)
172 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
173 !
174 haction='get attributs'
175 !status=nf90_inq_natts(KCDF_ID,NGATTS)
176 status=nf90_inquire_variable(kcdf_id,idvar,natts=ngatts)
177 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
178 !write(0,*) 'number of attributes = ',NGATTS
179 allocate(hname(1:ngatts))
180 !
181 DO jloop=1,ngatts
182  status=nf90_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
183  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
184  !write(0,*) 'attributes names = ', hname(JLOOP)
185  if (trim(hname(jloop))=='missing_value') then
186  !write(0,*) 'missing value search '
187  haction='get missing value'
188  status=nf90_get_att(kcdf_id,idvar,"missing_value",pmissvalue)
189  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
190  !write(0,*) 'missing value = ',PMISSVALUE
191  endif
192 ENDDO
193 !
194 !
195 DO jloop=1,ndims
196  haction='get variable dimensions name'
197  status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop),name=nvardimnam(jloop))
198  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
199  haction='get variable dimensions length'
200  status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop),len=nvardimlen(jloop))
201  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
202  !write(0,*) 'variable dimension ',JLOOP,' named ',NVARDIMNAM(JLOOP),&
203  ! &'has a length of',NVARDIMLEN(JLOOP)
204 ENDDO
205 !
206 ALLOCATE(zvalu2d(1:nvardimlen(1),1:nvardimlen(2)))
207 zvalu2d=0.
208 IF (kvartype>=5) then
209  haction='get variable values (2D)'
210  status=nf90_get_var(kcdf_id,idvar,zvalu2d(:,:))
211  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
212 ENDIF
213 pvalu2d(:,:)=zvalu2d(:,:)
214 !
215  CALL get1dcdf(kcdf_id,nvardimid(1),zmiss1,pdim1)
216  CALL get1dcdf(kcdf_id,nvardimid(2),zmiss2,pdim2)
217 hdim1name=nvardimnam(1)
218 hdim2name=nvardimnam(2)
219 IF (ALLOCATED(zvalu2d )) DEALLOCATE(zvalu2d)
220 IF (lhook) CALL dr_hook('MODE_READ_CDF:GET2DCDF',1,zhook_handle)
221 !
222 END SUBROUTINE get2dcdf
223 !--------------------------------------------------------------------
224 !-------------------------------------------------------------------
225 ! ####################
226  SUBROUTINE read_latlonval_cdf(HFILENAME,HNCVARNAME,PLON,PLAT,PVAL)
227 ! ####################
228 !
229 USE netcdf
230 !
231 IMPLICIT NONE
232 !
233  CHARACTER(LEN=28), INTENT(IN) :: HFILENAME ! Name of the field file.
234  CHARACTER(LEN=28), INTENT(IN) :: HNCVARNAME ! Name of variable to read in netcdf file
235 REAL, DIMENSION(:), INTENT(OUT) :: PLON,PLAT ! Longitudes/latitudes innetcdf file
236 REAL, DIMENSION(:), INTENT(OUT) :: PVAL ! value to get
237 !
238 integer :: status
239 integer :: kcdf_id
240 integer :: NBVARS
241 character(len=80) :: HACTION
242 character(len=80),DIMENSION(:),ALLOCATABLE :: VARNAME
243 integer ::JLOOP1,JDIM1,JDIM2,JLOOP
244 integer ::ID_VARTOGET,ID_VARTOGET1,ID_VARTOGET2
245 integer ::NVARDIMS
246 integer ::NLEN
247 integer,dimension(1) :: IDIMID
248 integer,DIMENSION(1:2) :: NLEN2D,IDIMID2D
249 integer,DIMENSION(:),ALLOCATABLE :: NVARDIMID,NVARDIMLEN
250 character(len=80),DIMENSION(:),ALLOCATABLE :: NVARDIMNAM
251 real,DIMENSION(:),ALLOCATABLE :: ZVALU
252 real,DIMENSION(:,:),ALLOCATABLE :: ZVALU2D
253 real :: ZMISS
254 real,DIMENSION(:),ALLOCATABLE :: ZDIM1
255 real,DIMENSION(:),ALLOCATABLE :: ZDIM2
256 character(len=80) :: YDIM1NAME,YDIM2NAME
257 integer :: ILONFOUND,ILATFOUND, IARG
258 REAL(KIND=JPRB) :: ZHOOK_HANDLE
259 !
260 !
261 !* 1. Open the netcdf file
262 ! --------------------
263 IF (lhook) CALL dr_hook('MODE_READ_CDF:READ_LATLONVAL_CDF',0,zhook_handle)
264 status=-9999
265 kcdf_id=-9999
266 haction='open netcdf'
267 status=nf90_open(hfilename,nf90_nowrite,kcdf_id)
268 !write(0,*) 'status=',status
269 !write(0,*) 'identifiant de ',HFILENAME,'=',kcdf_id
270 if (status/=nf90_noerr) then
271  CALL handle_err_cdf(status,haction)
272 !else
273 ! write(0,*) 'netcdf file opened: ',HFILENAME
274 endif
275 !
276 !-----------
277 !
278 !* 2. get the number of variables in netcdf file
279 ! ------------------------------------------
280 haction='get number of variables'
281 status=nf90_inquire(kcdf_id,nvariables=nbvars)
282 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
283 !write(0,*) 'nb vars', NBVARS
284 ALLOCATE(varname(nbvars))
285 !
286 !-----------
287 !
288 !* 3. get the variables names in netcdf file
289 ! --------------------------------------
290 id_vartoget1=0
291 id_vartoget2=0
292 DO jloop1=1,nbvars
293  haction='get variables names'
294  status=nf90_inquire_variable(kcdf_id,jloop1,name=varname(jloop1))
295  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
296  !write(0,*) 'var',JLOOP1,' name: ',VARNAME(JLOOP1)
297  if (varname(jloop1)==hncvarname) then
298  !write(0,*) 'var',JLOOP1,' corresponding to variable required'
299  id_vartoget1=jloop1
300  endif
301  if (varname(jloop1)/=hncvarname) then
302  if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
303  (scan(trim(varname(jloop1)),trim(hncvarname))==1)) then
304  !write(0,*) 'var',JLOOP1,VARNAME(JLOOP1),' could correspond to variable required ?'
305  !write(0,*) HNCVARNAME,' is variable required; only ',VARNAME(JLOOP1),' found'
306  id_vartoget2=jloop1
307  endif
308  endif
309 ENDDO
310 if (id_vartoget1/=0) then
311  id_vartoget=id_vartoget1
312 else
313  id_vartoget=id_vartoget2
314 endif
315 if (id_vartoget==0) then
316  haction='close netcdf'
317  status=nf90_close(kcdf_id)
318  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
319  CALL abor1_sfx('MODE_READ_NETCDF: READ_LATLONVAL_CDF')
320 endif
321 !-----------
322 !
323 !* 4. get the variable in netcdf file
324 ! -------------------------------
325 !
326 ! 4.1 get the variable dimensions number
327 ! -----------------------------------
328 !
329 haction='get variable dimensions number'
330 status=nf90_inquire_variable(kcdf_id,id_vartoget,ndims=nvardims)
331 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
332 !
333 ! 4.2 get the variable dimensions length and values
334 ! ----------------------------------------------
335 SELECT CASE (nvardims)
336 !CAS 1D
337  CASE (1)
338  !write(0,*) 'variable dimensions number = ',NVARDIMS
339  haction='get variable dimensions length'
340  status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=idimid)
341  status=nf90_inquire_dimension(kcdf_id,idimid(1),len=nlen)
342  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
343  ALLOCATE(zvalu(nlen))
344  !write(0,*) 'call GET1DCDF'
345  CALL get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
346  pval(:)=zvalu(:)
347 !CAS 2D
348  CASE (2)
349  !write(0,*) 'variable dimensions number = ',NVARDIMS
350  status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=idimid2d)
351  DO jloop=1,nvardims
352  haction='get variable dimensions length'
353  status=nf90_inquire_dimension(kcdf_id,idimid2d(jloop),len=nlen2d(jloop))
354  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
355  ENDDO
356  ALLOCATE(zvalu2d(nlen2d(1),nlen2d(2)))
357  ALLOCATE(zdim1(nlen2d(1)))
358  ALLOCATE(zdim2(nlen2d(2)))
359  !write(0,*) 'call GET2DCDF'
360  CALL get2dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,ydim2name,&
361  zmiss,zvalu2d)
362  !write(0,*) 'YDIM1NAME: ',YDIM1NAME
363  !write(0,*) 'YDIM2NAME: ',YDIM2NAME
364  if ((ydim1name=='lon').OR.(ydim1name=='longitude')) ilonfound=1
365  if ((ydim2name=='lon').OR.(ydim2name=='longitude')) ilonfound=2
366  if ((ydim1name=='lat').OR.(ydim1name=='latitude')) ilatfound=1
367  if ((ydim2name=='lat').OR.(ydim2name=='latitude')) ilatfound=2
368  iarg=0
369 !
370 ! 4.3 complete arrays
371 ! ---------------
372  IF ((ilonfound==1).AND.(ilatfound==2)) then
373  !write(0,*) 'ILONFOUND',ILONFOUND,'ILATFOUND',ILATFOUND
374  DO jdim1=1,SIZE(zdim1)
375  DO jdim2=1,SIZE(zdim2)
376  iarg=iarg+1
377  pval(iarg)=zvalu2d(jdim1,jdim2)
378  plon(iarg)=zdim1(jdim1)
379  plat(iarg)=zdim2(jdim2)
380  ENDDO
381  ENDDO
382  ELSEIF ((ilonfound==2).AND.(ilatfound==1)) then
383  !write(0,*) 'ILONFOUND',ILONFOUND,'ILATFOUND',ILATFOUND
384  DO jdim1=1,SIZE(zdim1)
385  DO jdim2=1,SIZE(zdim2)
386  iarg=iarg+1
387  pval(iarg)=zvalu2d(jdim1,jdim2)
388  plat(iarg)=zdim1(jdim1)
389  plon(iarg)=zdim2(jdim2)
390  ENDDO
391  ENDDO
392  ELSE
393  write(0,*) '*****WARNING*****: incompatible dimensions to lat/lon/value arrays'
394  ENDIF
395 !
396 END SELECT
397 !
398 !
399 !-----------
400 !* 10. Close the netcdf file
401 ! ---------------------
402 haction='close netcdf'
403 status=nf90_close(kcdf_id)
404 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
405 !write(0,*) 'OK: netcdf file closed: ',HFILENAME
406 !
407 !-----------
408 !* 11. Deallocate
409 ! ----------
410 IF (ALLOCATED(varname )) DEALLOCATE(varname)
411 IF (ALLOCATED(zvalu )) DEALLOCATE(zvalu )
412 IF (ALLOCATED(zvalu2d )) DEALLOCATE(zvalu2d)
413 IF (ALLOCATED(zdim1 )) DEALLOCATE(zdim1 )
414 IF (ALLOCATED(zdim2 )) DEALLOCATE(zdim2 )
415 !
416 IF (ALLOCATED(nvardimid )) DEALLOCATE(nvardimid )
417 IF (ALLOCATED(nvardimnam )) DEALLOCATE(nvardimnam)
418 IF (ALLOCATED(nvardimlen )) DEALLOCATE(nvardimlen)
419 IF (lhook) CALL dr_hook('MODE_READ_CDF:READ_LATLONVAL_CDF',1,zhook_handle)
420 END SUBROUTINE read_latlonval_cdf
421 !------------------------------------------------------------------------------
422 !==============================================================================
423 ! ####################
424  SUBROUTINE read_dim_cdf(HFILENAME,HNCVARNAME,KDIM)
425 ! ####################
426 !
427 USE netcdf
428 !
429 IMPLICIT NONE
430 !
431  CHARACTER(LEN=28), INTENT(IN) :: HFILENAME ! Name of the field file.
432  CHARACTER(LEN=28), INTENT(IN) :: HNCVARNAME ! Name of variable to read in netcdf file
433 INTEGER, INTENT(OUT):: KDIM ! value of dimension to get
434 !
435 integer :: status
436 integer :: kcdf_id
437 integer :: NBVARS
438 character(len=80) :: HACTION
439 character(len=80),DIMENSION(:),ALLOCATABLE :: VARNAME
440 integer ::JLOOP1,JLOOP
441 integer ::ID_VARTOGET,ID_VARTOGET1,ID_VARTOGET2
442 integer ::NVARDIMS
443 integer, dimension(1) :: NDIMID
444 integer,DIMENSION(2) ::NLEN2D, NDIMID2D
445 REAL(KIND=JPRB) :: ZHOOK_HANDLE
446 !
447 !
448 !* 1. Open the netcdf file
449 ! --------------------
450 IF (lhook) CALL dr_hook('MODE_READ_CDF:READ_DIM_CDF',0,zhook_handle)
451 haction='open netcdf'
452 status=nf90_open(hfilename,nf90_nowrite,kcdf_id)
453 if (status/=nf90_noerr) then
454  CALL handle_err_cdf(status,haction)
455 !else
456 endif
457 !
458 !-----------
459 !
460 !* 2. get the number of variables in netcdf file
461 ! ------------------------------------------
462 haction='get number of variables'
463 status=nf90_inquire(kcdf_id,nvariables=nbvars)
464 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
465 !write(0,*) 'nb vars', NBVARS
466 ALLOCATE(varname(nbvars))
467 !
468 !-----------
469 !
470 !* 3. get the variables names in netcdf file
471 ! --------------------------------------
472 id_vartoget1=0
473 id_vartoget2=0
474 DO jloop1=1,nbvars
475  haction='get variables names'
476  status=nf90_inquire_variable(kcdf_id,jloop1,name=varname(jloop1))
477  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
478  !write(0,*) 'var',JLOOP1,' name: ',VARNAME(JLOOP1)
479  if (varname(jloop1)==hncvarname) then
480  !write(0,*) 'var',JLOOP1,' corresponding to variable required'
481  id_vartoget1=jloop1
482  endif
483  if (varname(jloop1)/=hncvarname) then
484  if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
485  (scan(trim(varname(jloop1)),trim(hncvarname))==1)) then
486  !write(0,*) 'var',JLOOP1,VARNAME(JLOOP1),' could correspond to variable required ?'
487  !write(0,*) HNCVARNAME,' is variable required; only ',VARNAME(JLOOP1),' found'
488  id_vartoget2=jloop1
489  endif
490  endif
491 ENDDO
492 if (id_vartoget1/=0) then
493  id_vartoget=id_vartoget1
494 else
495  id_vartoget=id_vartoget2
496 endif
497 if (id_vartoget==0) then
498  haction='close netcdf'
499  status=nf90_close(kcdf_id)
500  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
501  CALL abor1_sfx('MODE_READ_CDF: READ_DIM_CDF')
502 endif
503 !-----------
504 !
505 !* 4. get the total dimension of HNCVARNAME
506 ! -------------------------------------
507 !
508 ! 4.1 get the variable dimensions number
509 ! -----------------------------------
510 !
511 haction='get variable dimensions number'
512 status=nf90_inquire_variable(kcdf_id,id_vartoget,ndims=nvardims)
513 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
514 !write(0,*) 'variable dimensions number = ',NVARDIMS
515 !
516 ! 4.2 get the variable dimensions length
517 ! ----------------------------------
518 SELECT CASE (nvardims)
519 !CAS 1D
520  CASE (1)
521  haction='get variable dimensions length'
522  status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=ndimid)
523  status=nf90_inquire_dimension(kcdf_id,ndimid(1),len=kdim)
524  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
525 !
526 !CAS 2D
527  CASE (2)
528  kdim=1
529  status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=ndimid2d)
530  DO jloop=1,nvardims
531  haction='get variable dimensions length'
532  status=nf90_inquire_dimension(kcdf_id,ndimid2d(jloop),len=nlen2d(jloop))
533  if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
534  kdim=kdim*nlen2d(jloop)
535  ENDDO
536 END SELECT
537 !-----------
538 !* 10. Close the netcdf file
539 ! ---------------------
540 haction='close netcdf'
541 status=nf90_close(kcdf_id)
542 if (status/=nf90_noerr) CALL handle_err_cdf(status,haction)
543 !write(0,*) 'OK: netcdf file closed: ',HFILENAME
544 !
545 !-----------
546 !* 11. Deallocate
547 ! ----------
548 IF (ALLOCATED(varname )) DEALLOCATE(varname)
549 IF (lhook) CALL dr_hook('MODE_READ_CDF:READ_DIM_CDF',1,zhook_handle)
550 !
551 END SUBROUTINE read_dim_cdf
552 !------------------------------------------------------------------------------
553 !==============================================================================
554 END MODULE mode_read_cdf
static const char * trim(const char *name, int *n)
Definition: drhook.c:2383
subroutine handle_err_cdf(status, line)
quick &counting sorts only inumt inumt name
subroutine read_latlonval_cdf(HFILENAME, HNCVARNAME, PLON, PLAT, PVAL)
subroutine abor1_sfx(YTEXT)
Definition: abor1_sfx.F90:7
integer, parameter jprb
Definition: parkind1.F90:32
subroutine get2dcdf(KCDF_ID, IDVAR, PDIM1, HDIM1NAME, PDIM2, HDIM2NAME, PMISSVALUE, PVALU2D)
logical lhook
Definition: yomhook.F90:15
subroutine get1dcdf(KCDF_ID, IDVAR, PMISSVALUE, PVALU1D)
subroutine read_dim_cdf(HFILENAME, HNCVARNAME, KDIM)