27 INTEGER,
INTENT(IN) :: status
28 CHARACTER(LEN=80),
INTENT(IN) :: line
29 REAL(KIND=JPRB) :: ZHOOK_HANDLE
32 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:HANDLE_ERR_MER',0,zhook_handle)
33 IF (status /= nf90_noerr)
THEN 34 CALL abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: HANDLE_ERR_MER')
36 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:HANDLE_ERR_MER',1,zhook_handle)
41 SUBROUTINE get1dcdf(KCDF_ID,IDVAR,PMISSVALUE,PVALU1D)
48 INTEGER,
INTENT(IN) :: KCDF_ID
49 INTEGER,
INTENT(IN) :: IDVAR
50 REAL,
INTENT(OUT) :: PMISSVALUE
51 REAL,
DIMENSION(:),
INTENT(OUT) :: PVALU1D
54 character(len=80) :: HACTION
55 integer,
save :: NDIMS=1
57 integer,
DIMENSION(:),
ALLOCATABLE :: NVARDIMID,NVARDIMLEN
58 character(len=80),
DIMENSION(:),
ALLOCATABLE :: NVARDIMNAM
61 integer,
dimension(1) :: NDIMID
62 character(len=80),
DIMENSION(:),
ALLOCATABLE :: HNAME
63 REAL,
DIMENSION(:),
ALLOCATABLE :: ZVALU1D
64 REAL(KIND=JPRB) :: ZHOOK_HANDLE
67 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET1DCDF',0,zhook_handle)
69 ALLOCATE(nvardimid(ndims))
70 ALLOCATE(nvardimlen(ndims))
71 ALLOCATE(nvardimnam(ndims))
76 haction=
'get variable type' 77 status=nf90_inquire_variable(kcdf_id,idvar,xtype=kvartype)
81 status=nf90_inquire_variable(kcdf_id,idvar,dimids=ndimid)
82 haction=
'get variable dimensions name' 83 status=nf90_inquire_dimension(kcdf_id,ndimid(ndims),
name=nvardimnam(ndims))
86 haction=
'get variable dimensions length' 87 status=nf90_inquire_dimension(kcdf_id,ndimid(ndims),len=nvardimlen(ndims))
92 haction=
'get attributs' 93 status=nf90_inquire_variable(kcdf_id,idvar,natts=ngatts)
96 allocate(hname(1:ngatts))
98 ALLOCATE(zvalu1d(1:nvardimlen(ndims)))
101 IF (kvartype>=5)
then 102 haction=
'get variable values (1D)' 103 status=nf90_get_var(kcdf_id,idvar,zvalu1d(:))
107 pvalu1d(:)=zvalu1d(:)
109 IF (
ALLOCATED(zvalu1d ))
DEALLOCATE(zvalu1d)
110 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET1DCDF',1,zhook_handle)
116 SUBROUTINE get2dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
125 INTEGER,
INTENT(IN) :: KCDF_ID
126 INTEGER,
INTENT(IN) :: IDVAR
127 REAL,
DIMENSION(:),
INTENT(OUT) :: PDIM1,PDIM2
128 CHARACTER(len=80),
INTENT(OUT) :: HDIM1NAME,HDIM2NAME
129 REAL,
INTENT(OUT) :: PMISSVALUE
130 REAL,
DIMENSION(:,:),
INTENT(OUT) :: PVALU2D
133 character(len=80) :: HACTION
134 integer,
save :: NDIMS=2
136 integer,
DIMENSION(:),
ALLOCATABLE :: NVARDIMID,NVARDIMLEN
137 character(len=80),
DIMENSION(:),
ALLOCATABLE :: NVARDIMNAM
138 integer :: JLOOP2, JLOOP, J1, J2
140 character(len=80),
DIMENSION(:),
ALLOCATABLE :: HNAME
141 real :: ZMISS1,ZMISS2
143 REAL,
DIMENSION(:,:),
ALLOCATABLE :: ZVALU2D
144 INTEGER,
DIMENSION(:,:),
ALLOCATABLE :: IVALU2D
145 REAL(KIND=JPRB) :: ZHOOK_HANDLE
148 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET2DCDF',0,zhook_handle)
150 ALLOCATE(nvardimid(ndims))
151 ALLOCATE(nvardimlen(ndims))
152 ALLOCATE(nvardimnam(ndims))
157 haction=
'get variable type' 158 status=nf90_inquire_variable(kcdf_id,idvar,xtype=kvartype)
162 haction=
'get variable dimensions identifiant' 163 status=nf90_inquire_variable(kcdf_id,idvar,dimids=nvardimid)
166 haction=
'get attributs' 167 status=nf90_inquire_variable(kcdf_id,idvar,natts=ngatts)
170 allocate(hname(1:ngatts))
175 status=nf90_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
178 if (
trim(hname(jloop))==
'missing_value')
then 180 haction=
'get missing value' 181 status=nf90_get_att(kcdf_id,idvar,
"missing_value",pmissvalue)
185 if (
trim(hname(jloop))==
'_FillValue')
then 187 haction=
'get _FillValue' 188 status=nf90_get_att(kcdf_id,idvar,
"_FillValue",pmissvalue)
193 if (
trim(hname(jloop))==
'scale_factor')
then 195 haction=
'get scale factor' 196 status=nf90_get_att(kcdf_id,idvar,
"scale_factor",zscfa)
200 if (
trim(hname(jloop))==
'add_offset')
then 203 status=nf90_get_att(kcdf_id,idvar,
"add_offset",zoffs)
211 haction=
'get variable dimensions name' 212 status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop2),
name=nvardimnam(jloop2))
214 haction=
'get variable dimensions length' 215 status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop2),len=nvardimlen(jloop2))
221 IF (kvartype>=5)
then 222 ALLOCATE(zvalu2d(1:nvardimlen(1),1:nvardimlen(2)))
224 haction=
'get variable values (2D)' 225 status=nf90_get_var(kcdf_id,idvar,zvalu2d(:,:))
228 ALLOCATE(ivalu2d(1:nvardimlen(1),1:nvardimlen(2)))
230 haction=
'get variable values (2D)' 231 status=nf90_get_var(kcdf_id,idvar,ivalu2d(:,:))
235 DO j1=1,nvardimlen(1)
236 DO j2=1,nvardimlen(2)
237 IF (kvartype>=5)
THEN 238 IF (zvalu2d(j1,j2)/=pmissvalue) pvalu2d(j1,j2)=zvalu2d(j1,j2)*zscfa+zoffs
240 IF (zvalu2d(j1,j2)/=pmissvalue) pvalu2d(j1,j2)=ivalu2d(j1,j2)*zscfa+zoffs
245 CALL get1dcdf(kcdf_id,nvardimid(1),zmiss1,pdim1)
246 CALL get1dcdf(kcdf_id,nvardimid(2),zmiss2,pdim2)
247 hdim1name=nvardimnam(1)
248 hdim2name=nvardimnam(2)
249 IF (
ALLOCATED(zvalu2d ))
DEALLOCATE(zvalu2d)
250 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET2DCDF',1,zhook_handle)
256 SUBROUTINE get3dcdf(KCDF_ID,IDVAR,PDIM1,HDIM1NAME,PDIM2,HDIM2NAME,&
257 PDIM3,HDIM3NAME,PMISSVALUE,PVALU3D)
265 INTEGER,
INTENT(IN) :: KCDF_ID
266 INTEGER,
INTENT(IN) :: IDVAR
267 REAL,
DIMENSION(:),
INTENT(OUT) :: PDIM1,PDIM2,PDIM3
268 CHARACTER(len=80),
INTENT(OUT) :: HDIM1NAME,HDIM2NAME,HDIM3NAME
269 REAL,
INTENT(OUT) :: PMISSVALUE
270 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: PVALU3D
273 character(len=80) :: HACTION
274 integer,
save :: NDIMS=3
276 integer,
DIMENSION(:),
ALLOCATABLE :: NVARDIMID,NVARDIMLEN
277 character(len=80),
DIMENSION(:),
ALLOCATABLE :: NVARDIMNAM
278 integer :: JLOOP2, JLOOP
281 character(len=80),
DIMENSION(:),
ALLOCATABLE :: HNAME
282 real :: ZMISS1,ZMISS2,ZMISS3
284 REAL,
DIMENSION(:,:,:),
ALLOCATABLE :: ZVALU3D
285 INTEGER,
DIMENSION(:,:,:),
ALLOCATABLE :: IVALU3D
286 REAL(KIND=JPRB) :: ZHOOK_HANDLE
289 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET3DCDF',0,zhook_handle)
291 ALLOCATE(nvardimid(ndims))
292 ALLOCATE(nvardimlen(ndims))
293 ALLOCATE(nvardimnam(ndims))
298 haction=
'get variable type' 299 status=nf90_inquire_variable(kcdf_id,idvar,xtype=kvartype)
303 haction=
'get variable dimensions identifiant' 304 status=nf90_inquire_variable(kcdf_id,idvar,dimids=nvardimid)
308 haction=
'get attributs' 309 status=nf90_inquire_variable(kcdf_id,idvar,natts=ngatts)
312 allocate(hname(1:ngatts))
317 status=nf90_inq_attname(kcdf_id,idvar,jloop,hname(jloop))
320 if (
trim(hname(jloop))==
'missing_value')
then 322 haction=
'get missing value' 323 status=nf90_get_att(kcdf_id,idvar,
"missing_value",pmissvalue)
327 if (
trim(hname(jloop))==
'_FillValue')
then 329 haction=
'get _FillValue' 330 status=nf90_get_att(kcdf_id,idvar,
"_FillValue",pmissvalue)
335 if (
trim(hname(jloop))==
'scale_factor')
then 337 haction=
'get scale factor' 338 status=nf90_get_att(kcdf_id,idvar,
"scale_factor",zscfa)
342 if (
trim(hname(jloop))==
'add_offset')
then 345 status=nf90_get_att(kcdf_id,idvar,
"add_offset",zoffs)
353 haction=
'get variable dimensions name' 354 status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop2),
name=nvardimnam(jloop2))
356 haction=
'get variable dimensions length' 357 status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop2),len=nvardimlen(jloop2))
363 IF (kvartype>=5)
then 364 ALLOCATE(zvalu3d(1:nvardimlen(1),1:nvardimlen(2),1:nvardimlen(3)))
366 haction=
'get variable values (3D)' 367 status=nf90_get_var(kcdf_id,idvar,zvalu3d(:,:,:))
370 ALLOCATE(ivalu3d(1:nvardimlen(1),1:nvardimlen(2),1:nvardimlen(3)))
372 haction=
'get variable values (3D)' 373 status=nf90_get_var(kcdf_id,idvar,ivalu3d(:,:,:))
378 DO j1=1,nvardimlen(1)
379 DO j2=1,nvardimlen(2)
380 DO j3=1,nvardimlen(3)
381 IF (kvartype>=5)
THEN 382 IF (zvalu3d(j1,j2,j3)/=pmissvalue) pvalu3d(j1,j2,j3)=zvalu3d(j1,j2,j3)*zscfa+zoffs
384 IF (ivalu3d(j1,j2,j3)/=pmissvalue) pvalu3d(j1,j2,j3)=ivalu3d(j1,j2,j3)*zscfa+zoffs
390 CALL get1dcdf(kcdf_id,nvardimid(1),zmiss1,pdim1)
391 CALL get1dcdf(kcdf_id,nvardimid(2),zmiss2,pdim2)
392 CALL get1dcdf(kcdf_id,nvardimid(3),zmiss3,pdim3)
393 hdim1name=nvardimnam(1)
394 hdim2name=nvardimnam(2)
395 hdim3name=nvardimnam(3)
396 IF (
ALLOCATED(zvalu3d ))
DEALLOCATE(zvalu3d)
397 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:GET3DCDF',1,zhook_handle)
412 CHARACTER(LEN=28),
INTENT(IN) :: HFILENAME
413 CHARACTER(LEN=28),
INTENT(IN) :: HNCVARNAME
414 INTEGER,
INTENT(OUT):: KDIM
419 character(len=80) :: HACTION
420 character(len=80),
DIMENSION(:),
ALLOCATABLE :: YVARNAME
421 integer ::JLOOP1,JLOOP
422 integer ::ID_VARTOGET,ID_VARTOGET1,ID_VARTOGET2
424 INTEGER,
DIMENSION(1) :: NDIMID
425 integer,
DIMENSION(2) ::NLEN2D, NDIMID2D
426 REAL(KIND=JPRB) :: ZHOOK_HANDLE
431 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_DIM_CDF',0,zhook_handle)
432 haction=
'open netcdf' 433 status=nf90_open(hfilename,nf90_nowrite,kcdf_id)
435 if (status/=nf90_noerr)
then 445 haction=
'get number of variables' 446 status=nf90_inquire(kcdf_id,nvariables=nbvars)
449 ALLOCATE(yvarname(nbvars))
458 haction=
'get variables names' 459 status=nf90_inquire_variable(kcdf_id,jloop1,
name=yvarname(jloop1))
462 if (yvarname(jloop1)==hncvarname)
then 466 if (yvarname(jloop1)/=hncvarname)
then 467 if((lgt(
trim(yvarname(jloop1)),
trim(hncvarname))).AND.&
468 (scan(
trim(yvarname(jloop1)),
trim(hncvarname))==1))
then 475 if (id_vartoget1/=0)
then 476 id_vartoget=id_vartoget1
478 id_vartoget=id_vartoget2
480 if (id_vartoget==0)
then 481 haction=
'close netcdf' 482 status=nf90_close(kcdf_id)
484 CALL abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_DIM_CDF')
494 haction=
'get variable dimensions number' 495 status=nf90_inquire_variable(kcdf_id,id_vartoget,ndims=nvardims)
501 SELECT CASE (nvardims)
504 haction=
'get variable dimensions length' 505 status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=ndimid)
506 status=nf90_inquire_dimension(kcdf_id,ndimid(1),len=kdim)
512 status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=ndimid2d)
514 haction=
'get variable dimensions length' 515 status=nf90_inquire_dimension(kcdf_id,ndimid2d(jloop),len=nlen2d(jloop))
517 kdim=kdim*nlen2d(jloop)
523 haction=
'close netcdf' 524 status=nf90_close(kcdf_id)
531 IF (
ALLOCATED(yvarname ))
DEALLOCATE(yvarname)
532 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_DIM_CDF',1,zhook_handle)
549 USE modi_horibl_surf_init
550 USE modi_horibl_surf_coef
560 CHARACTER(LEN=28),
INTENT(IN) :: HFILENAME
561 CHARACTER(LEN=28),
INTENT(IN) :: HNCVARNAME
566 character(len=80) :: HACTION
567 character(len=80),
DIMENSION(:),
ALLOCATABLE :: YVARNAME
568 integer,
DIMENSION(:),
ALLOCATABLE :: NVARDIMID
569 integer ::JLOOP1,JLOOP
570 integer ::ID_VARTOGET,ID_VARTOGET1,ID_VARTOGET2
572 integer,
DIMENSION(3) ::INDIMLEN
573 character(LEN=80),
DIMENSION(3) :: NDIMNAM
576 INTEGER :: IINLA, INO
577 real :: ZZLAMISS,ZZLOMISS
579 REAL(KIND=JPRB) :: ZHOOK_HANDLE
582 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',0,zhook_handle)
595 haction=
'open netcdf' 596 status=nf90_open(hfilename,nf90_nowrite,kcdf_id)
598 if (status/=nf90_noerr)
then 608 haction=
'get number of variables' 609 status=nf90_inquire(kcdf_id,nvariables=inbvars)
612 ALLOCATE(yvarname(inbvars))
621 haction=
'get variables names' 622 status=nf90_inquire_variable(kcdf_id,jloop1,
name=yvarname(jloop1))
625 if (yvarname(jloop1)==hncvarname)
then 629 if (yvarname(jloop1)/=hncvarname)
then 630 if((lgt(
trim(yvarname(jloop1)),
trim(hncvarname))).AND.&
631 (scan(
trim(yvarname(jloop1)),
trim(hncvarname))==1))
then 639 if (id_vartoget1/=0)
then 640 id_vartoget=id_vartoget1
642 id_vartoget=id_vartoget2
645 if (id_vartoget==0)
then 646 haction=
'close netcdf' 647 status=nf90_close(kcdf_id)
655 CALL mpi_bcast(id_vartoget,kind(id_vartoget)/4,mpi_integer,
npio,
ncomm,infompi)
659 if (id_vartoget==0)
then 660 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',1,zhook_handle)
676 haction=
'get variable dimensions number' 677 status=nf90_inquire_variable(kcdf_id,id_vartoget,ndims=invardims)
680 ALLOCATE(nvardimid(invardims))
681 haction=
'get variable dimensions identifiant' 682 status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=nvardimid)
687 SELECT CASE (invardims)
690 haction=
'get variable dimensions length' 691 status=nf90_inquire_dimension(kcdf_id,nvardimid(1),len=idim)
697 haction=
'get variable dimensions length' 698 status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop),len=indimlen(jloop))
700 haction=
'get variable dimensions names' 701 status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop),
name=ndimnam(jloop))
703 if ((ndimnam(jloop)==
'lat').OR.(ndimnam(jloop)==
'latitude'))
then 709 if ((ndimnam(jloop)==
'lon').OR.(ndimnam(jloop)==
'longitude'))
then 710 inlon=indimlen(jloop)
714 if (ndimnam(jloop)==
'depth')
nindepth=indimlen(jloop)
721 haction=
'close netcdf' 722 status=nf90_close(kcdf_id)
731 CALL mpi_bcast(inlon,kind(inlon)/4,mpi_integer,
npio,
ncomm,infompi)
761 IF (
ALLOCATED(
no))
DEALLOCATE(
no)
762 IF (
ALLOCATED(
xla))
DEALLOCATE(
xla)
763 IF (
ALLOCATED(
xola))
DEALLOCATE(
xola)
764 IF (
ALLOCATED(
xolo))
DEALLOCATE(
xolo)
778 IF (
ALLOCATED(
np))
DEALLOCATE(
np)
781 ALLOCATE(
xloph(ino,12))
783 IF (
lglobs) iinla = iinla + 2
784 IF (
lglobn) iinla = iinla + 2
790 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:PREP_NETCDF_GRID',1,zhook_handle)
803 CHARACTER(LEN=28),
INTENT(IN) :: HFILENAME
804 CHARACTER(LEN=28),
INTENT(IN) :: HNCVARNAME
805 REAL,
DIMENSION(:),
INTENT(OUT) :: PVAL
810 character(len=80) :: HACTION
811 character(len=80),
DIMENSION(:),
ALLOCATABLE :: VARNAME
813 integer ::ID_VARTOGET,ID_VARTOGET1,ID_VARTOGET2
815 INTEGER,
DIMENSION(1) :: NDIMID
817 real,
DIMENSION(:),
ALLOCATABLE :: ZVALU
819 REAL(KIND=JPRB) :: ZHOOK_HANDLE
824 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_Z1D_CDF',0,zhook_handle)
827 haction=
'open netcdf' 828 status=nf90_open(hfilename,nf90_nowrite,kcdf_id)
829 if (status/=nf90_noerr)
then 835 haction=
'get number of variables' 836 status=nf90_inquire(kcdf_id,nvariables=nbvars)
839 ALLOCATE(varname(nbvars))
846 haction=
'get variables names' 847 status=nf90_inquire_variable(kcdf_id,jloop1,
name=varname(jloop1))
849 if (varname(jloop1)==hncvarname)
then 852 if (varname(jloop1)/=hncvarname)
then 853 if((lgt(
trim(varname(jloop1)),
trim(hncvarname))).AND.&
854 (scan(
trim(varname(jloop1)),
trim(hncvarname))==1))
then 859 if (id_vartoget1/=0)
then 860 id_vartoget=id_vartoget1
862 id_vartoget=id_vartoget2
864 if (id_vartoget==0)
then 865 haction=
'close netcdf' 866 status=nf90_close(kcdf_id)
868 CALL abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_Z1D_CDF')
875 haction=
'get variable dimensions number' 876 status=nf90_inquire_variable(kcdf_id,id_vartoget,ndims=nvardims)
881 SELECT CASE (nvardims)
884 haction=
'get variable dimensions length' 885 status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=ndimid)
886 status=nf90_inquire_dimension(kcdf_id,ndimid(1),len=nlen)
888 ALLOCATE(zvalu(nlen))
890 CALL get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
894 write(0,*)
'YOU ARE TRYING TO READ A 2D FIELD FOR :',
trim(hncvarname)
895 CALL abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_Z1D_CDF')
900 haction=
'close netcdf' 901 status=nf90_close(kcdf_id)
906 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
907 IF (
ALLOCATED(zvalu ))
DEALLOCATE(zvalu )
909 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_Z1D_CDF',1,zhook_handle)
921 CHARACTER(LEN=28),
INTENT(IN) :: HFILENAME
922 CHARACTER(LEN=28),
INTENT(IN) :: HNCVARNAME
923 REAL,
DIMENSION(:),
INTENT(OUT) :: PLON,PLAT
924 REAL,
DIMENSION(:),
INTENT(OUT) :: PVAL
929 character(len=80) :: HACTION
930 character(len=80),
DIMENSION(:),
ALLOCATABLE :: VARNAME
931 integer ::JLOOP1,JDIM1,JDIM2,JLOOP
932 integer ::ID_VARTOGET,ID_VARTOGET1,ID_VARTOGET2
935 integer,
dimension(1) :: NDIMID
936 integer,
DIMENSION(2) ::NLEN2D, NDIMID2D
937 integer,
DIMENSION(:),
ALLOCATABLE :: NVARDIMID,NVARDIMLEN
938 character(len=80),
DIMENSION(:),
ALLOCATABLE :: NVARDIMNAM
939 real,
DIMENSION(:),
ALLOCATABLE :: ZVALU
940 real,
DIMENSION(:,:),
ALLOCATABLE :: ZVALU2D
942 real,
DIMENSION(:),
ALLOCATABLE :: ZDIM1
943 real,
DIMENSION(:),
ALLOCATABLE :: ZDIM2
944 character(len=80) :: YDIM1NAME,YDIM2NAME
945 integer :: ILONFOUND,ILATFOUND, IARG
946 REAL(KIND=JPRB) :: ZHOOK_HANDLE
953 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONVAL_CDF',0,zhook_handle)
956 haction=
'open netcdf' 957 status=nf90_open(hfilename,nf90_nowrite,kcdf_id)
960 if (status/=nf90_noerr)
then 970 haction=
'get number of variables' 971 status=nf90_inquire(kcdf_id,nvariables=nbvars)
974 ALLOCATE(varname(nbvars))
983 haction=
'get variables names' 984 status=nf90_inquire_variable(kcdf_id,jloop1,
name=varname(jloop1))
987 if (varname(jloop1)==hncvarname)
then 991 if (varname(jloop1)/=hncvarname)
then 992 if((lgt(
trim(varname(jloop1)),
trim(hncvarname))).AND.&
993 (scan(
trim(varname(jloop1)),
trim(hncvarname))==1))
then 1000 if (id_vartoget1/=0)
then 1001 id_vartoget=id_vartoget1
1003 id_vartoget=id_vartoget2
1005 if (id_vartoget==0)
then 1006 haction=
'close netcdf' 1007 status=nf90_close(kcdf_id)
1009 CALL abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_LATLONVAL_CDF')
1019 haction=
'get variable dimensions number' 1020 status=nf90_inquire_variable(kcdf_id,id_vartoget,ndims=nvardims)
1026 SELECT CASE (nvardims)
1029 haction=
'get variable dimensions length' 1030 status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=ndimid)
1031 status=nf90_inquire_dimension(kcdf_id,ndimid(1),len=nlen)
1033 ALLOCATE(zvalu(nlen))
1035 CALL get1dcdf(kcdf_id,id_vartoget,zmiss,zvalu)
1039 status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=ndimid2d)
1041 haction=
'get variable dimensions length' 1042 status=nf90_inquire_dimension(kcdf_id,ndimid2d(jloop),len=nlen2d(jloop))
1045 ALLOCATE(zvalu2d(nlen2d(1),nlen2d(2)))
1046 ALLOCATE(zdim1(nlen2d(1)))
1047 ALLOCATE(zdim2(nlen2d(2)))
1049 CALL get2dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,ydim2name,&
1053 if ((ydim1name==
'lon').OR.(ydim1name==
'longitude')) ilonfound=1
1054 if ((ydim2name==
'lon').OR.(ydim2name==
'longitude')) ilonfound=2
1055 if ((ydim1name==
'lat').OR.(ydim1name==
'latitude')) ilatfound=1
1056 if ((ydim2name==
'lat').OR.(ydim2name==
'latitude')) ilatfound=2
1061 IF ((ilonfound==1).AND.(ilatfound==2))
then 1063 DO jdim1=1,
SIZE(zdim1)
1064 DO jdim2=1,
SIZE(zdim2)
1066 pval(iarg)=zvalu2d(jdim1,jdim2)
1067 plon(iarg)=zdim1(jdim1)
1068 plat(iarg)=zdim2(jdim2)
1071 ELSEIF ((ilonfound==2).AND.(ilatfound==1))
then 1073 DO jdim1=1,
SIZE(zdim1)
1074 DO jdim2=1,
SIZE(zdim2)
1076 pval(iarg)=zvalu2d(jdim1,jdim2)
1077 plat(iarg)=zdim1(jdim1)
1078 plon(iarg)=zdim2(jdim2)
1082 write(0,*)
'*****WARNING*****: incompatible dimensions to lat/lon/value arrays' 1091 haction=
'close netcdf' 1092 status=nf90_close(kcdf_id)
1099 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
1100 IF (
ALLOCATED(zvalu ))
DEALLOCATE(zvalu )
1101 IF (
ALLOCATED(zvalu2d ))
DEALLOCATE(zvalu2d)
1102 IF (
ALLOCATED(zdim1 ))
DEALLOCATE(zdim1 )
1103 IF (
ALLOCATED(zdim2 ))
DEALLOCATE(zdim2 )
1106 IF (
ALLOCATED(nvardimid ))
DEALLOCATE(nvardimid )
1107 IF (
ALLOCATED(nvardimnam ))
DEALLOCATE(nvardimnam)
1108 IF (
ALLOCATED(nvardimlen ))
DEALLOCATE(nvardimlen)
1109 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONVAL_CDF',1,zhook_handle)
1121 CHARACTER(LEN=28),
INTENT(IN) :: HFILENAME
1122 CHARACTER(LEN=28),
INTENT(IN) :: HNCVARNAME
1123 REAL,
DIMENSION(:),
INTENT(OUT) :: PLON,PLAT
1124 REAL,
DIMENSION(:),
INTENT(OUT) :: PDEP
1125 REAL,
DIMENSION(:,:),
INTENT(OUT) :: PVAL
1130 character(len=80) :: HACTION
1131 character(len=80),
DIMENSION(:),
ALLOCATABLE :: VARNAME
1132 integer ::JLOOP1,JDIM1,JDIM2,JDIM3,JLOOP
1134 integer ::ID_VARTOGET,ID_VARTOGET1,ID_VARTOGET2
1136 integer,
DIMENSION(3) ::NLEN3D
1137 integer,
DIMENSION(:),
ALLOCATABLE :: NVARDIMID,NVARDIMLEN
1138 character(len=80),
DIMENSION(:),
ALLOCATABLE :: NVARDIMNAM
1139 real,
DIMENSION(:,:,:),
ALLOCATABLE :: ZVALU3D
1141 real,
DIMENSION(:),
ALLOCATABLE :: ZDIM1
1142 real,
DIMENSION(:),
ALLOCATABLE :: ZDIM2
1143 real,
DIMENSION(:),
ALLOCATABLE :: ZDIM3
1144 character(len=80) :: YDIM1NAME,YDIM2NAME,YDIM3NAME
1145 integer :: ILONFOUND,ILATFOUND,IDEPFOUND
1147 REAL(KIND=JPRB) :: ZHOOK_HANDLE
1154 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONDEPVAL_CDF',0,zhook_handle)
1155 haction=
'open netcdf' 1156 status=nf90_open(hfilename,nf90_nowrite,kcdf_id)
1164 haction=
'get number of variables' 1165 status=nf90_inquire(kcdf_id,nvariables=nbvars)
1168 ALLOCATE(varname(nbvars))
1177 haction=
'get variables names' 1178 status=nf90_inquire_variable(kcdf_id,jloop1,
name=varname(jloop1))
1181 if (varname(jloop1)==hncvarname)
then 1185 if (varname(jloop1)/=hncvarname)
then 1186 if((lgt(
trim(varname(jloop1)),
trim(hncvarname))).AND.&
1187 (scan(
trim(varname(jloop1)),
trim(hncvarname))==1))
then 1194 if (id_vartoget1/=0)
then 1195 id_vartoget=id_vartoget1
1197 id_vartoget=id_vartoget2
1199 if (id_vartoget==0)
then 1200 haction=
'close netcdf' 1201 status=nf90_close(kcdf_id)
1203 CALL abor1_sfx(
'MODE_READ_NETCDF_MERCATOR: READ_LATLONDEPVAL_CDF')
1213 haction=
'get variable dimensions number' 1214 status=nf90_inquire_variable(kcdf_id,id_vartoget,ndims=nvardims)
1217 ALLOCATE(nvardimid(nvardims))
1218 haction=
'get variable dimensions identifiant' 1219 status=nf90_inquire_variable(kcdf_id,id_vartoget,dimids=nvardimid)
1225 SELECT CASE (nvardims)
1228 write(0,*)
'********************************************' 1229 write(0,*)
'* number of dimension to low: ',nvardims,
' *' 1230 write(0,*)
'* you need a 3-dimension variable *' 1231 write(0,*)
'********************************************' 1235 haction=
'get variable dimensions length' 1236 status=nf90_inquire_dimension(kcdf_id,nvardimid(jloop),len=nlen3d(jloop))
1239 ALLOCATE(zvalu3d(nlen3d(1),nlen3d(2),nlen3d(3)))
1240 ALLOCATE(zdim1(nlen3d(1)))
1241 ALLOCATE(zdim2(nlen3d(2)))
1242 ALLOCATE(zdim3(nlen3d(3)))
1244 CALL get3dcdf(kcdf_id,id_vartoget,zdim1,ydim1name,zdim2,ydim2name,&
1245 zdim3,ydim3name,zmiss,zvalu3d)
1249 if ((ydim1name==
'lon').OR.(ydim1name==
'longitude')) ilonfound=1
1250 if ((ydim2name==
'lon').OR.(ydim2name==
'longitude')) ilonfound=2
1251 if ((ydim3name==
'lon').OR.(ydim3name==
'longitude')) ilonfound=3
1252 if ((ydim1name==
'lat').OR.(ydim1name==
'latitude')) ilatfound=1
1253 if ((ydim2name==
'lat').OR.(ydim2name==
'latitude')) ilatfound=2
1254 if ((ydim3name==
'lat').OR.(ydim3name==
'latitude')) ilatfound=3
1255 if (ydim1name==
'depth') idepfound=1
1256 if (ydim2name==
'depth') idepfound=2
1257 if (ydim3name==
'depth') idepfound=3
1263 IF ((ilonfound==1).AND.(ilatfound==2).AND.(idepfound==3))
then 1267 DO jdim2=1,
SIZE(zdim2)
1268 DO jdim1=1,
SIZE(zdim1)
1270 plon(iarg)=zdim1(jdim1)
1271 plat(iarg)=zdim2(jdim2)
1272 DO jdim3=1,
SIZE(zdim3)
1273 pval(iarg,jdim3)=zvalu3d(jdim1,jdim2,jdim3)
1279 ELSEIF ((ilonfound==2).AND.(ilatfound==1).AND.(idepfound==3))
then 1281 DO jdim1=1,
SIZE(zdim1)
1282 DO jdim2=1,
SIZE(zdim2)
1284 plon(iarg)=zdim2(jdim2)
1285 plat(iarg)=zdim1(jdim1)
1286 DO jdim3=1,
SIZE(zdim3)
1287 pval(iarg,jdim3)=zvalu3d(jdim1,jdim2,jdim3)
1292 ELSEIF ((ilonfound==1).AND.(ilatfound==3).AND.(idepfound==2))
then 1294 DO jdim3=1,
SIZE(zdim3)
1295 DO jdim1=1,
SIZE(zdim1)
1297 plon(iarg)=zdim1(jdim1)
1298 plat(iarg)=zdim3(jdim3)
1299 DO jdim2=1,
SIZE(zdim2)
1300 pval(iarg,jdim2)=zvalu3d(jdim1,jdim2,jdim3)
1305 ELSEIF ((ilatfound==1).AND.(ilonfound==3).AND.(idepfound==2))
then 1307 DO jdim1=1,
SIZE(zdim1)
1308 DO jdim3=1,
SIZE(zdim3)
1310 plon(iarg)=zdim3(jdim3)
1311 plat(iarg)=zdim1(jdim1)
1312 DO jdim2=1,
SIZE(zdim2)
1313 pval(iarg,jdim2)=zvalu3d(jdim1,jdim2,jdim3)
1318 ELSEIF ((ilonfound==2).AND.(ilatfound==3).AND.(idepfound==1))
then 1320 DO jdim3=1,
SIZE(zdim3)
1321 DO jdim2=1,
SIZE(zdim2)
1323 plon(iarg)=zdim2(jdim2)
1324 plat(iarg)=zdim3(jdim3)
1325 DO jdim1=1,
SIZE(zdim1)
1326 pval(iarg,jdim1)=zvalu3d(jdim1,jdim2,jdim3)
1331 ELSEIF ((ilatfound==2).AND.(ilonfound==3).AND.(idepfound==1))
then 1333 DO jdim2=1,
SIZE(zdim2)
1334 DO jdim3=1,
SIZE(zdim3)
1336 plon(iarg)=zdim3(jdim3)
1337 plat(iarg)=zdim2(jdim2)
1338 DO jdim1=1,
SIZE(zdim1)
1339 pval(iarg,jdim1)=zvalu3d(jdim1,jdim2,jdim3)
1345 write(0,*)
'*****WARNING*****: incompatible dimensions to lat/lon/value arrays' 1353 haction=
'close netcdf' 1355 status=nf90_close(kcdf_id)
1362 IF (
ALLOCATED(varname ))
DEALLOCATE(varname)
1363 IF (
ALLOCATED(zvalu3d ))
DEALLOCATE(zvalu3d)
1364 IF (
ALLOCATED(zdim1 ))
DEALLOCATE(zdim1 )
1365 IF (
ALLOCATED(zdim2 ))
DEALLOCATE(zdim2 )
1366 IF (
ALLOCATED(zdim3 ))
DEALLOCATE(zdim3 )
1369 IF (
ALLOCATED(nvardimid ))
DEALLOCATE(nvardimid )
1370 IF (
ALLOCATED(nvardimnam ))
DEALLOCATE(nvardimnam)
1371 IF (
ALLOCATED(nvardimlen ))
DEALLOCATE(nvardimlen)
1372 201
FORMAT(4(3x,f10.4))
1373 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_LATLONDEPVAL_CDF',1,zhook_handle)
1390 CHARACTER(LEN=28),
INTENT(IN) :: HFILENAME
1391 CHARACTER(LEN=28),
INTENT(IN) :: HNCVARNAME
1392 REAL,
POINTER,
DIMENSION(:) :: PFIELD
1394 REAL,
DIMENSION(:),
ALLOCATABLE :: ZLATI
1395 REAL,
DIMENSION(:),
ALLOCATABLE :: ZLONG
1396 REAL,
TARGET,
DIMENSION(:,:),
ALLOCATABLE :: ZVALUE
1397 REAL,
DIMENSION(:),
ALLOCATABLE :: ZDEPTH
1398 REAL,
TARGET,
DIMENSION(:),
ALLOCATABLE :: ZVAL
1402 REAL(KIND=JPRB) :: ZHOOK_HANDLE
1406 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_SST',0,zhook_handle)
1419 IF (
trim(hncvarname) ==
'temperature')
THEN 1420 WHERE (zval(:)/=zundef .AND. zval(:)<100.) pfield(:)=pfield(:)+
xtt 1435 pfield(:)=zvalue(:,1)
1436 IF (
trim(hncvarname) ==
'temperature')
THEN 1437 WHERE (zvalue(:,1)/=zundef .AND. zvalue(:,1)<100.) pfield(:)=pfield(:)+
xtt 1444 IF (
ALLOCATED(zvalue ))
DEALLOCATE(zvalue )
1445 IF (
ALLOCATED(zlong ))
DEALLOCATE(zlong )
1446 IF (
ALLOCATED(zlati ))
DEALLOCATE(zlati )
1447 IF (
ALLOCATED(zdepth ))
DEALLOCATE(zdepth )
1448 IF (
ALLOCATED(zval ))
DEALLOCATE(zval )
1450 202
FORMAT(3(3x,f10.4))
1451 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_SST',1,zhook_handle)
1467 CHARACTER(LEN=28),
INTENT(IN) :: HFILENAME
1468 CHARACTER(LEN=28),
INTENT(IN) :: HNCVARNAME
1469 REAL,
POINTER,
DIMENSION(:) :: PFIELD
1471 REAL,
DIMENSION(:),
ALLOCATABLE :: ZLATI
1472 REAL,
DIMENSION(:),
ALLOCATABLE :: ZLONG
1473 REAL,
TARGET,
DIMENSION(:),
ALLOCATABLE:: ZVALUE
1475 REAL(KIND=JPRB) :: ZHOOK_HANDLE
1479 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_ZS_SEA',0,zhook_handle)
1502 IF (
ALLOCATED(zvalue ))
DEALLOCATE(zvalue )
1503 IF (
ALLOCATED(zlong ))
DEALLOCATE(zlong )
1504 IF (
ALLOCATED(zlati ))
DEALLOCATE(zlati )
1506 202
FORMAT(3(3x,f10.4))
1507 IF (
lhook)
CALL dr_hook(
'MODE_READ_NETCDF_MERCATOR:READ_NETCDF_ZS_SEA',1,zhook_handle)
subroutine horibl_surf_coef(KOLEN, OINTERP, OGLOBLON, PILO1, PILO2, POLO, KO, KINLO, KP, PLOP)
integer, dimension(:), allocatable ninloh
static const char * trim(const char *name, int *n)
real, dimension(:), allocatable xilatarray
subroutine read_z1d_cdf(HFILENAME, HNCVARNAME, PVAL)
subroutine read_latlondepval_cdf(HFILENAME, HNCVARNAME, PLON, PLAT, PDEP, PVAL)
real, dimension(:), allocatable xola
subroutine get3dcdf(KCDF_ID, IDVAR, PDIM1, HDIM1NAME, PDIM2, HDIM2NAME, PDIM3, HDIM3NAME, PMISSVALUE, PVALU3D)
integer, dimension(:,:), allocatable np
real, dimension(:), allocatable xlon_out
subroutine prep_netcdf_grid(HFILENAME, HNCVARNAME)
quick &counting sorts only inumt inumt name
character(len=6) cinterp_type
subroutine get1dcdf(KCDF_ID, IDVAR, PMISSVALUE, PVALU1D)
subroutine get2dcdf(KCDF_ID, IDVAR, PDIM1, HDIM1NAME, PDIM2, HDIM2NAME, PMISSVALUE, PVALU2D)
subroutine handle_err_mer(status, line)
integer, dimension(:,:), allocatable no
subroutine abor1_sfx(YTEXT)
logical, dimension(:), allocatable linterp
real, dimension(:), allocatable xilonarray
real, dimension(:), allocatable xlat_out
integer, parameter nundef
subroutine read_netcdf_zs_sea(HFILENAME, HNCVARNAME, PFIELD)
real, dimension(:), allocatable xolo
subroutine read_dim_cdf(HFILENAME, HNCVARNAME, KDIM)
subroutine horibl_surf_init(PILA1, PILO1, PILA2, PILO2, KINLA, KINLO, KOLEN, PXOUT, PYOUT, OINTERP, OGLOBLON, OGLOBN, OGLOBS, KO, KINLO_OUT, POLA, POLO, PILO1_OUT, PILO2_OUT, PLA, PILATARRAY)
integer, dimension(:), allocatable ninlon
subroutine read_latlonval_cdf(HFILENAME, HNCVARNAME, PLON, PLAT, PVAL)
real, dimension(:,:), allocatable xla
real, dimension(:,:), allocatable xloph
subroutine wlog_mpi(HLOG, PLOG, KLOG, KLOG2, OLOG)
subroutine read_netcdf_sst(HFILENAME, HNCVARNAME, PFIELD)