11 USE yomhook
,ONLY : lhook, dr_hook
12 USE parkind1
,ONLY : jprb
21 INTEGER,
INTENT(IN) :: status
22 CHARACTER(*),
INTENT(IN) :: line
23 REAL(KIND=JPRB) :: zhook_handle
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))
31 IF (lhook) CALL dr_hook(
'MODE_READ_CDF:HANDLE_ERR_CDF',1,zhook_handle)
36 SUBROUTINE get1dcdf(KCDF_ID,IDVAR,PMISSVALUE,PVALU1D)
41 INTEGER,
INTENT(IN) :: kcdf_id
42 INTEGER,
INTENT(IN) :: idvar
43 REAL,
INTENT(OUT) :: pmissvalue
44 REAL,
DIMENSION(:),
INTENT(OUT) :: pvalu1d
47 character(len=80) :: haction
48 integer,
save :: ndims=1
50 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid,nvardimlen
51 character(len=80),
DIMENSION(:),
ALLOCATABLE :: nvardimnam
54 character(len=80),
DIMENSION(:),
ALLOCATABLE :: hname
55 REAL,
DIMENSION(:),
ALLOCATABLE :: zvalu1d
56 REAL(KIND=JPRB) :: zhook_handle
60 IF (lhook) CALL dr_hook(
'MODE_READ_CDF:GET1DCDF',0,zhook_handle)
62 ALLOCATE(nvardimid(ndims))
63 ALLOCATE(nvardimlen(ndims))
64 ALLOCATE(nvardimnam(ndims))
69 haction=
'get variable type'
70 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
74 haction=
'get variable dimensions identifiant'
75 status=nf_inq_vardimid(kcdf_id,idvar,nvardimid(ndims))
79 haction=
'get variable dimensions name'
80 status=nf_inq_dimname(kcdf_id,nvardimid(ndims),nvardimnam(ndims))
83 haction=
'get variable dimensions length'
84 status=nf_inq_dimlen(kcdf_id,nvardimid(ndims),nvardimlen(ndims))
89 haction=
'get attributs'
91 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
94 allocate(hname(1:ngatts))
97 status=nf_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
100 if (trim(hname(jloop))==
'missing_value')
then
102 haction=
'get missing value'
103 status=nf_get_att_double(kcdf_id,idvar,
"missing_value",pmissvalue)
109 ALLOCATE(zvalu1d(1:nvardimlen(ndims)))
112 IF (kvartype>=5)
then
113 haction=
'get variable values (1D)'
114 status=nf_get_var_double(kcdf_id,idvar,zvalu1d(:))
117 pvalu1d(:)=zvalu1d(:)
118 IF (
ALLOCATED(zvalu1d ))
DEALLOCATE(zvalu1d)
119 IF (lhook) CALL dr_hook(
'MODE_READ_CDF:GET1DCDF',1,zhook_handle)
125 SUBROUTINE get2dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
131 INTEGER,
INTENT(IN) :: kcdf_id
132 INTEGER,
INTENT(IN) :: idvar
133 REAL,
DIMENSION(:),
INTENT(OUT) :: pdim1,pdim2
134 CHARACTER(len=80),
INTENT(OUT) :: hdim1name,hdim2name
135 REAL,
INTENT(OUT) :: pmissvalue
136 REAL,
DIMENSION(:,:),
INTENT(OUT) :: pvalu2d
139 character(len=80) :: haction
140 integer,
save :: ndims=2
142 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid,nvardimlen
143 character(len=80),
DIMENSION(:),
ALLOCATABLE :: nvardimnam
144 integer :: jloop2, jloop
146 character(len=80),
DIMENSION(:),
ALLOCATABLE :: hname
147 real :: zmiss1,zmiss2
148 REAL,
DIMENSION(:,:),
ALLOCATABLE :: zvalu2d
149 REAL(KIND=JPRB) :: zhook_handle
153 IF (lhook) CALL dr_hook(
'MODE_READ_CDF:GET2DCDF',0,zhook_handle)
155 ALLOCATE(nvardimid(ndims))
156 ALLOCATE(nvardimlen(ndims))
157 ALLOCATE(nvardimnam(ndims))
162 haction=
'get variable type'
163 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
167 haction=
'get variable dimensions identifiant'
168 status=nf_inq_vardimid(kcdf_id,idvar,nvardimid)
171 haction=
'get attributs'
173 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
176 allocate(hname(1:ngatts))
179 status=nf_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
182 if (trim(hname(jloop))==
'missing_value')
then
184 haction=
'get missing value'
185 status=nf_get_att_double(kcdf_id,idvar,
"missing_value",pmissvalue)
193 haction=
'get variable dimensions name'
194 status=nf_inq_dimname(kcdf_id,jloop2,nvardimnam(jloop2))
196 haction=
'get variable dimensions length'
197 status=nf_inq_dimlen(kcdf_id,jloop2,nvardimlen(jloop2))
203 ALLOCATE(zvalu2d(1:nvardimlen(1),1:nvardimlen(2)))
205 IF (kvartype>=5)
then
206 haction=
'get variable values (2D)'
207 status=nf_get_var_double(kcdf_id,idvar,zvalu2d(:,:))
210 pvalu2d(:,:)=zvalu2d(:,:)
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)
228 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
229 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
230 REAL,
DIMENSION(:),
INTENT(OUT) :: plon,plat
231 REAL,
DIMENSION(:),
INTENT(OUT) :: pval
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
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
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
258 IF (lhook) CALL dr_hook(
'MODE_READ_CDF:READ_LATLONVAL_CDF',0,zhook_handle)
261 haction=
'open netcdf'
262 status=nf_open(hfilename,nf_nowrite,kcdf_id)
265 if (status/=nf_noerr)
then
275 haction=
'get number of variables'
276 status=nf_inq_nvars(kcdf_id,nbvars)
279 ALLOCATE(varname(nbvars))
288 haction=
'get variables names'
289 status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
292 if (varname(jloop1)==hncvarname)
then
296 if (varname(jloop1)/=hncvarname)
then
297 if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
298 (scan(trim(varname(jloop1)),trim(hncvarname))==1))
then
305 if (id_vartoget1/=0)
then
306 id_vartoget=id_vartoget1
308 id_vartoget=id_vartoget2
310 if (id_vartoget==0)
then
311 haction=
'close netcdf'
312 status=nf_close(kcdf_id)
314 CALL
abor1_sfx(
'MODE_READ_NETCDF: READ_LATLONVAL_CDF')
324 haction=
'get variable dimensions number'
325 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
330 SELECT CASE (nvardims)
334 haction=
'get variable dimensions length'
335 status=nf_inq_dimlen(kcdf_id,nvardims,nlen)
337 ALLOCATE(zvalu(nlen))
339 CALL
get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
345 haction=
'get variable dimensions length'
346 status=nf_inq_dimlen(kcdf_id,jloop,nlen2d(jloop))
349 ALLOCATE(zvalu2d(nlen2d(1),nlen2d(2)))
350 ALLOCATE(zdim1(nlen2d(1)))
351 ALLOCATE(zdim2(nlen2d(2)))
353 CALL
get2dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,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
365 IF ((ilonfound==1).AND.(ilatfound==2))
then
367 DO jdim1=1,
SIZE(zdim1)
368 DO jdim2=1,
SIZE(zdim2)
370 pval(iarg)=zvalu2d(jdim1,jdim2)
371 plon(iarg)=zdim1(jdim1)
372 plat(iarg)=zdim2(jdim2)
375 ELSEIF ((ilonfound==2).AND.(ilatfound==1))
then
377 DO jdim1=1,
SIZE(zdim1)
378 DO jdim2=1,
SIZE(zdim2)
380 pval(iarg)=zvalu2d(jdim1,jdim2)
381 plat(iarg)=zdim1(jdim1)
382 plon(iarg)=zdim2(jdim2)
386 write(0,*)
'*****WARNING*****: incompatible dimensions to lat/lon/value arrays'
395 haction=
'close netcdf'
396 status=nf_close(kcdf_id)
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 )
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)
422 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
423 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
424 INTEGER,
INTENT(OUT):: kdim
429 character(len=80) :: haction
430 character(len=80),
DIMENSION(:),
ALLOCATABLE :: varname
431 integer ::jloop1,jloop
432 integer ::id_vartoget,id_vartoget1,id_vartoget2
434 integer,
DIMENSION(2) ::nlen2d
435 REAL(KIND=JPRB) :: zhook_handle
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
453 haction=
'get number of variables'
454 status=nf_inq_nvars(kcdf_id,nbvars)
457 ALLOCATE(varname(nbvars))
466 haction=
'get variables names'
467 status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
470 if (varname(jloop1)==hncvarname)
then
474 if (varname(jloop1)/=hncvarname)
then
475 if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
476 (scan(trim(varname(jloop1)),trim(hncvarname))==1))
then
483 if (id_vartoget1/=0)
then
484 id_vartoget=id_vartoget1
486 id_vartoget=id_vartoget2
488 if (id_vartoget==0)
then
489 haction=
'close netcdf'
490 status=nf_close(kcdf_id)
492 CALL
abor1_sfx(
'MODE_READ_CDF: READ_DIM_CDF')
502 haction=
'get variable dimensions number'
503 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
509 SELECT CASE (nvardims)
512 haction=
'get variable dimensions length'
513 status=nf_inq_dimlen(kcdf_id,nvardims,kdim)
520 haction=
'get variable dimensions length'
521 status=nf_inq_dimlen(kcdf_id,jloop,nlen2d(jloop))
523 kdim=kdim*nlen2d(jloop)
529 haction=
'close netcdf'
530 status=nf_close(kcdf_id)
537 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
538 IF (lhook) CALL dr_hook(
'MODE_READ_CDF:READ_DIM_CDF',1,zhook_handle)
subroutine get1dcdf(KCDF_ID, IDVAR, PMISSVALUE, PVALU1D)
subroutine handle_err_cdf(status, line)
subroutine abor1_sfx(YTEXT)
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)