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