15 USE yomhook
,ONLY : lhook, dr_hook
16 USE parkind1
,ONLY : jprb
25 INTEGER,
INTENT(IN) :: status
26 CHARACTER(LEN=80),
INTENT(IN) :: line
27 REAL(KIND=JPRB) :: zhook_handle
31 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:HANDLE_ERR_MER',0,zhook_handle)
32 IF (status /= nf_noerr)
THEN
33 CALL
abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: HANDLE_ERR_MER')
35 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:HANDLE_ERR_MER',1,zhook_handle)
40 SUBROUTINE get1dcdf(KCDF_ID,IDVAR,PMISSVALUE,PVALU1D)
45 INTEGER,
INTENT(IN) :: kcdf_id
46 INTEGER,
INTENT(IN) :: idvar
47 REAL,
INTENT(OUT) :: pmissvalue
48 REAL,
DIMENSION(:),
INTENT(OUT) :: pvalu1d
51 character(len=80) :: haction
52 integer,
save :: ndims=1
54 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid,nvardimlen
55 character(len=80),
DIMENSION(:),
ALLOCATABLE :: nvardimnam
58 character(len=80),
DIMENSION(:),
ALLOCATABLE :: hname
59 REAL,
DIMENSION(:),
ALLOCATABLE :: zvalu1d
60 REAL(KIND=JPRB) :: zhook_handle
64 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET1DCDF',0,zhook_handle)
66 ALLOCATE(nvardimid(ndims))
67 ALLOCATE(nvardimlen(ndims))
68 ALLOCATE(nvardimnam(ndims))
73 haction=
'get variable type'
74 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
78 haction=
'get variable dimensions name'
79 status=nf_inq_dimname(kcdf_id,idvar,nvardimnam(ndims))
82 haction=
'get variable dimensions length'
83 status=nf_inq_dimlen(kcdf_id,idvar,nvardimlen(ndims))
88 haction=
'get attributs'
89 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
92 allocate(hname(1:ngatts))
94 ALLOCATE(zvalu1d(1:nvardimlen(ndims)))
98 haction=
'get variable values (1D)'
99 status=nf_get_var_double(kcdf_id,idvar,zvalu1d(:))
103 pvalu1d(:)=zvalu1d(:)
105 IF (
ALLOCATED(zvalu1d ))
DEALLOCATE(zvalu1d)
106 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET1DCDF',1,zhook_handle)
112 SUBROUTINE get2dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
119 INTEGER,
INTENT(IN) :: kcdf_id
120 INTEGER,
INTENT(IN) :: idvar
121 REAL,
DIMENSION(:),
INTENT(OUT) :: pdim1,pdim2
122 CHARACTER(len=80),
INTENT(OUT) :: hdim1name,hdim2name
123 REAL,
INTENT(OUT) :: pmissvalue
124 REAL,
DIMENSION(:,:),
INTENT(OUT) :: pvalu2d
127 character(len=80) :: haction
128 integer,
save :: ndims=2
130 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid,nvardimlen
131 character(len=80),
DIMENSION(:),
ALLOCATABLE :: nvardimnam
132 integer :: jloop2, jloop, j1, j2
134 character(len=80),
DIMENSION(:),
ALLOCATABLE :: hname
135 real :: zmiss1,zmiss2
137 REAL,
DIMENSION(:,:),
ALLOCATABLE :: zvalu2d
138 INTEGER,
DIMENSION(:,:),
ALLOCATABLE :: ivalu2d
139 REAL(KIND=JPRB) :: zhook_handle
143 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET2DCDF',0,zhook_handle)
145 ALLOCATE(nvardimid(ndims))
146 ALLOCATE(nvardimlen(ndims))
147 ALLOCATE(nvardimnam(ndims))
152 haction=
'get variable type'
153 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
157 haction=
'get variable dimensions identifiant'
158 status=nf_inq_vardimid(kcdf_id,idvar,nvardimid)
161 haction=
'get attributs'
162 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
165 allocate(hname(1:ngatts))
170 status=nf_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
173 if (trim(hname(jloop))==
'missing_value')
then
175 haction=
'get missing value'
176 status=nf_get_att_double(kcdf_id,idvar,
"missing_value",pmissvalue)
180 if (trim(hname(jloop))==
'_FillValue')
then
182 haction=
'get _FillValue'
183 status=nf_get_att_double(kcdf_id,idvar,
"_FillValue",pmissvalue)
188 if (trim(hname(jloop))==
'scale_factor')
then
190 haction=
'get scale factor'
191 status=nf_get_att_double(kcdf_id,idvar,
"scale_factor",zscfa)
195 if (trim(hname(jloop))==
'add_offset')
then
198 status=nf_get_att_double(kcdf_id,idvar,
"add_offset",zoffs)
206 haction=
'get variable dimensions name'
207 status=nf_inq_dimname(kcdf_id,nvardimid(jloop2),nvardimnam(jloop2))
209 haction=
'get variable dimensions length'
210 status=nf_inq_dimlen(kcdf_id,nvardimid(jloop2),nvardimlen(jloop2))
216 IF (kvartype>=5)
then
217 ALLOCATE(zvalu2d(1:nvardimlen(1),1:nvardimlen(2)))
219 haction=
'get variable values (2D)'
220 status=nf_get_var_double(kcdf_id,idvar,zvalu2d(:,:))
223 ALLOCATE(ivalu2d(1:nvardimlen(1),1:nvardimlen(2)))
225 haction=
'get variable values (2D)'
226 status=nf_get_var_int(kcdf_id,idvar,ivalu2d(:,:))
230 DO j1=1,nvardimlen(1)
231 DO j2=1,nvardimlen(2)
232 IF (kvartype>=5)
THEN
233 IF (zvalu2d(j1,j2)/=pmissvalue) pvalu2d(j1,j2)=zvalu2d(j1,j2)*zscfa+zoffs
235 IF (zvalu2d(j1,j2)/=pmissvalue) pvalu2d(j1,j2)=ivalu2d(j1,j2)*zscfa+zoffs
240 CALL
get1dcdf(kcdf_id,nvardimid(1),zmiss1,pdim1)
241 CALL
get1dcdf(kcdf_id,nvardimid(2),zmiss2,pdim2)
242 hdim1name=nvardimnam(1)
243 hdim2name=nvardimnam(2)
244 IF (
ALLOCATED(zvalu2d ))
DEALLOCATE(zvalu2d)
245 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET2DCDF',1,zhook_handle)
251 SUBROUTINE get3dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
252 pdim3,hdim3name,pmissvalue,pvalu3d)
258 INTEGER,
INTENT(IN) :: kcdf_id
259 INTEGER,
INTENT(IN) :: idvar
260 REAL,
DIMENSION(:),
INTENT(OUT) :: pdim1,pdim2,pdim3
261 CHARACTER(len=80),
INTENT(OUT) :: hdim1name,hdim2name,hdim3name
262 REAL,
INTENT(OUT) :: pmissvalue
263 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: pvalu3d
266 character(len=80) :: haction
267 integer,
save :: ndims=3
269 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid,nvardimlen
270 character(len=80),
DIMENSION(:),
ALLOCATABLE :: nvardimnam
271 integer :: jloop2, jloop
274 character(len=80),
DIMENSION(:),
ALLOCATABLE :: hname
275 real :: zmiss1,zmiss2,zmiss3
277 REAL,
DIMENSION(:,:,:),
ALLOCATABLE :: zvalu3d
278 INTEGER,
DIMENSION(:,:,:),
ALLOCATABLE :: ivalu3d
279 REAL(KIND=JPRB) :: zhook_handle
283 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET3DCDF',0,zhook_handle)
285 ALLOCATE(nvardimid(ndims))
286 ALLOCATE(nvardimlen(ndims))
287 ALLOCATE(nvardimnam(ndims))
292 haction=
'get variable type'
293 status=nf_inq_vartype(kcdf_id,idvar,kvartype)
297 haction=
'get variable dimensions identifiant'
298 status=nf_inq_vardimid(kcdf_id,idvar,nvardimid)
302 haction=
'get attributs'
303 status=nf_inq_varnatts(kcdf_id,idvar,ngatts)
306 allocate(hname(1:ngatts))
311 status=nf_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
314 if (trim(hname(jloop))==
'missing_value')
then
316 haction=
'get missing value'
317 status=nf_get_att_double(kcdf_id,idvar,
"missing_value",pmissvalue)
321 if (trim(hname(jloop))==
'_FillValue')
then
323 haction=
'get _FillValue'
324 status=nf_get_att_double(kcdf_id,idvar,
"_FillValue",pmissvalue)
329 if (trim(hname(jloop))==
'scale_factor')
then
331 haction=
'get scale factor'
332 status=nf_get_att_double(kcdf_id,idvar,
"scale_factor",zscfa)
336 if (trim(hname(jloop))==
'add_offset')
then
339 status=nf_get_att_double(kcdf_id,idvar,
"add_offset",zoffs)
347 haction=
'get variable dimensions name'
348 status=nf_inq_dimname(kcdf_id,nvardimid(jloop2),nvardimnam(jloop2))
350 haction=
'get variable dimensions length'
351 status=nf_inq_dimlen(kcdf_id,nvardimid(jloop2),nvardimlen(jloop2))
357 IF (kvartype>=5)
then
358 ALLOCATE(zvalu3d(1:nvardimlen(1),1:nvardimlen(2),1:nvardimlen(3)))
360 haction=
'get variable values (3D)'
361 status=nf_get_var_double(kcdf_id,idvar,zvalu3d(:,:,:))
364 ALLOCATE(ivalu3d(1:nvardimlen(1),1:nvardimlen(2),1:nvardimlen(3)))
366 haction=
'get variable values (3D)'
367 status=nf_get_var_int(kcdf_id,idvar,ivalu3d(:,:,:))
371 pvalu3d(:,:,:)=xundef
372 DO j1=1,nvardimlen(1)
373 DO j2=1,nvardimlen(2)
374 DO j3=1,nvardimlen(3)
375 IF (kvartype>=5)
THEN
376 IF (zvalu3d(j1,j2,j3)/=pmissvalue) pvalu3d(j1,j2,j3)=zvalu3d(j1,j2,j3)*zscfa+zoffs
378 IF (ivalu3d(j1,j2,j3)/=pmissvalue) pvalu3d(j1,j2,j3)=ivalu3d(j1,j2,j3)*zscfa+zoffs
384 CALL
get1dcdf(kcdf_id,nvardimid(1),zmiss1,pdim1)
385 CALL
get1dcdf(kcdf_id,nvardimid(2),zmiss2,pdim2)
386 CALL
get1dcdf(kcdf_id,nvardimid(3),zmiss3,pdim3)
387 hdim1name=nvardimnam(1)
388 hdim2name=nvardimnam(2)
389 hdim3name=nvardimnam(3)
390 IF (
ALLOCATED(zvalu3d ))
DEALLOCATE(zvalu3d)
391 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET3DCDF',1,zhook_handle)
404 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
405 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
406 INTEGER,
INTENT(OUT):: kdim
411 character(len=80) :: haction
412 character(len=80),
DIMENSION(:),
ALLOCATABLE :: varname
413 integer ::jloop1,jloop
414 integer ::id_vartoget,id_vartoget1,id_vartoget2
416 integer,
DIMENSION(2) ::nlen2d
417 REAL(KIND=JPRB) :: zhook_handle
423 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_DIM_CDF',0,zhook_handle)
424 haction=
'open netcdf'
425 status=nf_open(hfilename,nf_nowrite,kcdf_id)
427 if (status/=nf_noerr)
then
437 haction=
'get number of variables'
438 status=nf_inq_nvars(kcdf_id,nbvars)
441 ALLOCATE(varname(nbvars))
450 haction=
'get variables names'
451 status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
454 if (varname(jloop1)==hncvarname)
then
458 if (varname(jloop1)/=hncvarname)
then
459 if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
460 (scan(trim(varname(jloop1)),trim(hncvarname))==1))
then
467 if (id_vartoget1/=0)
then
468 id_vartoget=id_vartoget1
470 id_vartoget=id_vartoget2
472 if (id_vartoget==0)
then
473 haction=
'close netcdf'
474 status=nf_close(kcdf_id)
476 CALL
abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_DIM_CDF')
486 haction=
'get variable dimensions number'
487 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
493 SELECT CASE (nvardims)
496 haction=
'get variable dimensions length'
497 status=nf_inq_dimlen(kcdf_id,id_vartoget,kdim)
504 haction=
'get variable dimensions length'
505 status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen2d(jloop))
507 kdim=kdim*nlen2d(jloop)
513 haction=
'close netcdf'
514 status=nf_close(kcdf_id)
521 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
522 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_DIM_CDF',1,zhook_handle)
536 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
537 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
542 character(len=80) :: haction
543 character(len=80),
DIMENSION(:),
ALLOCATABLE :: varname
544 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid
545 integer ::jloop1,jloop
546 integer ::id_vartoget,id_vartoget1,id_vartoget2
548 integer,
DIMENSION(3) ::ndimlen
549 character(LEN=80),
DIMENSION(3) :: ndimnam
552 real :: zzlamiss,zzlomiss
553 REAL(KIND=JPRB) :: zhook_handle
557 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',0,zhook_handle)
568 haction=
'open netcdf'
569 status=nf_open(hfilename,nf_nowrite,kcdf_id)
571 if (status/=nf_noerr)
then
581 haction=
'get number of variables'
582 status=nf_inq_nvars(kcdf_id,nbvars)
585 ALLOCATE(varname(nbvars))
594 haction=
'get variables names'
595 status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
598 if (varname(jloop1)==hncvarname)
then
602 if (varname(jloop1)/=hncvarname)
then
603 if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
604 (scan(trim(varname(jloop1)),trim(hncvarname))==1))
then
611 if (id_vartoget1/=0)
then
612 id_vartoget=id_vartoget1
614 id_vartoget=id_vartoget2
616 if (id_vartoget==0)
then
617 haction=
'close netcdf'
618 status=nf_close(kcdf_id)
620 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',1,zhook_handle)
632 haction=
'get variable dimensions number'
633 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
636 ALLOCATE(nvardimid(nvardims))
637 haction=
'get variable dimensions identifiant'
638 status=nf_inq_vardimid(kcdf_id,id_vartoget,nvardimid)
643 SELECT CASE (nvardims)
646 haction=
'get variable dimensions length'
647 status=nf_inq_dimlen(kcdf_id,id_vartoget,idim)
653 haction=
'get variable dimensions length'
654 status=nf_inq_dimlen(kcdf_id,nvardimid(jloop),ndimlen(jloop))
656 haction=
'get variable dimensions names'
657 status=nf_inq_dimname(kcdf_id,nvardimid(jloop),ndimnam(jloop))
659 if ((ndimnam(jloop)==
'lat').OR.(ndimnam(jloop)==
'latitude'))
then
660 ninlat=ndimlen(jloop)
661 if (.not.
allocated(xilatarray))
allocate(xilatarray(ndimlen(jloop)))
662 if (.not.
allocated(ninlon))
allocate(ninlon(ninlat))
663 CALL
get1dcdf(kcdf_id,nvardimid(jloop),zzlamiss,xilatarray(:))
665 if ((ndimnam(jloop)==
'lon').OR.(ndimnam(jloop)==
'longitude'))
then
667 if (.not.
allocated(xilonarray))
allocate(xilonarray(ndimlen(jloop)))
668 CALL
get1dcdf(kcdf_id,nvardimid(jloop),zzlomiss,xilonarray(:))
670 if (ndimnam(jloop)==
'depth') nindepth=ndimlen(jloop)
677 haction=
'close netcdf'
678 status=nf_close(kcdf_id)
685 nilength = nilength + ninlon(jloop1)
689 xilat2=xilatarray(
SIZE(xilatarray))
690 xilon2=xilonarray(
SIZE(xilonarray))
694 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
695 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',1,zhook_handle)
706 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
707 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
708 REAL,
DIMENSION(:),
INTENT(OUT) :: pval
713 character(len=80) :: haction
714 character(len=80),
DIMENSION(:),
ALLOCATABLE :: varname
716 integer ::id_vartoget,id_vartoget1,id_vartoget2
719 real,
DIMENSION(:),
ALLOCATABLE :: zvalu
721 REAL(KIND=JPRB) :: zhook_handle
727 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_Z1D_CDF',0,zhook_handle)
730 haction=
'open netcdf'
731 status=nf_open(hfilename,nf_nowrite,kcdf_id)
732 if (status/=nf_noerr)
then
738 haction=
'get number of variables'
739 status=nf_inq_nvars(kcdf_id,nbvars)
742 ALLOCATE(varname(nbvars))
749 haction=
'get variables names'
750 status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
752 if (varname(jloop1)==hncvarname)
then
755 if (varname(jloop1)/=hncvarname)
then
756 if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
757 (scan(trim(varname(jloop1)),trim(hncvarname))==1))
then
762 if (id_vartoget1/=0)
then
763 id_vartoget=id_vartoget1
765 id_vartoget=id_vartoget2
767 if (id_vartoget==0)
then
768 haction=
'close netcdf'
769 status=nf_close(kcdf_id)
771 CALL
abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_Z1D_CDF')
778 haction=
'get variable dimensions number'
779 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
784 SELECT CASE (nvardims)
787 haction=
'get variable dimensions length'
788 status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen)
790 ALLOCATE(zvalu(nlen))
792 CALL
get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
796 write(0,*)
'YOU ARE TRYING TO READ A 2D FIELD FOR :', trim(hncvarname)
797 CALL
abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_Z1D_CDF')
802 haction=
'close netcdf'
803 status=nf_close(kcdf_id)
808 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
809 IF (
ALLOCATED(zvalu ))
DEALLOCATE(zvalu )
811 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_Z1D_CDF',1,zhook_handle)
821 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
822 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
823 REAL,
DIMENSION(:),
INTENT(OUT) :: plon,plat
824 REAL,
DIMENSION(:),
INTENT(OUT) :: pval
829 character(len=80) :: haction
830 character(len=80),
DIMENSION(:),
ALLOCATABLE :: varname
831 integer ::jloop1,jdim1,jdim2,jloop
832 integer ::id_vartoget,id_vartoget1,id_vartoget2
835 integer,
DIMENSION(2) ::nlen2d
836 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid,nvardimlen
837 character(len=80),
DIMENSION(:),
ALLOCATABLE :: nvardimnam
838 real,
DIMENSION(:),
ALLOCATABLE :: zvalu
839 real,
DIMENSION(:,:),
ALLOCATABLE :: zvalu2d
841 real,
DIMENSION(:),
ALLOCATABLE :: zdim1
842 real,
DIMENSION(:),
ALLOCATABLE :: zdim2
843 character(len=80) :: ydim1name,ydim2name
844 integer :: ilonfound,ilatfound, iarg
845 REAL(KIND=JPRB) :: zhook_handle
853 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONVAL_CDF',0,zhook_handle)
856 haction=
'open netcdf'
857 status=nf_open(hfilename,nf_nowrite,kcdf_id)
860 if (status/=nf_noerr)
then
870 haction=
'get number of variables'
871 status=nf_inq_nvars(kcdf_id,nbvars)
874 ALLOCATE(varname(nbvars))
883 haction=
'get variables names'
884 status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
887 if (varname(jloop1)==hncvarname)
then
891 if (varname(jloop1)/=hncvarname)
then
892 if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
893 (scan(trim(varname(jloop1)),trim(hncvarname))==1))
then
900 if (id_vartoget1/=0)
then
901 id_vartoget=id_vartoget1
903 id_vartoget=id_vartoget2
905 if (id_vartoget==0)
then
906 haction=
'close netcdf'
907 status=nf_close(kcdf_id)
909 CALL
abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_LATLONVAL_CDF')
919 haction=
'get variable dimensions number'
920 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
926 SELECT CASE (nvardims)
929 haction=
'get variable dimensions length'
930 status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen)
932 ALLOCATE(zvalu(nlen))
934 CALL
get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
939 haction=
'get variable dimensions length'
940 status=nf_inq_dimlen(kcdf_id,id_vartoget,nlen2d(jloop))
943 ALLOCATE(zvalu2d(nlen2d(1),nlen2d(2)))
944 ALLOCATE(zdim1(nlen2d(1)))
945 ALLOCATE(zdim2(nlen2d(2)))
947 CALL
get2dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,ydim2name,&
951 if ((ydim1name==
'lon').OR.(ydim1name==
'longitude')) ilonfound=1
952 if ((ydim2name==
'lon').OR.(ydim2name==
'longitude')) ilonfound=2
953 if ((ydim1name==
'lat').OR.(ydim1name==
'latitude')) ilatfound=1
954 if ((ydim2name==
'lat').OR.(ydim2name==
'latitude')) ilatfound=2
959 IF ((ilonfound==1).AND.(ilatfound==2))
then
961 DO jdim1=1,
SIZE(zdim1)
962 DO jdim2=1,
SIZE(zdim2)
964 pval(iarg)=zvalu2d(jdim1,jdim2)
965 plon(iarg)=zdim1(jdim1)
966 plat(iarg)=zdim2(jdim2)
969 ELSEIF ((ilonfound==2).AND.(ilatfound==1))
then
971 DO jdim1=1,
SIZE(zdim1)
972 DO jdim2=1,
SIZE(zdim2)
974 pval(iarg)=zvalu2d(jdim1,jdim2)
975 plat(iarg)=zdim1(jdim1)
976 plon(iarg)=zdim2(jdim2)
980 write(0,*)
'*****WARNING*****: incompatible dimensions to lat/lon/value arrays'
989 haction=
'close netcdf'
990 status=nf_close(kcdf_id)
997 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
998 IF (
ALLOCATED(zvalu ))
DEALLOCATE(zvalu )
999 IF (
ALLOCATED(zvalu2d ))
DEALLOCATE(zvalu2d)
1000 IF (
ALLOCATED(zdim1 ))
DEALLOCATE(zdim1 )
1001 IF (
ALLOCATED(zdim2 ))
DEALLOCATE(zdim2 )
1004 IF (
ALLOCATED(nvardimid ))
DEALLOCATE(nvardimid )
1005 IF (
ALLOCATED(nvardimnam ))
DEALLOCATE(nvardimnam)
1006 IF (
ALLOCATED(nvardimlen ))
DEALLOCATE(nvardimlen)
1007 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONVAL_CDF',1,zhook_handle)
1017 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
1018 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
1019 REAL,
DIMENSION(:),
INTENT(OUT) :: plon,plat
1020 REAL,
DIMENSION(:),
INTENT(OUT) :: pdep
1021 REAL,
DIMENSION(:,:),
INTENT(OUT) :: pval
1026 character(len=80) :: haction
1027 character(len=80),
DIMENSION(:),
ALLOCATABLE :: varname
1028 integer ::jloop1,jdim1,jdim2,jdim3,jloop
1030 integer ::id_vartoget,id_vartoget1,id_vartoget2
1032 integer,
DIMENSION(3) ::nlen3d
1033 integer,
DIMENSION(:),
ALLOCATABLE :: nvardimid,nvardimlen
1034 character(len=80),
DIMENSION(:),
ALLOCATABLE :: nvardimnam
1035 real,
DIMENSION(:,:,:),
ALLOCATABLE :: zvalu3d
1037 real,
DIMENSION(:),
ALLOCATABLE :: zdim1
1038 real,
DIMENSION(:),
ALLOCATABLE :: zdim2
1039 real,
DIMENSION(:),
ALLOCATABLE :: zdim3
1040 character(len=80) :: ydim1name,ydim2name,ydim3name
1041 integer :: ilonfound,ilatfound,idepfound
1043 REAL(KIND=JPRB) :: zhook_handle
1045 include
'netcdf.inc'
1051 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONDEPVAL_CDF',0,zhook_handle)
1052 haction=
'open netcdf'
1053 status=nf_open(hfilename,nf_nowrite,kcdf_id)
1061 haction=
'get number of variables'
1062 status=nf_inq_nvars(kcdf_id,nbvars)
1065 ALLOCATE(varname(nbvars))
1074 haction=
'get variables names'
1075 status=nf_inq_varname(kcdf_id,jloop1,varname(jloop1))
1078 if (varname(jloop1)==hncvarname)
then
1082 if (varname(jloop1)/=hncvarname)
then
1083 if((lgt(trim(varname(jloop1)),trim(hncvarname))).AND.&
1084 (scan(trim(varname(jloop1)),trim(hncvarname))==1))
then
1091 if (id_vartoget1/=0)
then
1092 id_vartoget=id_vartoget1
1094 id_vartoget=id_vartoget2
1096 if (id_vartoget==0)
then
1097 haction=
'close netcdf'
1098 status=nf_close(kcdf_id)
1100 CALL
abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_LATLONDEPVAL_CDF')
1110 haction=
'get variable dimensions number'
1111 status=nf_inq_varndims(kcdf_id,id_vartoget,nvardims)
1114 ALLOCATE(nvardimid(nvardims))
1115 haction=
'get variable dimensions identifiant'
1116 status=nf_inq_vardimid(kcdf_id,id_vartoget,nvardimid)
1122 SELECT CASE (nvardims)
1125 write(0,*)
'********************************************'
1126 write(0,*)
'* number of dimension to low: ',nvardims,
' *'
1127 write(0,*)
'* you need a 3-dimension variable *'
1128 write(0,*)
'********************************************'
1132 haction=
'get variable dimensions length'
1133 status=nf_inq_dimlen(kcdf_id,nvardimid(jloop),nlen3d(jloop))
1136 ALLOCATE(zvalu3d(nlen3d(1),nlen3d(2),nlen3d(3)))
1137 ALLOCATE(zdim1(nlen3d(1)))
1138 ALLOCATE(zdim2(nlen3d(2)))
1139 ALLOCATE(zdim3(nlen3d(3)))
1141 CALL
get3dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,ydim2name,&
1142 zdim3,ydim3name,zmiss,zvalu3d)
1146 if ((ydim1name==
'lon').OR.(ydim1name==
'longitude')) ilonfound=1
1147 if ((ydim2name==
'lon').OR.(ydim2name==
'longitude')) ilonfound=2
1148 if ((ydim3name==
'lon').OR.(ydim3name==
'longitude')) ilonfound=3
1149 if ((ydim1name==
'lat').OR.(ydim1name==
'latitude')) ilatfound=1
1150 if ((ydim2name==
'lat').OR.(ydim2name==
'latitude')) ilatfound=2
1151 if ((ydim3name==
'lat').OR.(ydim3name==
'latitude')) ilatfound=3
1152 if (ydim1name==
'depth') idepfound=1
1153 if (ydim2name==
'depth') idepfound=2
1154 if (ydim3name==
'depth') idepfound=3
1160 IF ((ilonfound==1).AND.(ilatfound==2).AND.(idepfound==3))
then
1164 DO jdim2=1,
SIZE(zdim2)
1165 DO jdim1=1,
SIZE(zdim1)
1167 plon(iarg)=zdim1(jdim1)
1168 plat(iarg)=zdim2(jdim2)
1169 DO jdim3=1,
SIZE(zdim3)
1170 pval(iarg,jdim3)=zvalu3d(jdim1,jdim2,jdim3)
1176 ELSEIF ((ilonfound==2).AND.(ilatfound==1).AND.(idepfound==3))
then
1178 DO jdim1=1,
SIZE(zdim1)
1179 DO jdim2=1,
SIZE(zdim2)
1181 plon(iarg)=zdim2(jdim2)
1182 plat(iarg)=zdim1(jdim1)
1183 DO jdim3=1,
SIZE(zdim3)
1184 pval(iarg,jdim3)=zvalu3d(jdim1,jdim2,jdim3)
1189 ELSEIF ((ilonfound==1).AND.(ilatfound==3).AND.(idepfound==2))
then
1191 DO jdim3=1,
SIZE(zdim3)
1192 DO jdim1=1,
SIZE(zdim1)
1194 plon(iarg)=zdim1(jdim1)
1195 plat(iarg)=zdim3(jdim3)
1196 DO jdim2=1,
SIZE(zdim2)
1197 pval(iarg,jdim2)=zvalu3d(jdim1,jdim2,jdim3)
1202 ELSEIF ((ilatfound==1).AND.(ilonfound==3).AND.(idepfound==2))
then
1204 DO jdim1=1,
SIZE(zdim1)
1205 DO jdim3=1,
SIZE(zdim3)
1207 plon(iarg)=zdim3(jdim3)
1208 plat(iarg)=zdim1(jdim1)
1209 DO jdim2=1,
SIZE(zdim2)
1210 pval(iarg,jdim2)=zvalu3d(jdim1,jdim2,jdim3)
1215 ELSEIF ((ilonfound==2).AND.(ilatfound==3).AND.(idepfound==1))
then
1217 DO jdim3=1,
SIZE(zdim3)
1218 DO jdim2=1,
SIZE(zdim2)
1220 plon(iarg)=zdim2(jdim2)
1221 plat(iarg)=zdim3(jdim3)
1222 DO jdim1=1,
SIZE(zdim1)
1223 pval(iarg,jdim1)=zvalu3d(jdim1,jdim2,jdim3)
1228 ELSEIF ((ilatfound==2).AND.(ilonfound==3).AND.(idepfound==1))
then
1230 DO jdim2=1,
SIZE(zdim2)
1231 DO jdim3=1,
SIZE(zdim3)
1233 plon(iarg)=zdim3(jdim3)
1234 plat(iarg)=zdim2(jdim2)
1235 DO jdim1=1,
SIZE(zdim1)
1236 pval(iarg,jdim1)=zvalu3d(jdim1,jdim2,jdim3)
1242 write(0,*)
'*****WARNING*****: incompatible dimensions to lat/lon/value arrays'
1250 haction=
'close netcdf'
1252 status=nf_close(kcdf_id)
1259 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
1260 IF (
ALLOCATED(zvalu3d ))
DEALLOCATE(zvalu3d)
1261 IF (
ALLOCATED(zdim1 ))
DEALLOCATE(zdim1 )
1262 IF (
ALLOCATED(zdim2 ))
DEALLOCATE(zdim2 )
1263 IF (
ALLOCATED(zdim3 ))
DEALLOCATE(zdim3 )
1266 IF (
ALLOCATED(nvardimid ))
DEALLOCATE(nvardimid )
1267 IF (
ALLOCATED(nvardimnam ))
DEALLOCATE(nvardimnam)
1268 IF (
ALLOCATED(nvardimlen ))
DEALLOCATE(nvardimlen)
1269 201
FORMAT(4(3x,f10.4))
1270 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONDEPVAL_CDF',1,zhook_handle)
1285 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
1286 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
1287 REAL,
POINTER,
DIMENSION(:) :: pfield
1289 REAL,
DIMENSION(:),
ALLOCATABLE :: zlati
1290 REAL,
DIMENSION(:),
ALLOCATABLE :: zlong
1291 REAL,
TARGET,
DIMENSION(:,:),
ALLOCATABLE :: zvalue
1292 REAL,
DIMENSION(:),
ALLOCATABLE :: zdepth
1293 REAL,
TARGET,
DIMENSION(:),
ALLOCATABLE :: zval
1297 REAL(KIND=JPRB) :: zhook_handle
1300 include
'netcdf.inc'
1302 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_SST',0,zhook_handle)
1303 IF (nilength<0)
then
1306 cinterp_type=
'UNIF '
1307 ELSEIF (nindepth<0)
THEN
1308 ALLOCATE(zlati(nilength) )
1309 ALLOCATE(zlong(nilength) )
1310 ALLOCATE(zval(nilength) )
1312 ALLOCATE(pfield(nilength))
1315 IF (trim(hncvarname) ==
'temperature')
THEN
1316 WHERE (zval(:)/=zundef .AND. zval(:)<100.) pfield(:)=pfield(:)+xtt
1318 cinterp_type=
'HORIBL'
1322 ALLOCATE(zvalue(nilength,nindepth))
1323 ALLOCATE(zlati(nilength) )
1324 ALLOCATE(zlong(nilength) )
1325 ALLOCATE(zdepth(nindepth))
1329 ALLOCATE(pfield(nilength))
1331 pfield(:)=zvalue(:,1)
1332 IF (trim(hncvarname) ==
'temperature')
THEN
1333 WHERE (zvalue(:,1)/=zundef .AND. zvalue(:,1)<100.) pfield(:)=pfield(:)+xtt
1335 cinterp_type=
'HORIBL'
1340 IF (
ALLOCATED(zvalue ))
DEALLOCATE(zvalue )
1341 IF (
ALLOCATED(zlong ))
DEALLOCATE(zlong )
1342 IF (
ALLOCATED(zlati ))
DEALLOCATE(zlati )
1343 IF (
ALLOCATED(zdepth ))
DEALLOCATE(zdepth )
1344 IF (
ALLOCATED(zval ))
DEALLOCATE(zval )
1346 202
FORMAT(3(3x,f10.4))
1347 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_SST',1,zhook_handle)
1361 CHARACTER(LEN=28),
INTENT(IN) :: hfilename
1362 CHARACTER(LEN=28),
INTENT(IN) :: hncvarname
1363 REAL,
POINTER,
DIMENSION(:) :: pfield
1365 REAL,
DIMENSION(:),
ALLOCATABLE :: zlati
1366 REAL,
DIMENSION(:),
ALLOCATABLE :: zlong
1367 REAL,
TARGET,
DIMENSION(:),
ALLOCATABLE:: zvalue
1369 REAL(KIND=JPRB) :: zhook_handle
1372 include
'netcdf.inc'
1374 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_ZS_SEA',0,zhook_handle)
1379 cinterp_type=
'UNIF '
1380 elseif(nilength>0)
then
1381 ALLOCATE(zvalue(nilength))
1382 ALLOCATE(zlati(nilength) )
1383 ALLOCATE(zlong(nilength) )
1386 ALLOCATE(pfield(nilength))
1388 cinterp_type=
'HORIBL'
1394 cinterp_type=
'UNIF '
1397 IF (
ALLOCATED(zvalue ))
DEALLOCATE(zvalue )
1398 IF (
ALLOCATED(zlong ))
DEALLOCATE(zlong )
1399 IF (
ALLOCATED(zlati ))
DEALLOCATE(zlati )
1401 202
FORMAT(3(3x,f10.4))
1402 IF (lhook) CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_ZS_SEA',1,zhook_handle)
subroutine get1dcdf(KCDF_ID, IDVAR, PMISSVALUE, PVALU1D)
subroutine get3dcdf(KCDF_ID, IDVAR, PDIM1, HDIM1NAME, PDIM2, HDIM2NAME, PDIM3, HDIM3NAME, PMISSVALUE, PVALU3D)
subroutine read_netcdf_sst(HFILENAME, HNCVARNAME, PFIELD)
subroutine abor1_sfx(YTEXT)
subroutine read_latlondepval_cdf(HFILENAME, HNCVARNAME, PLON, PLAT, PDEP, PVAL)
subroutine read_dim_cdf(HFILENAME, HNCVARNAME, KDIM)
subroutine prep_netcdf_grid(HFILENAME, HNCVARNAME)
subroutine read_latlonval_cdf(HFILENAME, HNCVARNAME, PLON, PLAT, PVAL)
subroutine read_z1d_cdf(HFILENAME, HNCVARNAME, PVAL)
subroutine handle_err_mer(status, line)
subroutine get2dcdf(KCDF_ID, IDVAR, PDIM1, HDIM1NAME, PDIM2, HDIM2NAME, PMISSVALUE, PVALU2D)
subroutine read_netcdf_zs_sea(HFILENAME, HNCVARNAME, PFIELD)