srvctl.f90 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411
  1. program srvctl
  2. ! SERVICE to GrADS format converter
  3. ! INPUT (on commandline):
  4. ! arg#1: name of inputfile with data in SERVICE-format
  5. ! The program assumes that the records of the file in SERVICE format
  6. ! are in the order produced by the PUMABURNER or AFTERBURNER program,
  7. ! i. e. the sequence and the number of codes for each timestep must not
  8. ! vary in the inputfile, different levels of one code at one timestep
  9. ! must be ordered and in sequence. Usually this is not a problem.
  10. !
  11. ! The order is (because GrADS assumes this): level - code - time
  12. !
  13. ! Author: Simon Blessing (blessing(at)dkrz.de), 27.11.2000
  14. ! Last maintenance: Oct. 27, 2004, got rid of complaint by sunf90
  15. ! Jan. 1st, 2005, added recunit detection
  16. !
  17. implicit none
  18. type levlst
  19. integer :: level
  20. type (levlst), pointer :: next
  21. end type levlst
  22. type codelst
  23. integer :: code
  24. integer :: nlev
  25. type (codelst), pointer :: next
  26. end type codelst
  27. type (codelst), pointer :: cl_first, cl_help, codelist
  28. type (levlst), pointer :: ll_first, levlist, ll_help
  29. integer :: nlon, nlat, nlev,j, deltat,fdate,fhour,flev, nlevmax
  30. integer :: dd,mm,mi,yyyy
  31. integer :: ddl,mml,mil,yyyyl,lhour,ldate
  32. integer :: tdef=1
  33. integer :: nvar=1
  34. logical :: l_vsect=.false.
  35. real, allocatable :: si(:) !sinus of gaussian latitude
  36. real, parameter :: PI = 3.14159265359 ! Pi
  37. character*2 :: ut
  38. character*16 :: vform = '(a1,i3,i6,a,3i6)'
  39. integer,parameter :: maxdpm(12) =(/ 31,29,31,30,31,30,31,31,30,31,30,31/)
  40. character*3,parameter :: month(12) = (/'jan','feb','mar','apr','may','jun',&
  41. 'jul','aug','sep','oct','nov','dec'/)
  42. integer irec, i, narg
  43. integer ih(8) !header of the SERVICE file
  44. REAL, ALLOCATABLE :: x(:) !data array
  45. Character*255 filein, filedes !filenames
  46. ! external iargc
  47. !==========================================================================
  48. narg=iargc()
  49. if (narg < 1) call usage
  50. call getarg(1,filein)
  51. if ( filein == '' .or. filein(1:1) == '-' ) call usage
  52. i = len_trim(filein)
  53. if (i > 4 .and. filein(i-3:i) == '.srv') i = i - 4
  54. filedes = filein(1:i) // '.ctl'
  55. if (narg > 1) call getarg(2,filedes)
  56. print*,'reading binary service file <',trim(filein),'>'
  57. print*,'writing grads control file <',trim(filedes),'>'
  58. !
  59. ! Open the input file
  60. !
  61. open(1,File=filein,FORM='UNFORMATTED',status='OLD')
  62. !
  63. ! Read the first header to know the length of the data array
  64. !
  65. irec=0
  66. read (1,End=999) ih
  67. nlon = ih(5)
  68. nlat = ih(6)
  69. nlev = ih(2) !very first guess
  70. if (ih(2)==-1) l_vsect=.true.
  71. flev = ih(2) !remember first level in file
  72. fdate = ih(3) !first date in dataset
  73. fhour = ih(4) !first hour in dataset
  74. deltat = 1 !first guess
  75. nlevmax = 1 !first guess
  76. if ( NLAT .lt. 32 .or. NLAT .gt. 1000) then
  77. print*,ih
  78. print*,'Unusual number of Latitudes in header of inputfile:', &
  79. & ' "" ',filein(1:len_trim(filein)),'"" : ',NLAT
  80. if (nlat.gt.1000) then
  81. print*,'Exiting.'
  82. stop
  83. end if
  84. end if
  85. !
  86. ALLOCATE( x(ih(5)*ih(6)) )
  87. !
  88. ! Open the descriptor file (formatted, sequential)
  89. !
  90. Open (3,File=filedes,FORM='FORMATTED',ACCESS='SEQUENTIAL')
  91. !
  92. ! Calculate Gaussian latitudes
  93. !
  94. if (l_vsect) then
  95. print*,'Presuming ',nlon,' gaussian latitudes.'
  96. allocate (si(nlon))
  97. call inigau(nlon,si)
  98. else
  99. print*,'Presuming ',nlat,' gaussian latitudes.'
  100. allocate (si(nlat))
  101. call inigau(nlat,si)
  102. end if
  103. !
  104. ! Write to descriptor file what we know so far
  105. !
  106. write(3,'(2a)')'DSET ^',trim(filein)
  107. ! write(3,'(2a)')'TITLE ',name
  108. if (kind(x)==8) then
  109. write(3,'(a)')'UNDEF 9e+99'
  110. else
  111. write(3,'(a)')'UNDEF 9e+09'
  112. end if
  113. ! write(3,'(a)')'OPTIONS BIG_ENDIAN'
  114. if (l_vsect) then
  115. write(3,'(a,i6,a,f8.4,a,f8.4)')'XDEF ',1, &
  116. & ' LEVELS ',1.
  117. write(3,'(a)')'OPTIONS YREV'
  118. write(3,'(a,i6,a/(8f10.4))')'YDEF ',nlon, &
  119. & ' LEVELS ',-180./PI*asin(si(:))
  120. else
  121. write(3,'(a,i6,a,f8.4,a,f8.4)')'XDEF ',nlon, &
  122. & ' LINEAR ',0.,' ',(360./float(nlon))
  123. write(3,'(a)')'OPTIONS YREV'
  124. if (nlat.gt.1) then
  125. write(3,'(a,i6,a/(8f10.4))')'YDEF ',nlat, &
  126. & ' LEVELS ',-180./PI*asin(si(:))
  127. else
  128. write(3,'(a,i6,a,f8.4,a,f8.4)')'YDEF ',nlat, &
  129. & ' LINEAR ',0.,' ',180.
  130. end if
  131. end if
  132. write(3,'(A)') 'OPTIONS SEQUENTIAL'
  133. write(3,'(A)') 'XYHEADER 40'
  134. !
  135. ! Start list about codes and levels
  136. !
  137. nvar=1
  138. allocate(cl_first)
  139. cl_first%code=ih(1)
  140. cl_first%nlev=min(ih(2),1) !special treatment of level zero
  141. nullify(cl_first%next) !last in list so far; Lahey and VAST need this
  142. allocate(ll_first)
  143. ll_first%level=ih(2)
  144. nullify(ll_first%next) !last in list so far
  145. !
  146. ! Do the real work (read & write data)
  147. !
  148. irec = 1 !start record count
  149. do
  150. !
  151. ! Read the input data
  152. !
  153. read (1,END=999) x
  154. !
  155. ! Do for the next record
  156. !
  157. read (1,End=999) ih
  158. irec = irec+1 !count records
  159. !
  160. ! count time
  161. !
  162. if ( ih(1)==cl_first%code.and.ih(2)==flev ) then
  163. !first code first level again?
  164. tdef=tdef+1 !so time must have advanced!
  165. end if
  166. !
  167. ! get first non-zero level
  168. !
  169. if (ll_first%level==0 .and. ih(2).ne.0) then
  170. ll_first%level=ih(2)
  171. end if
  172. !
  173. ! Work on codelist
  174. !
  175. codelist=>cl_first
  176. do while (ih(1).ne.codelist%code)
  177. if (.not. associated(codelist%next)) then !add new code to list
  178. cl_help=>codelist
  179. allocate(codelist)
  180. cl_help%next=>codelist
  181. nvar=nvar+1 !count codes
  182. codelist%code=ih(1)
  183. codelist%nlev=min(ih(2),1) !special treatment of level zero
  184. nullify(codelist%next) !this is last in list
  185. else
  186. codelist=>codelist%next !try next code in list
  187. end if
  188. end do
  189. ! count levels for that code
  190. if ( ih(2).ne.ll_first%level .and. ih(2).ne.0 .and. tdef==1 ) &
  191. & then
  192. !level information should be complete in
  193. !first time level so we are only looking there
  194. codelist%nlev=codelist%nlev+1 !code repeats at first~
  195. ! timelevel => must be new level and can't be zero
  196. levlist=>ll_first
  197. ! maintain list of levels
  198. do while (levlist%level.ne.ih(2))
  199. if (.not. associated(levlist%next)) then
  200. ll_help=>levlist
  201. allocate(levlist)
  202. ll_help%next=>levlist
  203. nullify(levlist%next)
  204. levlist%level=ih(2)
  205. nlevmax=nlevmax+1
  206. else
  207. levlist=>levlist%next !try next level in list
  208. end if
  209. end do
  210. end if
  211. end do
  212. 999 Continue
  213. if (irec==0) then
  214. print*,"*** No records processed. Input file was empty! ***"
  215. stop
  216. end if
  217. !
  218. ! write to descriptorfile levels and codes
  219. !
  220. if (l_vsect) then
  221. write(3,'(a,i4,a,i8,a)')'ZDEF ',nlat,' LINEAR ', &
  222. & abs(ll_first%level),' 1'
  223. else
  224. if (nlevmax==1) then
  225. write(3,'(a,i8,a)')'ZDEF 1 LINEAR ',ll_first%level,' 1'
  226. else
  227. write(3,'(a,i6,a)')'ZDEF ',nlevmax,' LEVELS '
  228. levlist=>ll_first
  229. do while (associated(levlist))
  230. write(3,'(i10)')levlist%level
  231. levlist=>levlist%next
  232. end do
  233. endif
  234. end if
  235. ! make a date
  236. ! fdate and ih(3) are of the form yyyymmdd
  237. ! extract first time:
  238. if (fhour.gt.99) then
  239. mi=mod(fhour,100)
  240. fhour=max(min(fhour/100,23),0)
  241. else
  242. fhour=max(min(fhour,23),0)
  243. mi=0
  244. end if
  245. ! extract first date
  246. mm=min(12,max(mod(fdate,10000)/100,1))
  247. dd=min(maxdpm(mm),max(mod(fdate,100),1))
  248. yyyy=max(1,min((fdate-100*mm-dd)/10000,9999)) !There is no year zero
  249. ! extract last time
  250. if (ih(4).gt.99) then
  251. mil=mod(ih(4),100)
  252. lhour=max(min(ih(4)/100,23),0)
  253. else
  254. lhour=max(min(ih(4),23),0)
  255. mil=0
  256. end if
  257. ! extract last date
  258. ddl=min(31,max(mod(ih(3),100),1))
  259. mml=min(12,max(mod(ih(3),10000)/100,1))
  260. yyyyl=max(1,min((ih(3)-100*mml-ddl)/10000,9999)) ! There is no year zero
  261. ! convert everything to minutes
  262. fdate=mi +60*(fhour+24*(dd +30*(mm +12*(yyyy )))) ! a modeler's month
  263. ldate=mil+60*(lhour+24*(ddl+30*(mml+12*(yyyyl)))) ! has 30 days :-)
  264. ! calculate time increment
  265. deltat=max(1,(ldate-fdate)/(max(tdef-1,1))) ! just let's be robust
  266. ! find appropriate unit for time increment
  267. ut='mn' ! fallback
  268. if (mod(deltat,60).eq.0) then
  269. deltat=deltat/60
  270. ut='hr'
  271. if (mod(deltat,24).eq.0) then
  272. deltat=deltat/24
  273. ut='dy'
  274. if (mod(deltat,30).eq.0) then
  275. deltat=deltat/30
  276. ut='mo'
  277. if (mod(deltat,12).eq.0) then
  278. deltat=deltat/12
  279. ut='yr'
  280. if (deltat.eq.0) then
  281. deltat=6
  282. ut='hr'
  283. print*,'Could not extract proper time increment.'
  284. print*,'Assuming 6-hourly data...'
  285. end if
  286. end if
  287. end if
  288. end if
  289. end if
  290. ! Was this complicated?
  291. ! We need a date string of the form hh:mmZddmmmyyyy
  292. write(3,'(a4,i6,a8,i2.2,a1,i2.2,a1,i2.2,a3,i4.4,i9,a2)') &
  293. & 'TDEF',tdef,' LINEAR ',fhour,':',mi,'Z',dd,month(mm), &
  294. & yyyy,deltat,ut
  295. ! write codelist
  296. write(3,'(a,i6)')'VARS ',nvar
  297. codelist=>cl_first
  298. do while(associated(codelist))
  299. write(vform(6:6),'(i1)') &
  300. & (min(9,1+int(log10(max(1.,real(codelist%code))))))!varying format
  301. if (l_vsect) codelist%nlev=nlat
  302. write(3,vform)'c',codelist%code, &
  303. & codelist%nlev,' 99 ',codelist%code,ih(7),ih(8)
  304. if (l_vsect) then
  305. print*,"level -1 detected, assuming vertical section with"
  306. print*,nlon,"latitudes and ",nlat," linear levels."
  307. else
  308. if (codelist%nlev.ne.0 .and. codelist%nlev.ne.nlevmax) then
  309. print*,'***CAUTION!*** code ',codelist%code, &
  310. & ' is not defined on all levels. '
  311. print*,'GrADS will assume that they are the first few', &
  312. & ' from the descriptor file: ',filedes
  313. end if
  314. end if
  315. codelist=>codelist%next
  316. end do
  317. write(3,'(a)')'ENDVARS'
  318. Write (*,*) 'Processed ', irec, ' record(s)'
  319. End
  320. subroutine usage
  321. implicit none
  322. print*,'Usage: srvctl <service file>'
  323. print*,'Example: srvctl z500.srv will create z500.ctl'
  324. stop
  325. end subroutine usage
  326. subroutine inigau(kl,pga)
  327. implicit none
  328. integer :: kl ! Number of Gaussian latitudes
  329. integer :: j0,j1 ! Loop index
  330. real :: pga(kl) ! Gaussian abscissas
  331. real (kind=8) :: pz0(kl) ! Gaussian abscissas
  332. real (kind=8), parameter :: pi = 3.14159265358979_8
  333. real (kind=8) :: z0,z1,z2,z3
  334. real (kind=8) :: ql,qld
  335. ! 1. Compute Gaussian abscissas
  336. z0 = pi / (2*kl+1)
  337. z1 = 1.0_8 / (kl*kl*8)
  338. do j0 = 1 , kl/2
  339. z2 = z0 * (2*j0 - 0.5_8)
  340. z2 = cos(z2 + z1 / tan(z2))
  341. do j1 = 1 , 50
  342. z3 = ql(kl,z2) / qld(kl,z2)
  343. z2 = z2 - z3
  344. if (abs(z3) < 1.0e-16) exit
  345. enddo ! j1
  346. ! print *,j1,' iter for lat= ',j0,' eps= ',z3
  347. pz0(j0) = z2
  348. pz0(kl-j0+1) = -z2
  349. enddo ! j0
  350. pga(:) = pz0(:) ! Double precision to native
  351. return
  352. end
  353. real (kind=8) function ql(k,p)
  354. implicit none
  355. integer :: k
  356. real (kind=8) :: p
  357. real (kind=8) :: z0,z1,z2,z3,z4
  358. integer :: j
  359. z0 = acos(p)
  360. z1 = 1.0
  361. z2 = 0.0
  362. do j = k , 0 , -2
  363. z3 = z1 * cos(z0 * j)
  364. z2 = z2 + z3
  365. z4 = (k-j+1) * (k+j) * 0.5
  366. z1 = z1 * z4 / (z4 + (j-1))
  367. enddo ! j
  368. if (mod(k,2) == 0) z2 = z2 - 0.5 * z3
  369. z0 = sqrt(2.0_8)
  370. do j = 1 ,k
  371. z0 = z0 * sqrt(1.0_8 - 0.25_8 / (j*j))
  372. enddo ! j
  373. ql = z0 * z2
  374. return
  375. end
  376. real (kind=8) function qld(k,p)
  377. implicit none
  378. integer :: k
  379. real (kind=8) :: p
  380. real (kind=8) :: z0,z1,z2,z3,z4
  381. real (kind=8) :: ql
  382. integer :: j
  383. if (k == 0) then
  384. qld = 0.0
  385. return
  386. endif
  387. z0 = 1.0 / (p*p - 1.0)
  388. z1 = k
  389. z2 = 2.0 * z1
  390. z3 = sqrt((z2 + 1.0) / (z2 - 1.0))
  391. z4 = p * ql(k,p) - z3 * ql(k-1,p)
  392. qld = z0 * z1 * z4
  393. return
  394. end