[2018-03-27_15:13]
[/work2/am/2017.KYUSHU-HOKUBU.HEAVY.RAIN/JCOPE2M.HEAT.BUDGET/TEST.PLOT.SNAPSHOT]
[am@localhost]
$ test.plot.snapshot.sh
SOURCE FILE:
-rw-rw-r-- 1 am am 18K 3月 26 15:38 test.plot.snapshot.f90
Compiling test.plot.snapshot.f90 ...
ifort -convert big_endian -assume byterecl test.plot.snapshot.f90 -o test.plot.snapshot.exe
Done Compile.
-rwxrwxr-x 1 am am 878K 3月 27 15:16 test.plot.snapshot.exe
test.plot.snapshot.exe is running ...
### isyy=2017 ismm=07 isdd=31
### ieyy=2017 iemm=07 iedd=31
### xlons= 107.95834 ylats= 10.45833
### intday= 1
### stepday= 1DY
### depthfile=/misc/raid112MANDA/KakenA.H28/JCOPE2M/201706-07/data/basic.dat
### READ START
### read runz 123456
### READING INPUT DATA ...
# read(U_) um80 2017 7 31 180
# read(V_) vm80 2017 7 31 180
# read(T_) tm80 2017 7 31 180
# read(S_) sm80 2017 7 31 180
# NUM,yymmdd = 0 2017 7 31
### VERTICAL INTERPOLATION ...
### DONE.
### DESTAGGER (Arakawa C-gird to A-grid) ...
### DONE.
### OUTPUT BINARY FILE ...
### OUTPUT = /misc/raid112MANDA/KakenA.H28/JCOPE2M/201706-07/bin/JCOPE2M20170731.dat
### DONE.
### CREATE GrADS CONTROL FILE ...
### CTL FILE = JCOPE2M.ctl
### DONE.
Done test.plot.snapshot.exe
合計 14G
-rw-rw-r-- 1 am am 230M 3月 26 15:41 JCOPE2M20170601.dat
-rw-rw-r-- 1 am am 230M 3月 26 15:41 JCOPE2M20170602.dat
-rw-rw-r-- 1 am am 230M 3月 26 15:41 JCOPE2M20170603.dat
-rw-rw-r-- 1 am am 230M 3月 26 15:42 JCOPE2M20170604.dat
....
-rw-rw-r-- 1 am am 230M 3月 26 15:53 JCOPE2M20170727.dat
-rw-rw-r-- 1 am am 230M 3月 26 15:53 JCOPE2M20170728.dat
-rw-rw-r-- 1 am am 230M 3月 26 15:53 JCOPE2M20170729.dat
-rw-rw-r-- 1 am am 230M 3月 26 15:53 JCOPE2M20170730.dat
-rw-rw-r-- 1 am am 230M 3月 27 15:16 JCOPE2M20170731.dat
[2018-03-27_15:26]
[/work2/am/2017.KYUSHU-HOKUBU.HEAVY.RAIN/JCOPE2M.HEAT.BUDGET/TEST.PLOT.SNAPSHOT]
[am@localhost]
$ test.plot.snapshot.pl.sh
Tue_Mar_27_15:26:55_JST_2018
Grid Analysis and Display System (GrADS) Version 2.1.a3
Copyright (C) 1988-2015 by the Institute for Global Environment and Society (IGES)
GrADS comes with ABSOLUTELY NO WARRANTY
See file COPYRIGHT for more information
Config: v2.1.a3 little-endian readline grib2 netcdf hdf4-sds hdf5 opendap-grids,stn geotiff shapefile cairo
Issue 'q config' command for more detailed configuration information
GX Package Initialization: Size = 8.5 11
Running in Batch mode
sdate 00Z01JUN2017
edate 00Z01JUN2017
lev -10
OPEN JCOPE2M.ctl
Time = 00Z01JUN2017 to 00Z01JUN2017 Thu to Thu
OUTPUT: Fig/JCOPE2M.lev-10_20170601.eps
------------------------------
List of the following files:
------------------------------
test.plot.snapshot.sh
test.plot.snapshot.f90
test.plot.snapshot.pl.sh
test.plot.snapshot.gs
------------------------------
Machine info
------------------------------
localhost
/work2/am/2017.KYUSHU-HOKUBU.HEAVY.RAIN/JCOPE2M.HEAT.BUDGET/TEST.PLOT.SNAPSHOT
Tue Mar 27 15:28:06 JST 2018
======================
test.plot.snapshot.sh
======================
#!/bin/bash
#
# Description:
#
src=$(basename $0 .sh).f90
exe=$(basename $0 .sh).exe
nml=$(basename $0 .sh).nml
f90=ifort
opt="-convert big_endian -assume byterecl"
#opt="-CB -traceback -fpe0 -convert big_endian -assume byterecl"
prog_name=$(echo test_sample| sed -e 's/\./_/g')
datadir=/misc/raid112MANDA/KakenA.H28/JCOPE2M/201706-07/data/
odir=/misc/raid112MANDA/KakenA.H28/JCOPE2M/201706-07/bin
#odir=.
mkdir -vp $odir
isyy=2017
ismm=06
isdd=01
ieyy=2017
iemm=06
iedd=02
cat <<EOF>$nml
¶
datadir="$datadir"
odir="$odir"
isyy = $isyy
ismm = $ismm
isdd = $isdd
ieyy = $ieyy
iemm = $iemm
iedd = $iedd
&end
EOF
echo
echo SOURCE FILE:
echo
ls -lh ${src}
echo
echo Compiling ${src} ...
echo
echo ${f90} ${opt} ${src} -o ${exe}
echo
${f90} ${opt} ${src} -o ${exe}
if [ $? -ne 0 ]; then
echo
echo "=============================================="
echo
echo " COMPILE ERROR!!!"
echo
echo "=============================================="
echo
echo TERMINATED.
echo
exit 1
fi
echo "Done Compile."
echo
ls -lh ${exe}
echo
echo
echo ${exe} is running ...
echo
${exe} < $nml
if [ $? -ne 0 ]; then
echo
echo "=============================================="
echo
echo " ERROR in $exe: RUNTIME ERROR!!!"
echo
echo "=============================================="
echo
echo Terminated.
echo
exit 1
fi
echo
echo "Done ${exe}"
echo
ls -lh $odir |head -5;echo;echo " ...." ;echo; ls -lh $odir|tail -5;echo
----------------------
End of test.plot.snapshot.sh
----------------------
======================
test.plot.snapshot.f90
======================
!===============================================================
!
! PRINT JCOPE2M DATASET
!
! [2018-03-26_14:40]
! [/work2/am/2017.KYUSHU-HOKUBU.HEAVY.RAIN/JCOPE2M.HEAT.BUDGET/TEST.PLOT.SNAPSHOT]
! [am@localhost]
!
!===============================================================
program main
!------------------------------
parameter( im=866, jm=620, km=47 )
parameter( klev=28 )
!------------------------------
real &
z(im,jm,km),zz(im,jm,km),dz(im,jm,km)
real &
dxlon(im),dylat(jm),xlon(im),ylat(jm),xlonu(im),ylatu(jm) &
,xlonv(im),ylatv(jm)
real el(im,jm),wtsurf(im,jm)
real u(im,jm,km),v(im,jm,km),t(im,jm,km),s(im,jm,km)
real h0(klev)
real ulev(im,jm,klev),vlev(im,jm,klev) &
,tlev(im,jm,klev),slev(im,jm,klev)
real udsg(im,jm,klev),vdsg(im,jm,klev)
character(3) :: namemonth(12),cmon
data namemonth/'JAN','FEB','MAR','APR','MAY','JUN' &
,'JUL','AUG','SEP','OCT','NOV','DEC'/
character datadir*80,depthfile*200,outdir*80
character odir*300
character varname*100
character(500) :: ctlname,datname
character(2) :: cmonth,cday
character(4) :: cyear,&
fildscH,fildscEL,fildscWT,fildscU,fildscV,fildscT,fildscS
character(5) :: stepday
data h0/ &
-0., -5., -10., -50., -75., -100., -150., -200., -250., -300. &
, -400., -500., -600., -700., -800., -900., -1000., -1200. &
, -1400., -1600., -1800., -2000., -2500, -3000., -3500., -4000. &
, -5000., -6000. /
namelist /para/datadir,odir,isyy,ismm,isdd,ieyy,iemm,iedd
!------- SETTING START -------
read(5,nml=para)
xlons = 108 - 1./24.
ylats = 10.5- 1./24.
intday=1
stepday=' 1DY'
depthfile=trim(datadir)//'basic.dat'
print '(A,i4.4,A,i2.2,A,i2.2)',&
' ### isyy=',isyy,' ismm=',ismm,' isdd=',isdd
print '(A,i4.4,A,i2.2,A,i2.2)',&
' ### ieyy=',isyy,' iemm=',iemm,' iedd=',iedd
print '(A,f12.5,A,f12.5)',' ### xlons=',xlons,' ylats=',ylats
print '(A,i5)',' ### intday=',intday
print '(A,A)',' ### stepday=',stepday
print '(A,A)',' ### depthfile=',trim(depthfile)
print *
!------- SETTING END -------
!-------- READ START ---------
print *,'### READ START'
tbias = 10.
sbias = 35.
valmis = -99.9999
open(10,file=depthfile,form='unformatted')
read(10) fildscH,Z,ZZ,DZ,ichflg
close(10)
print *,'### read ',fildscH,ichflg
timejde = julday(ieyy,iemm,iedd)
timejd = julday(isyy,ismm,isdd)
nm = int(timejde -timejd)/intday +1
TLOOP: do n=1,nm
call caldat(int(timejd+0.0001),iyy,imm,idd)
write(cmonth,'(i2.2)') imm
write(cday,'(i2.2)') idd
write(cyear,'(i4.4)') iyy
print *, '### READING INPUT DATA ...'
!-----------------------------------------------------------
! write(6,*) '# read EL_ '
! open(10,file=datadir(1:lnblnk(datadir))
! & //'EL_'//cyear//cmonth//cday
! & ,form='unformatted')
! read(10) fildscEL,iyf,imf,idf,el,mall
! close(10)
! write(6,*) '### read END (EL_) ',fildscEL,iyf,imf,idf,mall
!-----------------------------------------------------------
! open(10,file=datadir(1:lnblnk(datadir)) &
! //'WTSURF_'//cyear//cmonth//cday,form='unformatted')
! read(10) fildscWT,iyf,imf,idf,wtsruf !,mall
! close(10)
! print *,'# read (WTSURF_) ',fildscWT,iyf,imf,idf,mall
! print *,wtsurf
!-----------------------------------------------------------
open(10,file=datadir(1:lnblnk(datadir)) &
//'U_'//cyear//cmonth//cday &
,form='unformatted')
read(10) &
fildscU,iyf,imf,idf &
,(((U(i,j,k),i=1,im),j=1,jm),k=1,km-1),mall
close(10)
print *, '# read(U_) ',fildscU,iyf,imf,idf,mall
!-----------------------------------------------------------
open(10,file=datadir(1:lnblnk(datadir)) &
//'V_'//cyear//cmonth//cday &
,form='unformatted')
read(10) &
fildscV,iyf,imf,idf &
,(((V(i,j,k),i=1,im),j=1,jm),k=1,km-1),mall
close(10)
print *, '# read(V_) ',fildscV,iyf,imf,idf,mall
!-----------------------------------------------------------
open(10,file=datadir(1:lnblnk(datadir)) &
//'T_'//cyear//cmonth//cday &
,form='unformatted')
read(10) &
fildscT,iyf,imf,idf &
,(((T(i,j,k),i=1,im),j=1,jm),k=1,km-1),mall
close(10)
print *,'# read(T_) ',fildscT,iyf,imf,idf,mall
!-----------------------------------------------------------
open(10,file=datadir(1:lnblnk(datadir)) &
//'S_'//cyear//cmonth//cday &
,form='unformatted')
read(10) &
fildscS,iyf,imf,idf &
,(((S(i,j,k),i=1,im),j=1,jm),k=1,km-1),mall
close(10)
print *, '# read(S_) ',fildscS,iyf,imf,idf,mall
!-----------------------------------------------------------
print *,'# NUM,yymmdd = ',numfil,iyy,imm,idd
!------- Make GrADS data START -------
do j=1,jm
do i=1,im
if( abs(z(i,j,km)) <= 1. ) then
el(i,j) = valmis
end if
enddo !i
enddo !j
print *, '### VERTICAL INTERPOLATION ...'
do k=1,klev
call gsigtoz3 &
( zz,u,ulev,h0(k),im,jm,km,klev,k,valmis )
call gsigtoz3 &
( zz,v,vlev,h0(k),im,jm,km,klev,k,valmis )
call gsigtoz3 &
( zz,t,tlev,h0(k),im,jm,km,klev,k,valmis )
call gsigtoz3 &
( zz,s,slev,h0(k),im,jm,km,klev,k,valmis )
enddo
print *, '### DONE.'
do k=1,klev
do j=1,jm
do i=1,im
if( abs(z(i,j,km)) <= 1. ) then
ulev(i,j,k) = valmis
vlev(i,j,k) = valmis
tlev(i,j,k) = valmis
slev(i,j,k) = valmis
else
tlev(i,j,k) = tlev(i,j,k) +tbias
if (tlev(i,j,k) < -20.0) tlev(i,j,k)=valmis
slev(i,j,k) = slev(i,j,k) +sbias
if (slev(i,j,k) < -20.0) slev(i,j,k)=valmis
end if
enddo
enddo
enddo
print *, '### DESTAGGER (Arakawa C-gird to A-grid) ...'
call destagger(im,jm,klev,ulev,vlev,valmis,udsg,vdsg)
print *, '### DONE.'
print *, '### OUTPUT BINARY FILE ...'
iu=20
call outbin(iu,odir,iyy,imm,idd,im,jm,klev,tlev,slev,udsg,vdsg)
print *, '### DONE.'
timejd = timejd +intday
end do TLOOP
call degset( dxlon,dylat,xlon,ylat,xlonu,ylatu,xlonv,ylatv &
,xlons,ylats,im,jm )
iu=20
itmax=1
cmon=namemonth(ismm)
print *, '### CREATE GrADS CONTROL FILE ...'
call create_ctl(iu,odir,im,jm,klev,nm,&
&valmis,isyy,cmon,isdd,stepday,xlon,ylat,h0)
print *, '### DONE.'
end program main
subroutine outbin(iu,odir,iyy,imm,idd,im,jm,klev,&
tlev,slev,udsg,vdsg)
integer,intent(in)::iu
character(len=*),intent(in)::odir
integer,intent(in)::iyy,imm,idd
integer,intent(in)::im,jm,klev
!real,intent(in),dimension(im,jm)::wtsurf
real,intent(in),dimension(im,jm,klev)::tlev,slev,udsg,vdsg
real,dimension(im,jm)::dum
character(len=1000)::ofle
is=1
ie=len_trim(odir)
write(ofle(is:ie),'(A)')trim(odir)
is=ie+1
ie=is
write(ofle(is:ie),'(A)')'/'
is=ie+1
ie=is+len('JCOPE2M')
write(ofle(is:ie),'(A)')'JCOPE2M'
is=ie
ie=is+3
write(ofle(is:ie),'(i4.4)')iyy
is=ie+1
ie=is+1
write(ofle(is:ie),'(i2.2)')imm
is=ie+1
ie=is+1
write(ofle(is:ie),'(i2.2)')idd
is=ie+1
ie=is+len('.dat')
write(ofle(is:ie),'(A)')'.dat'
print '(A,A)',' ### OUTPUT = ',trim(ofle)
open (iu,file=trim(ofle),form='unformatted',access='direct',&
&recl=im*jm*klev*4)
irec=1
write(iu,rec=irec) (((tlev(i,j,k),i=1,im),j=1,jm),k=klev,1,-1)
irec=irec+1
write(iu,rec=irec) (((slev(i,j,k),i=1,im),j=1,jm),k=klev,1,-1)
irec=irec+1
write(iu,rec=irec) (((udsg(i,j,k),i=1,im),j=1,jm),k=klev,1,-1)
irec=irec+1
write(iu,rec=irec) (((vdsg(i,j,k),i=1,im),j=1,jm),k=klev,1,-1)
close(iu)
end subroutine outbin
subroutine create_ctl(iu,odir,im,jm,klev,nm,&
&valmis,isyy,cmon,isdd,stepday,xlon,ylat,h0)
integer,intent(in)::iu
character(len=*),intent(in)::odir
integer,intent(in)::im,jm,klev,nm
real,intent(in):: xlon(im),ylat(jm),h0(klev)
character(len=*),intent(in)::cmon
character(len=*),intent(in)::stepday
character(len=90)::ctl
character(len=1000)::ofle
is=ie+1
ie=is+len('JCOPE2M.ctl')
write(ctl(is:ie),'(A)')'JCOPE2M.ctl'
print '(A,A)',' ### CTL FILE = ',trim(ctl)
is=1
ie=len_trim(odir)
write(ofle(is:ie),'(A)')trim(odir)
is=ie+1
ie=is
write(ofle(is:ie),'(A)')'/'
is=ie+1
ie=is+LNBLNK('JCOPE2M%y4%m2%d2.dat')-1
write(ofle(is:ie),'(A)')trim('JCOPE2M%y4%m2%d2.dat')
ie_ofle=ie
open(iu,file=trim(ctl),form='formatted')
write(iu,'(A,A)') 'dset ',ofle(1:ie_ofle)
write(iu,'(A)') 'options big_endian'
write(iu,'(A)') 'options template'
write(iu,'(A)') 'title JCOPE2M'
if( valmis > -1.e+3) then
write(iu,'(A,f10.4)') 'undef ',valmis
else
write(iu,'(A)') 'undef -1.e+21'
end if
write(iu,2000) 'xdef ',im,' levels ' ,(xlon(i), i=1,im)
write(iu,2000) 'ydef ',jm,' levels ' ,(ylat(j), j=1,jm)
write(iu,'(a,i5,a,/,6(1x,f10.3))') &
'zdef ',klev,' levels ',(h0(k), k=klev,1,-1)
write(iu,'(a,i4,a,i2,a3,i4,1x,a)') &
'tdef ',nm,' linear ',isdd,cmon,isyy,trim(stepday)
write(iu,'(A)') 'vars 4'
write(iu,'(A,i5,A)') 't ',klev,' 0 t'
write(iu,'(A,i5,A)') 's ',klev,' 0 s'
write(iu,'(A,i5,A)') 'u ',klev,' 0 u'
write(iu,'(A,i5,A)') 'v ',klev,' 0 v'
write(iu,'(A)') 'endvars'
close(iu)
1020 format(a,f8.4)
2000 format(a,i4,a,/,6(1x,f12.6))
3000 format(a,i4,a)
end subroutine create_ctl
subroutine destagger(im,jm,klev,ulev,vlev,valmis,udsg,vdsg)
integer,intent(in)::im,jm,klev
real,intent(in)::ulev(im,jm,klev),vlev(im,jm,klev)
real,intent(in)::valmis
real,intent(inout)::udsg(im,jm,klev),vdsg(im,jm,klev)
udsg=ulev
vdsg=vlev
! u
do k=1,klev
do j=1,jm
do i=1,im-1
udsg(i,j,k)=0.5*(ulev(i,j,k)+ulev(i+1,j,k))
end do !i
end do !j
end do !k
do k=1,klev
do j=1,jm
udsg(im,j,k)=ulev(im,j,k)
end do
end do
do k=1,klev
do j=1,jm
do i=1,im-1
if(ulev(i,j,k)==valmis.or.ulev(i+1,j,k)==valmis)then
udsg(i,j,k)=valmis
end if
end do
end do
end do
! v
do k=1,klev
do j=1,jm-1
do i=1,im
vdsg(i,j,k)=0.5*(vlev(i,j,k)+vlev(i,j+1,k))
end do !i
end do !j
end do !k
do k=1,klev
do i=1,im
vdsg(i,jm,k)=vlev(i,jm,k)
end do !i
end do !km
do k=1,klev
do j=1,jm
do i=1,im-1
if(vlev(i,j,k)==valmis.or.vlev(i,j+1,k)==valmis)then
vdsg(i,j,k)=valmis
end if
end do
end do
end do
end subroutine destagger
!*****************************************************************************
subroutine degset &
(dxlon,dylat,xlon,ylat,xlonu,ylatu,xlonv,ylatv,xlons,ylats,im,jm)
dimension dxlon(im),dylat(jm),xlon(im),ylat(jm)
dimension xlonu(im),ylatu(jm),xlonv(im),ylatv(jm)
do i=1,im
dxlon(i) = 1./12.
enddo
xlon(1) = xlons
do i=2,im
xlon(i) = xlon(i-1) +( dxlon(i-1)+dxlon(i) )/2.
enddo
xlonu(1) = xlons -dxlon(1)/2.
do i=2,im
xlonu(i) = xlonu(i-1) +dxlon(i-1)
enddo
xlonv(1) = xlons
do i=2,im
xlonv(i) = xlonv(i-1) +( dxlon(i-1)+dxlon(i) )/2.
enddo
! do i=1,im
! write(6,*) 'in degset i,xlonu = ',i,xlonu(i)
! end do
do j=1,jm
dylat(j) = 1./12.
enddo
ylat(1) = ylats
do j=2,jm
ylat(j) = ylat(j-1) +( dylat(j-1)+dylat(j) )/2.
enddo
ylatu(1) = ylats
do j=2,jm
ylatu(j) = ylatu(j-1) +( dylat(j-1)+dylat(j) )/2.
enddo
ylatv(1) = ylats -dylat(1)/2.
do j=2,jm
ylatv(j) = ylatv(j-1) +dylat(j-1)
enddo
! do j=1,jm
! write(6,*) 'in degset j,ylatv = ',j,ylatv(j)
! end do
return
end subroutine degset
!*****************************************************************************
SUBROUTINE GSIGTOZ3(ZZ,T,TLEV,ZLEV,IM,JM,KM &
,klev,kk,valmis)
!-------------------------------------------------------------------
! THIS ROUTINE LINERLY INTERPOLATES TLEV AT THE LEVEL,
! ZLEV, FROM T LOCATED ON SIGMA LEVELS, ZZ.
! ZLEV AND ZZ ARE NEGATIVE QUANTITIES.
! A SEARCH IS MADE TO FIND ZZ(I,J,K) AND ZZ(I,J,K+1) WHICH
! BRACKETS ZLEV; THEN THE INTERPOLATION IS MADE.
! IN THE REGION < 0 and > ZZ(I,J,1) AND, IN THE REGION,
! < ZZ(I,J,KM-1) and > -1, DATA IS EXTRAPOLATED.
! NOTE THAT A NEW MASK ,FSM, IS CREATED.
! THE VALUES OF H SUPPLIED TO THIS SUBROUTINE SHOULD BE
! APPROPRIATE TO THE VARIABLE T. FOR EXAMPLE, IF T IS THE
! X-COMPONENT OF VELOCITY, H SHOULD THE AVERAGE OF DEPTH
! AT I AND I-1.
!-------------------------------------------------------------------
DIMENSION ZZ(IM,JM,KM),T(IM,JM,KM),TLEV(IM,JM,KLEV)
do j=1,jm
do i=1,im
tlev(i,j,kk) = valmis
enddo
enddo
! c i = 227
! c j = 300
! c write(0,*) 'gsigoz3 zlev ',zlev
! c do k=1,km
! c write(0,*) 'zz t ',zz(i,j,k),t(i,j,k)
! c enddo
DO 200 J=1,JM
DO 200 I=1,IM
IF(ZLEV > ZZ(I,J,1) .AND. t(i,j,1) > valmis+1. .AND. &
t(i,j,2) > valmis+1. ) THEN
K=1
TLEV(I,J,kk)=T(I,J,K)+(ZLEV-ZZ(I,J,K)) &
*(T(I,J,K+1)-T(I,J,K))/(ZZ(I,J,K+1)-ZZ(I,J,K))
! c if( i.eq.227 .and. j.eq.300 ) then
! c write(0,*) 'gsigoz3 zlev>zz(1)'
! c write(0,*) 'kk k k+1 ',kk,k,k+1
! c write(0,*) 'zlev zz(k) zz(k+1) '
! c & ,zlev,zz(i,j,k),zz(i,j,k+1)
! c write(0,*) 'tlev t(k) t(k+1) '
! c & ,tlev(i,j,kk),t(i,j,k),t(i,j,k+1)
! c end if
GO TO 200
ENDIF
K=1
100 CONTINUE
IF(ZLEV <= ZZ(I,J,K) .AND. ZLEV >= ZZ(I,J,K+1) .AND. &
t(i,j,k) > valmis+1. .AND. t(i,j,k+1) > valmis+1. ) THEN
TLEV(I,J,kk)=T(I,J,K)+(ZLEV-ZZ(I,J,K)) &
*(T(I,J,K+1)-T(I,J,K))/(ZZ(I,J,K+1)-ZZ(I,J,K))
GO TO 200
! c if( i.eq.227 .and. j.eq.300 ) then
! c write(0,*) 'gsigoz3 zz(k)>zlev>zz(k+1)'
! c write(0,*) 'k kk k+1 ',k,kk,k+1
! c write(0,*) 'zz(k) zlev zz(k+1) '
! c & ,zz(i,j,k),zlev,zz(i,j,k+1)
! c write(0,*) 't(k) tlev t(k+1) '
! c & ,t(i,j,k),tlev(i,j,kk),t(i,j,k+1)
! c end if
ELSE
K=K+1
IF(K < (KM-1)) GO TO 100
ENDIF
IF(ZLEV <= ZZ(I,J,K) .AND. ZLEV >= ZZ(I,J,KM) .AND. &
t(i,j,k) > valmis+1. .AND. t(i,j,k-1) > valmis+1. ) THEN
TLEV(I,J,kk)=T(I,J,K-1)+(ZLEV-ZZ(I,J,K-1)) &
*(T(I,J,K)-T(I,J,K-1))/(ZZ(I,J,K)-ZZ(I,J,K-1))
! c if( i.eq.227 .and. j.eq.300 ) then
! c write(0,*) 'gsigoz3 zz(k)>zlev>zz(km)'
! c write(0,*) 'k-1 k kk ',k-1,k,kk
! c write(0,*) 'zz(k-1) zz(k) zlev '
! c & ,zz(i,j,k-1),zz(i,j,k),zlev
! c write(0,*) 't(k-1) t(k) tlev '
! c & ,t(i,j,k-1),t(i,j,k),tlev(i,j,kk)
! c end if
ENDIF
200 END DO
RETURN
END SUBROUTINE GSIGTOZ3
!*****************************************************************************
INTEGER FUNCTION JULDAY(IYYY,MONTH,DD)
! ********** SUBROUTINE DESCRIPTION:
! FINDS THE JULIAN DAY FROM A DATE.
! ********** ORIGINAL AUTHOR AND DATE:
! PRESS,FLANNERY,TEUKOLSKY,VETTERLING 1986.
! NUMERICAL RECIPES
! ********** REVISION HISTORY:
! ********** ARGUMENT DEFINITIONS:
INTEGER :: IYYY,MONTH,DD
! NAME IN/OUT DESCRIPTION
! IYYY I YEAR
! MONTH I MONTH (1 TO 12)
! DD I DAY OF MONTH
! JULDAY O JULIAN DAY
! ********** COMMON BLOCKS:
! NONE
! ********** LOCAL PARAMETER DEFINITIONS:
INTEGER :: IGREG
PARAMETER (IGREG = 15 + 31*(10 + 12*1582))
! ********** LOCAL VARIABLE DEFINITIONS:
INTEGER :: JY,JM,JA
! NAME DESCRIPTION
! ********** OTHER ROUTINES AND FUNCTIONS CALLED:
! INT - INTRINSIC TRUNCATE
!---+67--1----+----2----+----3----+----4----+----5----+----6----+----7--
IF (IYYY < 0) IYYY = IYYY + 1
IF (MONTH > 2) THEN
JY = IYYY
JM = MONTH + 1
ELSE
JY = IYYY - 1
JM = MONTH + 13
ENDIF
JULDAY = INT(365.25*JY) + INT(30.6001*JM) + DD + 1720995
IF (DD + 31*(MONTH + 12*IYYY) >= IGREG) THEN
JA = INT(0.01*JY)
JULDAY = JULDAY + 2 - JA + INT(0.25*JA)
ENDIF
RETURN
END FUNCTION JULDAY
!*****************************************************************************
SUBROUTINE CALDAT(JULIAN,IYYY,MONTH,DD)
! ********** SUBROUTINE DESCRIPTION:
! GIVEN THE JULIAN DAY, RETURNS THE YEAR, MONTH AND DAY OF MONTH.
! ********** ORIGINAL AUTHOR AND DATE:
! PRESS,FLANNERY,TEUKOLSKY,VETTERLING 1986.
! NUMERICAL RECIPES
! ********** REVISION HISTORY:
! ********** ARGUMENT DEFINITIONS:
INTEGER :: JULIAN,IYYY,MONTH,DD
! NAME IN/OUT DESCRIPTION
! JULIAN I THE JULIAN DAY
! IYYY O THE YEAR
! MONTH O THE MONTH (1 TO 12)
! DD O THE DAY OF THE MONTH
! ********** COMMON BLOCKS:
! NONE
! ********** LOCAL PARAMETER DEFINITIONS:
INTEGER :: IGREG
PARAMETER (IGREG=2299161)
! ********** LOCAL VARIABLE DEFINITIONS:
INTEGER :: JALPHA,JA,JB,JC,JD,JE
! NAME DESCRIPTION
! ********** OTHER ROUTINES AND FUNCTIONS CALLED:
!---+67--1----+----2----+----3----+----4----+----5----+----6----+----7--
IF (JULIAN >= IGREG) THEN
JALPHA = INT(((JULIAN - 1867216) - 0.25)/36524.25)
JA = JULIAN + 1 + JALPHA - INT(0.25*JALPHA)
ELSE
JA = JULIAN
ENDIF
JB = JA + 1524
JC = INT(6680. + ((JB - 2439870) - 122.1)/365.25)
JD = 365*JC + INT(0.25*JC)
JE = INT((JB - JD)/30.6001)
DD = JB - JD - INT(30.6001*JE)
MONTH = JE - 1
IF (MONTH > 12) MONTH = MONTH - 12
IYYY = JC - 4715
IF (MONTH > 2) IYYY = IYYY - 1
IF (IYYY <= 0) IYYY = IYYY - 1
RETURN
END SUBROUTINE CALDAT
----------------------
End of test.plot.snapshot.f90
----------------------
======================
test.plot.snapshot.pl.sh
======================
#!/bin/bash
export LANG=C
host=$(hostname)
cwd=$(pwd)
temp=$(date)
timestamp=$(echo $temp | sed -e 's/ /_/g')
echo $timestamp
gs=$(basename $0 .pl.sh).gs
grads -bcp "$gs -c JCOPE2M.ctl -l -10 \
-s 00Z01JUN2017 -e 00Z01JUN2017 \
-host $host \
-cwd $cwd \
-gs $gs \
-timestamp $timestamp \
-q"
exit
----------------------
End of test.plot.snapshot.pl.sh
----------------------
======================
test.plot.snapshot.gs
======================
#
# [2018-03-26_17:39]
# [/work2/am/2017.KYUSHU-HOKUBU.HEAVY.RAIN/JCOPE2M.HEAT.BUDGET/TEST.PLOT.SNAPSHOT]
# [am@localhost]
#
function sample(args)
ctl='ctl.ctl'
outdir='Fig'
*
* Decode options
*
i = 1
while( 1 )
arg = subwrd( args, i )
i = i + 1;
if( arg = '' ); break; endif
while( 1 )
if( arg = '-c' ); ctl = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-d' ); outdir = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-s' ); sdate = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-e' ); edate = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-l'); lev = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-host'); host = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-cwd'); cwd = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-gs'); gs = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-timestamp'); timestamp = subwrd(args,i) ;i=i+1; break; endif
if( arg = '-q'); quitflag=true; break; endif
say 'Syntax error : arg= 'arg
return
endwhile
endwhile
'!mkdir -vp 'outdir
say 'sdate 'sdate
say 'edate 'edate
say 'lev 'lev
say 'OPEN 'ctl
'open 'ctl
'set time 'sdate
'q dims'
lin=sublin(result,5)
wrd=subwrd(lin ,9)
t1=wrd
'set time 'edate
'q dims'
lin=sublin(result,5)
wrd=subwrd(lin ,9)
t2=wrd
'set t 't1' 't2
'q time'
say result
t=t1
while (t<=t2)
'cc'
'set t 't
'set lev 'lev
'set grid off'
'set grads off'
'set xlab off'
'set ylab off'
'set gxout shade2'
'set cmin 2'
'set cmax 32'
'set cint 2'
'rgbset2'
'd t'
'set gxout contour'
'set xlab off'
'set ylab off'
'set cmin 2'
'set cmax 32'
'set cint 2'
'set ccolor 0'
'set cthick 1'
'set clab off'
'd t'
'xcbar 1.5 7. 1.3 1.4 -line on -edge circle'
'set gxout vector'
'set xlab off'
'set ylab off'
'set ccolor 0'
'set cthick 10'
'set arrscl 0.2 1.5'
'set arrowhead -0.4'
'd skip(maskout(u, mag(u, v)-0.2),10,10);v'
'set gxout vector'
'set xlab off'
'set ylab off'
'set ccolor 1'
'set cthick 1'
'set arrscl 0.2 1.5'
'set arrowhead -0.4'
'd skip(maskout(u, mag(u, v)-0.2),10,10);v'
'set gxout contour'
'set xlab off'
'set ylab off'
'set cmin 4'
'set cmax 32'
'set cint 1'
'set clab off'
'set ccolor 0'
'set cthick 10'
'd s'
'set gxout contour'
'set xlab on'
'set ylab on'
'set cmin 4'
'set cmax 32'
'set cint 1'
'set ccolor 9'
'set cthick 1'
'set clab on'
'set clskip 2'
'd s'
'q dims'
line=sublin(result,5)
datetime=subwrd(line,6)
hh=substr(datetime,1,2)
dd=substr(datetime,4,2)
mmm=substr(datetime,6,3)
yyyy=substr(datetime,9,4)
if (mmm='JAN') ;mm=01 ;endif
if (mmm='FEB') ;mm=02 ;endif
if (mmm='MAR') ;mm=03 ;endif
if (mmm='APR') ;mm=04 ;endif
if (mmm='MAY') ;mm=05 ;endif
if (mmm='JUN') ;mm=06 ;endif
if (mmm='JUL') ;mm=07 ;endif
if (mmm='AUG') ;mm=08 ;endif
if (mmm='SEP') ;mm=09 ;endif
if (mmm='OCT') ;mm=10 ;endif
if (mmm='NOV') ;mm=11 ;endif
if (mmm='DEC') ;mm=12 ;endif
depth=-1*lev
title='JCOPE2M 'dd''mmm''yyyy' 'depth'm'
'draw title 'title
# Header
'set strsiz 0.1 0.12'
'set string 1 l 2'
'draw string 0.1 9.6 'timestamp' 'host
'draw string 0.1 9.45 'cwd
'draw string 0.1 9.3 'gs
out=outdir'/JCOPE2M.lev'lev'_'yyyy''mm''dd'.eps'
say 'OUTPUT: 'out
'gxprint 'out
t=t+1
endwhile ;*t
'allclose'
if (quitflag = 'true')
quit
endif
return
----------------------
End of test.plot.snapshot.gs
----------------------
[2018-03-27_15:33]
[/work2/am/2017.KYUSHU-HOKUBU.HEAVY.RAIN/JCOPE2M.HEAT.BUDGET/TEST.PLOT.SNAPSHOT/Fig]
[am@localhost]
$ convert -density 400 -delay 30 JCOPE2M.lev-10_20170*.eps JCOPE2M.lev-10_2017.gif