56 SUBROUTINE ncopen(KLISTING,ORW,OVERBOSE,HFILENAME,KNCID)
68 CHARACTER(LEN=NF90_MAX_NAME),
INTENT(IN) :: HFILENAME
70 LOGICAL,
INTENT(IN) :: ORW, OVERBOSE
72 INTEGER,
INTENT(IN) :: KLISTING
74 INTEGER,
INTENT(OUT) :: KNCID
78 CHARACTER(LEN=NF90_MAX_NAME) :: YFNAME
81 REAL(KIND=JPRB) :: ZHOOK_HANDLE
85 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCOPEN',0,zhook_handle)
87 yfname = hfilename(1:len_trim(hfilename))
90 ic = nf90_open(yfname,nf90_write,kncid)
92 ic = nf90_open(yfname,nf90_nowrite,kncid)
95 IF(ic/=nf90_noerr)
THEN 96 WRITE(klisting,*)
'NCOPEN for TRIP : Error opening file ',hfilename(1:len_trim(hfilename))
97 WRITE(klisting,*)nf90_strerror(ic)
98 CALL abort_trip(
'NCOPEN for TRIP : Error opening file')
100 WRITE(klisting,*)
'NCOPEN for TRIP : Opening file ',hfilename(1:len_trim(hfilename))
103 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCOPEN',1,zhook_handle)
110 SUBROUTINE ncclose(KLISTING,OVERBOSE,HFILENAME,KNCID)
122 CHARACTER(LEN=NF90_MAX_NAME),
INTENT(IN) :: HFILENAME
124 LOGICAL,
INTENT(IN) :: OVERBOSE
125 INTEGER,
INTENT(IN) :: KLISTING
126 INTEGER,
INTENT(IN) :: KNCID
131 REAL(KIND=JPRB) :: ZHOOK_HANDLE
135 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCCLOSE',0,zhook_handle)
137 ic = nf90_close(kncid)
138 IF(overbose)
WRITE(klisting,*)
'NCCLOSE for TRIP : Close file ',hfilename(1:len_trim(hfilename))
140 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCCLOSE',1,zhook_handle)
147 SUBROUTINE ncread_x(KLISTING,KNCID,HVNAME,PVECT,OVERBOSE)
159 CHARACTER(LEN=NF90_MAX_NAME),
INTENT(IN) :: HVNAME
161 INTEGER,
INTENT(IN) :: KLISTING
162 INTEGER,
INTENT(IN) :: KNCID
164 LOGICAL,
INTENT(IN) :: OVERBOSE
166 REAL,
DIMENSION(:),
INTENT(OUT) :: PVECT
170 CHARACTER(LEN=NF90_MAX_NAME) :: YVNAME
172 REAL*4,
DIMENSION(SIZE(PVECT)) :: ZVECT
173 INTEGER,
DIMENSION(SIZE(PVECT)) :: IVECT
178 INTEGER :: IC, ID, IVAR_TYPE
180 REAL(KIND=JPRB) :: ZHOOK_HANDLE
184 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCREAD_X',0,zhook_handle)
186 yvname = hvname(1:len_trim(hvname))
188 ic = nf90_inq_varid(kncid,yvname,id)
189 IF(ic/=nf90_noerr)
THEN 190 WRITE(klisting,*)
'NCREAD_X for TRIP : Error reading id for variable ',hvname(1:len_trim(hvname))
191 WRITE(klisting,*)nf90_strerror(ic)
192 CALL abort_trip(
'NCREAD_X for TRIP : Error reading id variable ')
195 ic = nf90_inquire_variable(kncid,id,xtype=ivar_type)
196 IF(ic/=nf90_noerr)
THEN 197 WRITE(klisting,*)
'NCREAD_X for TRIP : Error reading type for variable ',hvname(1:len_trim(hvname))
198 WRITE(klisting,*)nf90_strerror(ic)
199 CALL abort_trip(
'NCREAD_X for TRIP : Error reading type variable ')
202 IF(ivar_type==nf90_double)
THEN 203 ic=nf90_get_att(kncid,id,
'missing_value',zmiss)
204 IF(ic/=nf90_noerr) ic=nf90_get_att(kncid,id,
'_FillValue',zmiss)
205 ic=nf90_get_var(kncid,id,pvect)
206 WHERE(pvect(:)==zmiss)
209 ELSEIF(ivar_type==nf90_float)
THEN 210 ic=nf90_get_att(kncid,id,
'missing_value',zmiss4)
211 IF(ic/=nf90_noerr) ic=nf90_get_att(kncid,id,
'_FillValue',zmiss4)
212 ic=nf90_get_var(kncid,id,zvect)
213 WHERE(zvect(:)/=zmiss4)
218 ELSEIF(ivar_type==nf90_int)
THEN 219 ic=nf90_get_var(kncid,id,ivect)
220 pvect(:)=
REAL(ivect(:))
222 WRITE(klisting,*)
'NCREAD_X for TRIP : type of the variable must be integer, float or double ',hvname(1:len_trim(hvname))
223 CALL abort_trip(
'NCREAD_X for TRIP : type of the variable not good')
225 IF(ic/=nf90_noerr)
THEN 226 WRITE(klisting,*)
'NCREAD_X for TRIP : Error reading variable ',hvname(1:len_trim(hvname))
227 WRITE(klisting,*)nf90_strerror(ic)
228 CALL abort_trip(
'NCREAD_X for TRIP : Error reading variable')
230 WRITE(klisting,*)
'NCREAD_X for TRIP : Success in reading variable ',hvname(1:len_trim(hvname))
233 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCREAD_X',1,zhook_handle)
239 SUBROUTINE ncread_xy(KLISTING,KNCID,HVNAME,PVECT,OVERBOSE)
251 CHARACTER(LEN=NF90_MAX_NAME),
INTENT(IN) :: HVNAME
253 INTEGER,
INTENT(IN) :: KLISTING
254 INTEGER,
INTENT(IN) :: KNCID
256 LOGICAL,
INTENT(IN) :: OVERBOSE
258 REAL,
DIMENSION(:,:),
INTENT(OUT) :: PVECT
262 CHARACTER(LEN=NF90_MAX_NAME) :: YVNAME
264 REAL*4,
DIMENSION(SIZE(PVECT,1),SIZE(PVECT,2)) :: ZVECT
269 INTEGER :: IC, ID, IVAR_TYPE
271 REAL(KIND=JPRB) :: ZHOOK_HANDLE
275 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCREAD_XY',0,zhook_handle)
277 yvname = hvname(1:len_trim(hvname))
279 ic = nf90_inq_varid(kncid,yvname,id)
280 IF(ic/=nf90_noerr)
THEN 281 WRITE(klisting,*)
'NCREAD_XY for TRIP : Error reading id for variable ',hvname(1:len_trim(hvname))
282 WRITE(klisting,*)nf90_strerror(ic)
283 CALL abort_trip(
'NCREAD_XY for TRIP : Error reading id variable ')
286 ic = nf90_inquire_variable(kncid,id,xtype=ivar_type)
287 IF(ic/=nf90_noerr)
THEN 288 WRITE(klisting,*)
'NCREAD_XY for TRIP : Error reading type for variable ',hvname(1:len_trim(hvname))
289 WRITE(klisting,*)nf90_strerror(ic)
290 CALL abort_trip(
'NCREAD_XY for TRIP : Error reading type variable ')
293 IF(ivar_type==nf90_double)
THEN 294 ic=nf90_get_att(kncid,id,
'missing_value',zmiss)
295 IF(ic/=nf90_noerr) ic=nf90_get_att(kncid,id,
'_FillValue',zmiss)
296 ic=nf90_get_var(kncid,id,pvect)
297 WHERE(pvect(:,:)==zmiss)
300 ELSEIF(ivar_type==nf90_float)
THEN 301 ic=nf90_get_att(kncid,id,
'missing_value',zmiss4)
302 IF(ic/=nf90_noerr) ic=nf90_get_att(kncid,id,
'_FillValue',zmiss4)
303 ic=nf90_get_var(kncid,id,zvect)
304 WHERE(zvect(:,:)/=zmiss4)
305 pvect(:,:)=zvect(:,:)
310 WRITE(klisting,*)
'NCREAD_XY for TRIP : type of the variable must be float or double ',hvname(1:len_trim(hvname))
311 CALL abort_trip(
'NCREAD_XY for TRIP : type of the variable not good')
314 IF(ic/=nf90_noerr)
THEN 315 WRITE(klisting,*)
'NCREAD_XY for TRIP : Error reading variable ',hvname(1:len_trim(hvname))
316 WRITE(klisting,*)nf90_strerror(ic)
317 CALL abort_trip(
'NCREAD_XY for TRIP : Error reading variable')
319 WRITE(klisting,*)
'NCREAD_XY for TRIP : Success in reading variable ',hvname(1:len_trim(hvname))
322 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCREAD_XY',1,zhook_handle)
328 SUBROUTINE ncread_xyz(KLISTING,KNCID,HVNAME,PVECT,OVERBOSE)
340 CHARACTER(LEN=NF90_MAX_NAME),
INTENT(IN) :: HVNAME
342 INTEGER,
INTENT(IN) :: KLISTING
343 INTEGER,
INTENT(IN) :: KNCID
345 LOGICAL,
INTENT(IN) :: OVERBOSE
347 REAL,
DIMENSION(:,:,:),
INTENT(OUT) :: PVECT
351 CHARACTER(LEN=NF90_MAX_NAME) :: YVNAME
353 REAL*4,
DIMENSION(SIZE(PVECT,1),SIZE(PVECT,2),SIZE(PVECT,3)) :: ZVECT
358 INTEGER :: IC, ID, IVAR_TYPE
360 REAL(KIND=JPRB) :: ZHOOK_HANDLE
364 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCREAD_XYZ',0,zhook_handle)
365 yvname = hvname(1:len_trim(hvname))
367 ic = nf90_inq_varid(kncid,yvname,id)
369 ic = nf90_inq_varid(kncid,yvname,id)
370 IF(ic/=nf90_noerr)
THEN 371 WRITE(klisting,*)
'NCREAD_XYZ for TRIP : Error reading id for variable ',hvname(1:len_trim(hvname))
372 WRITE(klisting,*)nf90_strerror(ic)
373 CALL abort_trip(
'NCREAD_XYZ for TRIP : Error reading id variable ')
376 ic = nf90_inquire_variable(kncid,id,xtype=ivar_type)
377 IF(ic/=nf90_noerr)
THEN 378 WRITE(klisting,*)
'NCREAD_XYZ for TRIP : Error reading type for variable ',hvname(1:len_trim(hvname))
379 WRITE(klisting,*)nf90_strerror(ic)
380 CALL abort_trip(
'NCREAD_XYZ for TRIP : Error reading type variable ')
383 IF(ivar_type==nf90_double)
THEN 384 ic=nf90_get_att(kncid,id,
'missing_value',zmiss)
385 IF(ic/=nf90_noerr) ic=nf90_get_att(kncid,id,
'_FillValue',zmiss)
386 ic=nf90_get_var(kncid,id,pvect)
387 WHERE(pvect(:,:,:)==zmiss)
390 ELSEIF(ivar_type==nf90_float)
THEN 391 ic=nf90_get_att(kncid,id,
'missing_value',zmiss4)
392 IF(ic/=nf90_noerr) ic=nf90_get_att(kncid,id,
'_FillValue',zmiss4)
393 ic=nf90_get_var(kncid,id,zvect)
394 WHERE(zvect(:,:,:)/=zmiss4)
395 pvect(:,:,:)=zvect(:,:,:)
400 WRITE(klisting,*)
'NCREAD_XYZ for TRIP : type of the variable must be float or double ',hvname(1:len_trim(hvname))
401 CALL abort_trip(
'NCREAD_XYZ for TRIP : type of the variable not good')
403 IF(ic/=nf90_noerr)
THEN 404 WRITE(klisting,*)
'NCREAD_XYZ for TRIP : Error reading variable ',hvname(1:len_trim(hvname))
405 WRITE(klisting,*)nf90_strerror(ic)
406 CALL abort_trip(
'NCREAD_XYZ for TRIP : Error reading variable')
408 WRITE(klisting,*)
'NCREAD_XYZ for TRIP : Success in reading variable ',hvname(1:len_trim(hvname))
411 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCREAD_XYZ',1,zhook_handle)
418 SUBROUTINE nccreate(KLISTING,HFILENAME,HTITLE,HTIMEUNIT, &
419 HVNAME,HVLNAME,HUNIT,PLON,PLAT, &
420 PMISSVAL,OVERBOSE,KNCID,OTIME, &
421 KZLEN,OVARZDIM,ODOUBLE)
433 CHARACTER(LEN=NF90_MAX_NAME),
INTENT(IN) :: HFILENAME, HTITLE, HTIMEUNIT
435 CHARACTER(LEN=NF90_MAX_NAME),
DIMENSION(:),
INTENT(IN) :: HVNAME, HVLNAME, HUNIT
437 REAL,
DIMENSION(:),
INTENT(IN) :: PLON
438 REAL,
DIMENSION(:),
INTENT(IN) :: PLAT
440 LOGICAL,
INTENT(IN) :: OVERBOSE
442 REAL,
INTENT(IN) :: PMISSVAL
444 INTEGER,
INTENT(IN) :: KLISTING
446 INTEGER,
INTENT(OUT) :: KNCID
448 LOGICAL,
INTENT(IN) :: OTIME
450 INTEGER,
INTENT(IN),
OPTIONAL :: KZLEN
451 LOGICAL,
DIMENSION(:),
INTENT(IN),
OPTIONAL :: OVARZDIM
452 LOGICAL,
DIMENSION(:),
INTENT(IN),
OPTIONAL :: ODOUBLE
456 CHARACTER(LEN=NF90_MAX_NAME) :: YWORK
458 REAL*4,
DIMENSION(:),
ALLOCATABLE :: ZWORK
460 REAL*4,
DIMENSION(SIZE(PLON)) :: ZLON
461 REAL*4,
DIMENSION(SIZE(PLAT)) :: ZLAT
463 LOGICAL,
DIMENSION(SIZE(HVNAME)) :: LDOUBLE
467 INTEGER :: ILONDIM, ILATDIM, ILEVDIM, ITIMEDIM
468 INTEGER :: ILON_ID, ILAT_ID, ILEV_ID, ITIME_ID, VAR_ID
469 INTEGER :: IC, IWORK, INVAR
472 REAL(KIND=JPRB) :: ZHOOK_HANDLE
476 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCCREATE',0,zhook_handle)
478 IF(
PRESENT(odouble))
THEN 479 ldouble(:)=odouble(:)
486 ywork = hfilename(1:len_trim(hfilename))
488 ic = nf90_create(ywork,nf90_64bit_offset,kncid)
489 IF(ic/=nf90_noerr)
THEN 490 WRITE(klisting,*)
'NCCREATE for TRIP : Error create file ',hfilename(1:len_trim(hfilename))
491 WRITE(klisting,*)nf90_strerror(ic)
492 CALL abort_trip(
'NCCREATE for TRIP : Error create file')
494 WRITE(klisting,*)
'NCCREATE for TRIP : Success in creating file ',hfilename(1:len_trim(hfilename))
498 ywork = htitle(1:len_trim(htitle))
499 ic = nf90_put_att(kncid,nf90_global,
'title',ywork)
501 ic = nf90_put_att(kncid,nf90_global,
'Conventions',ywork)
505 ic = nf90_def_dim(kncid,
'longitude',iwork,ilondim)
507 ic = nf90_def_dim(kncid,
'latitude',iwork,ilatdim)
508 IF(
PRESENT(kzlen)) ic = nf90_def_dim(kncid,
'level',kzlen,ilevdim)
509 IF(otime) ic = nf90_def_dim(kncid,
'time',nf90_unlimited,itimedim)
513 ic = nf90_def_var(kncid,
'longitude',nf90_float,ilondim,ilon_id)
514 ic = nf90_def_var(kncid,
'latitude' ,nf90_float,ilatdim,ilat_id)
515 ywork =
'degrees_east' 516 ic = nf90_put_att(kncid,ilon_id,
'units',ywork)
517 ywork =
'degrees_north' 518 ic = nf90_put_att(kncid,ilat_id,
'units',ywork)
519 IF(
PRESENT(kzlen))
THEN 520 ic = nf90_def_var(kncid,
'level',nf90_float,ilevdim,ilev_id)
522 ic = nf90_put_att(kncid,ilev_id,
'units',ywork)
526 ywork = htimeunit(1:len_trim(htimeunit))
527 ic = nf90_def_var(kncid,
'time',nf90_int,itimedim,itime_id)
528 ic = nf90_put_att(kncid,itime_id,
'units',ywork)
533 zmissval =
REAL(pmissval)
538 ywork = hvname(jvar)(1:len_trim(hvname(jvar)))
539 IF(
PRESENT(kzlen))
THEN 541 IF(ovarzdim(jvar))
THEN 542 IF(ldouble(jvar))
THEN 543 ic = nf90_def_var(kncid,ywork,nf90_double,(/ilondim,ilatdim,ilevdim,itimedim/),var_id)
545 ic = nf90_def_var(kncid,ywork,nf90_float,(/ilondim,ilatdim,ilevdim,itimedim/),var_id)
548 IF(ldouble(jvar))
THEN 549 ic = nf90_def_var(kncid,ywork,nf90_double,(/ilondim,ilatdim,itimedim/),var_id)
551 ic = nf90_def_var(kncid,ywork,nf90_float,(/ilondim,ilatdim,itimedim/),var_id)
555 IF(ovarzdim(jvar))
THEN 556 IF(ldouble(jvar))
THEN 557 ic = nf90_def_var(kncid,ywork,nf90_double,(/ilondim,ilatdim,ilevdim/),var_id)
559 ic = nf90_def_var(kncid,ywork,nf90_float,(/ilondim,ilatdim,ilevdim/),var_id)
562 IF(ldouble(jvar))
THEN 563 ic = nf90_def_var(kncid,ywork,nf90_double,(/ilondim,ilatdim/),var_id)
565 ic = nf90_def_var(kncid,ywork,nf90_float,(/ilondim,ilatdim/),var_id)
571 IF(ldouble(jvar))
THEN 572 ic = nf90_def_var(kncid,ywork,nf90_double,(/ilondim,ilatdim,itimedim/),var_id)
574 ic = nf90_def_var(kncid,ywork,nf90_float,(/ilondim,ilatdim,itimedim/),var_id)
577 IF(ldouble(jvar))
THEN 578 ic = nf90_def_var(kncid,ywork,nf90_double,(/ilondim,ilatdim/),var_id)
580 ic = nf90_def_var(kncid,ywork,nf90_float,(/ilondim,ilatdim/),var_id)
585 ywork = hvlname(jvar)(1:len_trim(hvlname(jvar)))
586 ic = nf90_put_att(kncid,var_id,
'long_name',ywork)
587 ywork = hunit(jvar)(1:len_trim(hunit(jvar)))
588 ic = nf90_put_att(kncid,var_id,
'units',ywork)
589 IF(ldouble(jvar))
THEN 590 ic = nf90_put_att(kncid,var_id,
'_FillValue',pmissval)
592 ic = nf90_put_att(kncid,var_id,
'_FillValue',zmissval)
597 ic = nf90_enddef(kncid)
598 IF(ic/=nf90_noerr)
THEN 599 WRITE(klisting,*)
'NCCREATE for TRIP : Error file definition' 600 WRITE(klisting,*)nf90_strerror(ic)
601 CALL abort_trip(
'NCCREATE for TRIP : Error file definition')
603 WRITE(klisting,*)
'NCCREATE for TRIP : file definition is ok' 611 ic = nf90_put_var(kncid,ilon_id,zlon)
612 ic = nf90_put_var(kncid,ilat_id,zlat)
613 IF(
PRESENT(kzlen))
THEN 614 ALLOCATE(zwork(kzlen))
618 ic = nf90_put_var(kncid,ilev_id,zwork)
623 WRITE(klisting,*)
'NCCREATE ',hfilename(1:len_trim(hfilename)),
' for TRIP OK !' 626 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCCREATE',1,zhook_handle)
633 SUBROUTINE ncstore(KLISTING,KNCID,HVNAME,PWRITE, &
634 OVERBOSE,KTIMENUM,KTIMEVAL, &
635 KLEVEL,OVARZDIM,ODOUBLE )
647 CHARACTER(LEN=NF90_MAX_NAME),
INTENT(IN) :: HVNAME
649 REAL,
DIMENSION(:,:),
INTENT(IN) :: PWRITE
651 INTEGER,
INTENT(IN) :: KLISTING, KNCID
653 LOGICAL,
INTENT(IN) :: OVERBOSE
655 INTEGER,
INTENT(IN),
OPTIONAL :: KTIMENUM
656 INTEGER,
INTENT(IN),
OPTIONAL :: KTIMEVAL
658 INTEGER,
INTENT(IN),
OPTIONAL :: KLEVEL
659 LOGICAL,
INTENT(IN),
OPTIONAL :: OVARZDIM
660 LOGICAL,
INTENT(IN),
OPTIONAL :: ODOUBLE
664 CHARACTER(LEN=NF90_MAX_NAME) :: YWORK
666 REAL*4,
DIMENSION(SIZE(PWRITE,1),SIZE(PWRITE,2)) :: ZWRITE
668 INTEGER,
DIMENSION(4) :: ISTART, ICOUNT
670 INTEGER :: IUNLIMID, ITIMEID, ILENGHT, INDIM
671 INTEGER :: IC, IVAR_ID
674 REAL(KIND=JPRB) :: ZHOOK_HANDLE
678 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCSTORE',0,zhook_handle)
680 IF(
PRESENT(odouble))
THEN 686 IF(
PRESENT(klevel).AND..NOT.
PRESENT(ovarzdim))
THEN 687 WRITE(klisting,*)
'NCSTORE for TRIP : Error writing variable dimenssion ',hvname(1:len_trim(hvname))
688 WRITE(klisting,*)
'ILEVEL present but not LVARZDIM' 689 WRITE(klisting,*)nf90_strerror(ic)
690 CALL abort_trip(
'NCSTORE for TRIP : Error writing variable dimenssion')
693 IF(
PRESENT(ktimenum).AND.
PRESENT(ktimeval))
THEN 694 ic = nf90_inquire(kncid,unlimiteddimid=iunlimid)
696 ic = nf90_inquire_dimension(kncid,iunlimid,len=ilenght)
697 IF(ktimenum/=ilenght)
THEN 698 ic = nf90_inq_varid(kncid,
'time',itimeid)
699 ic = nf90_put_var(kncid,itimeid,ktimeval,start=(/ktimenum/))
701 WRITE(klisting,*)
'NCSTORE : re-writing of time variable number=',&
702 ktimenum,
' and value=',ktimeval
708 ywork = hvname(1:len_trim(hvname))
710 ic = nf90_inq_varid(kncid,ywork,ivar_id)
711 ic = nf90_inquire_variable(kncid,ivar_id,ndims=indim)
713 icount(1) =
SIZE(pwrite,1)
714 icount(2) =
SIZE(pwrite,2)
721 zwrite(:,:) = pwrite(:,:)
723 IF(
PRESENT(klevel).AND.ovarzdim)
THEN 725 IF(
PRESENT(ktimenum).AND.
PRESENT(ktimeval))
THEN 728 ic = nf90_put_var(kncid,ivar_id,pwrite,start=istart(1:4),
count=icount(1:4))
730 ic = nf90_put_var(kncid,ivar_id,zwrite,start=istart(1:4),
count=icount(1:4))
734 ic = nf90_put_var(kncid,ivar_id,pwrite,start=istart(1:3),
count=icount(1:3))
736 ic = nf90_put_var(kncid,ivar_id,zwrite,start=istart(1:3),
count=icount(1:3))
740 IF(
PRESENT(ktimenum).AND.
PRESENT(ktimeval))
THEN 743 ic = nf90_put_var(kncid,ivar_id,pwrite,start=istart(1:3),
count=icount(1:3))
745 ic = nf90_put_var(kncid,ivar_id,zwrite,start=istart(1:3),
count=icount(1:3))
749 ic = nf90_put_var(kncid,ivar_id,pwrite,start=istart(1:2),
count=icount(1:2))
751 ic = nf90_put_var(kncid,ivar_id,zwrite,start=istart(1:2),
count=icount(1:2))
756 IF(ic/=nf90_noerr)
THEN 757 WRITE(klisting,*)
'NCSTORE for TRIP : Error writing variable ',hvname(1:len_trim(hvname))
758 WRITE(klisting,*)nf90_strerror(ic)
759 CALL abort_trip(
'NCSTORE for TRIP : Error writing variable')
761 WRITE(klisting,*)
'NCSTORE for TRIP : Success in writing variable ',hvname(1:len_trim(hvname))
762 IF(
PRESENT(klevel))
WRITE(klisting,*)
' level: ',klevel
765 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCSTORE',1,zhook_handle)
772 SUBROUTINE ncdate(KLISTING,KNCID,KYEAR,KMONTH,KDAY,PTIME,OVERBOSE)
784 INTEGER,
INTENT(IN) :: KLISTING
785 INTEGER,
INTENT(IN) :: KNCID
786 INTEGER,
INTENT(IN) :: KYEAR
787 INTEGER,
INTENT(IN) :: KMONTH
788 INTEGER,
INTENT(IN) :: KDAY
789 REAL,
INTENT(IN) :: PTIME
791 LOGICAL,
INTENT(IN) :: OVERBOSE
795 CHARACTER(LEN=NF90_MAX_NAME) :: YDATE, YWORK
797 REAL*4,
DIMENSION(4) :: ZDATE
800 INTEGER :: IC, IWORK, ISIZE, IDATEID
802 REAL(KIND=JPRB) :: ZHOOK_HANDLE
806 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCDATE',0,zhook_handle)
810 zdate(1) =
REAL(kyear)
811 zdate(2) =
REAL(kmonth)
812 zdate(3) =
REAL(kday)
813 zdate(4) =
REAL(ptime)
817 ic = nf90_redef(kncid)
821 ic = nf90_def_dim(kncid,ydate,iwork,isize)
822 IF(ic/=nf90_noerr)
THEN 823 WRITE(klisting,*)
'NCDATE for TRIP : Error creating date dimension ' 824 WRITE(klisting,*)nf90_strerror(ic)
825 CALL abort_trip(
'NCDATE for TRIP : Error creating date dimension')
829 ic = nf90_def_var(kncid,ydate,nf90_float,isize,idateid)
830 ywork =
'current date: year month day time' 831 ic = nf90_put_att(kncid,idateid,
'long_name',ywork)
832 IF(ic/=nf90_noerr)
THEN 833 WRITE(klisting,*)
'NCDATE for TRIP : Error creating date attributs' 834 WRITE(klisting,*)nf90_strerror(ic)
835 CALL abort_trip(
'NCDATE for TRIP : Error creating date attributs')
837 ic = nf90_put_att(kncid,idateid,
'_FillValue',zmissval)
838 IF(ic/=nf90_noerr)
THEN 839 WRITE(klisting,*)
'NCDATE for TRIP : Error creating date attributs' 840 WRITE(klisting,*)nf90_strerror(ic)
841 CALL abort_trip(
'NCDATE for TRIP : Error creating date attributs')
844 ic = nf90_enddef(kncid)
846 ic = nf90_put_var(kncid,idateid,zdate)
849 WRITE(klisting,*)
'NCDATE for TRIP : Sucess in writting current date' 852 IF (
lhook)
CALL dr_hook(
'MODE_TRIP_NETCDF:NCDATE',1,zhook_handle)
subroutine ncstore(KLISTING, KNCID, HVNAME, PWRITE, OVERBOSE, KTIMENUM, KTIMEVAL, KLEVEL, OVARZDIM, ODOUBLE)
subroutine ncread_xy(KLISTING, KNCID, HVNAME, PVECT, OVERBOSE)
subroutine ncread_xyz(KLISTING, KNCID, HVNAME, PVECT, OVERBOSE)
subroutine ncclose(KLISTING, OVERBOSE, HFILENAME, KNCID)
subroutine ncdate(KLISTING, KNCID, KYEAR, KMONTH, KDAY, PTIME, OVERBOSE)
subroutine ncread_x(KLISTING, KNCID, HVNAME, PVECT, OVERBOSE)
subroutine nccreate(KLISTING, HFILENAME, HTITLE, HTIMEUNIT, HVNAME, HVLNAME, HUNIT, PLON, PLAT, PMISSVAL, OVERBOSE, KNCID, OTIME, KZLEN, OVARZDIM, ODOUBLE)
subroutine abort_trip(YTEXT)
subroutine ncopen(KLISTING, ORW, OVERBOSE, HFILENAME, KNCID)