tmm_mf_ecmwf_mars.F90 62 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240124112421243124412451246124712481249125012511252125312541255125612571258125912601261126212631264126512661267126812691270127112721273127412751276127712781279128012811282128312841285128612871288128912901291129212931294129512961297129812991300130113021303130413051306130713081309131013111312131313141315131613171318131913201321132213231324132513261327132813291330133113321333133413351336133713381339134013411342134313441345134613471348134913501351135213531354135513561357135813591360136113621363136413651366136713681369137013711372137313741375137613771378137913801381138213831384138513861387138813891390139113921393139413951396139713981399140014011402140314041405140614071408140914101411141214131414141514161417141814191420142114221423142414251426142714281429143014311432143314341435143614371438143914401441144214431444144514461447144814491450145114521453145414551456145714581459146014611462146314641465146614671468146914701471147214731474147514761477147814791480148114821483148414851486148714881489149014911492149314941495149614971498149915001501150215031504150515061507150815091510151115121513151415151516151715181519152015211522152315241525152615271528152915301531153215331534153515361537153815391540154115421543154415451546154715481549155015511552155315541555155615571558155915601561156215631564156515661567156815691570157115721573157415751576157715781579158015811582158315841585158615871588158915901591159215931594159515961597159815991600160116021603160416051606160716081609161016111612161316141615161616171618161916201621162216231624162516261627162816291630163116321633163416351636163716381639164016411642164316441645164616471648164916501651165216531654165516561657165816591660166116621663166416651666166716681669167016711672167316741675167616771678167916801681168216831684168516861687168816891690169116921693169416951696169716981699170017011702170317041705170617071708170917101711171217131714171517161717171817191720172117221723172417251726172717281729173017311732173317341735173617371738173917401741174217431744174517461747174817491750175117521753175417551756175717581759176017611762176317641765176617671768176917701771177217731774177517761777177817791780178117821783178417851786178717881789179017911792179317941795179617971798179918001801180218031804180518061807180818091810181118121813181418151816181718181819182018211822182318241825182618271828182918301831183218331834183518361837183818391840184118421843184418451846184718481849185018511852185318541855185618571858185918601861186218631864
  1. !###############################################################################
  2. !
  3. ! Input/output of meteofiles : grib version.
  4. !
  5. !### macro's ###################################################################
  6. !
  7. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  8. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  9. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  10. !
  11. #include "tmm.inc"
  12. !
  13. !###############################################################################
  14. module tmm_mf_ecmwf_mars
  15. use GO , only : gol, goErr, goPr
  16. use GO , only : TDate
  17. implicit none
  18. ! --- in/out ----------------------------
  19. private
  20. public :: TMeteoFile_ecmwf_mars
  21. public :: Init, Done
  22. public :: Get
  23. public :: ReadRecord
  24. ! --- const ------------------------------
  25. character(len=*), parameter :: mname = 'tmm_mf_ecmwf_mars'
  26. !--- type ---------------------------------
  27. type TMeteoFile_ecmwf_mars
  28. ! file name:
  29. character(len=256) :: dir
  30. character(len=256) :: fname
  31. ! time range covered by file:
  32. type(TDate) :: trange(2)
  33. ! other time keys for this file:
  34. character(len=16) :: treskey
  35. ! current time range covered by grib record:
  36. type(TDate) :: tref, t1, t2
  37. !
  38. ! file description
  39. !
  40. character(len=16) :: ec_class, ec_type, ec_levs
  41. integer :: ec_sh, ec_gg
  42. character(len=256) :: paramkeys
  43. !
  44. ! adhoc flag:
  45. logical :: slhf_for_convec
  46. end type TMeteoFile_ecmwf_mars
  47. ! --- interfaces -------------------------
  48. interface Init
  49. module procedure mf_Init
  50. end interface
  51. interface Done
  52. module procedure mf_Done
  53. end interface
  54. interface Get
  55. module procedure mf_Get
  56. end interface
  57. interface ReadRecord
  58. module procedure mf_ReadRecord
  59. end interface
  60. contains
  61. ! ==============================================================
  62. subroutine mf_Init( mf, dir, archivekeys, paramkey, &
  63. tday, t1, t2, status )
  64. use GO, only : TDate
  65. use GO, only : goVarValue
  66. ! --- in/out ----------------------------
  67. type(TMeteoFile_ecmwf_mars), intent(out) :: mf
  68. character(len=*), intent(in) :: dir
  69. character(len=*), intent(in) :: archivekeys
  70. character(len=*), intent(in) :: paramkey
  71. type(TDate), intent(in) :: tday, t1, t2
  72. integer, intent(out) :: status
  73. ! --- const --------------------------------------
  74. character(len=*), parameter :: rname = mname//'/mf_Init'
  75. ! --- local --------------------------------
  76. ! --- begin --------------------------------
  77. ! store
  78. mf%dir = dir
  79. !
  80. ! extract fields from archivekey :
  81. ! form=tmpp;class=od;type=fg;levs=ml60;sh=159;gg=80;tres=_fg006up4tr3
  82. !
  83. mf%ec_class = 'od'
  84. call goVarValue( archivekeys, ';', 'class', '=', mf%ec_class, status )
  85. IF_ERROR_RETURN(status=1)
  86. !
  87. mf%ec_type = 'fc'
  88. call goVarValue( archivekeys, ';', 'type', '=', mf%ec_type, status )
  89. IF_ERROR_RETURN(status=1)
  90. !
  91. mf%ec_levs = 'ml60'
  92. call goVarValue( archivekeys, ';', 'levs', '=', mf%ec_levs, status )
  93. IF_ERROR_RETURN(status=1)
  94. !
  95. mf%ec_sh = 159
  96. call goVarValue( archivekeys, ';', 'sh', '=', mf%ec_sh, status )
  97. IF_ERROR_RETURN(status=1)
  98. !
  99. mf%ec_gg = 80
  100. call goVarValue( archivekeys, ';', 'gg', '=', mf%ec_gg, status )
  101. IF_ERROR_RETURN(status=1)
  102. !
  103. mf%treskey = '_fc012up2tr3'
  104. call goVarValue( archivekeys, ';', 'tres', '=', mf%treskey, status )
  105. IF_ERROR_RETURN(status=1)
  106. ! adhoc settings immplemented for surface heat fluxes for convection; enable ?
  107. mf%slhf_for_convec = .false.
  108. call goVarValue( archivekeys, ';', 'slhf_for_convec', '=', mf%slhf_for_convec, status )
  109. IF_ERROR_RETURN(status=1)
  110. ! specials
  111. select case ( paramkey )
  112. case ( 'oro', 'lsm' )
  113. ! overwrite timeresolutionkey, used later on to set trange
  114. mf%treskey = 'const'
  115. ! tmm_convec tries to read oro using the default sourcekey,
  116. ! which probably contains type=fc ; force to use an for oro ...
  117. mf%ec_type = 'an'
  118. end select
  119. ! single parameter in a file:
  120. mf%paramkeys = '-'//trim(paramkey)//'-'
  121. ! extract time range:
  122. call GetGribTime( mf, mf%treskey, tday, t1, t2, status, trange=mf%trange )
  123. IF_NOTOK_RETURN(status=1)
  124. ! dummy filename, might be used in error diagnose
  125. write (mf%fname,'("ecmwf mars grib file for param ",a)') trim(paramkey)
  126. ! ok
  127. status = 0
  128. end subroutine mf_Init
  129. ! ***
  130. subroutine mf_Done( mf, status )
  131. ! --- in/out ------------------------------------
  132. type(TMeteoFile_ecmwf_mars), intent(inout) :: mf
  133. integer, intent(out) :: status
  134. ! --- const --------------------------------------
  135. character(len=*), parameter :: rname = mname//'/mf_Done'
  136. ! --- begin -------------------------------------
  137. ! nothing to be done ...
  138. ! ok
  139. status = 0
  140. end subroutine mf_Done
  141. ! ***
  142. subroutine mf_Get( mf, status, trange1, trange2, paramkeys )
  143. use GO, only : TDate
  144. ! --- in/out ----------------------------
  145. type(TMeteoFile_ecmwf_mars), intent(in) :: mf
  146. integer, intent(out) :: status
  147. type(TDate), intent(out), optional :: trange1, trange2
  148. character(len=*), intent(out), optional :: paramkeys
  149. ! --- const --------------------------------------
  150. character(len=*), parameter :: rname = mname//'/mf_Get'
  151. ! --- local --------------------------------
  152. ! --- begin --------------------------------
  153. ! time range:
  154. if ( present(trange1) ) trange1 = mf%trange(1)
  155. if ( present(trange2) ) trange2 = mf%trange(2)
  156. ! contents:
  157. if ( present(paramkeys) ) paramkeys = trim(mf%paramkeys)
  158. ! ok
  159. status = 0
  160. end subroutine mf_Get
  161. ! ***
  162. subroutine mf_ReadRecord( mf, paramkey, tday, t1, t2, nuv, nw, &
  163. gridtype, levi, &
  164. lli, ll, sp_ll, &
  165. ggi, gg, sp_gg, &
  166. shi, sh, lnsp_sh, &
  167. tmi, status )
  168. use parray , only : pa_Init, pa_Done
  169. use GO , only : TDate
  170. use Grid , only : TLevelInfo
  171. use Grid , only : TllGridInfo, TggGridInfo, TshGridInfo
  172. use tmm_info , only : TMeteoInfo
  173. ! --- in/out -------------------------------
  174. type(TMeteoFile_ecmwf_mars), intent(inout) :: mf
  175. character(len=*), intent(in) :: paramkey
  176. type(TDate), intent(in) :: tday, t1, t2
  177. character(len=1), intent(in) :: nuv
  178. character(len=1), intent(in) :: nw
  179. character(len=2), intent(out) :: gridtype
  180. type(TLevelInfo), intent(out) :: levi
  181. type(TllGridInfo), intent(inout) :: lli
  182. real, pointer :: ll(:,:,:)
  183. real, pointer :: sp_ll(:,:)
  184. type(TggGridInfo), intent(inout) :: ggi
  185. real, pointer :: gg(:,:)
  186. real, pointer :: sp_gg(:)
  187. type(TshGridInfo), intent(inout) :: shi
  188. complex, pointer :: sh(:,:)
  189. complex, pointer :: lnsp_sh(:)
  190. type(TMeteoInfo), intent(out) :: tmi
  191. integer, intent(out) :: status
  192. ! --- const --------------------------------------
  193. character(len=*), parameter :: rname = mname//'/mf_ReadRecord'
  194. ! --- local -------------------------------
  195. real, pointer :: ll2(:,:,:)
  196. real, pointer :: gg2(:,:)
  197. complex, pointer :: sh2(:,:)
  198. integer :: iveg
  199. logical :: unit_change
  200. real :: unit_fac
  201. ! --- begin ---------------------------------
  202. ! combined field ?
  203. select case ( paramkey )
  204. ! *** surface stress
  205. case ( 'sstr' )
  206. ! read first field:
  207. call mf_ReadRecord_1( mf, 'ewss', tday, t1, t2, nuv, nw, &
  208. gridtype, levi, &
  209. lli, ll, sp_ll, &
  210. ggi, gg, sp_gg, &
  211. shi, sh, lnsp_sh, &
  212. tmi, status )
  213. IF_NOTOK_RETURN(status=1)
  214. ! init pointer:
  215. call pa_Init( ll2 )
  216. call pa_Init( gg2 )
  217. call pa_Init( sh2 )
  218. ! read second field:
  219. call mf_ReadRecord_1( mf, 'nsss', tday, t1, t2, nuv, nw, &
  220. gridtype, levi, &
  221. lli, ll2, sp_ll, &
  222. ggi, gg2, sp_gg, &
  223. shi, sh2, lnsp_sh, &
  224. tmi, status )
  225. IF_NOTOK_RETURN(status=1)
  226. ! process:
  227. select case ( gridtype )
  228. case ( 'll' ) ; ll = sqrt( ll**2 + ll2**2 )
  229. case ( 'gg' ) ; gg = sqrt( gg**2 + gg2**2 )
  230. case ( 'sh' ) ; sh = sqrt( sh**2 + sh2**2 )
  231. case default
  232. write (gol,'("unsupported gridtype for surface stress :",a)') gridtype; call goErr
  233. TRACEBACK; status=1; return
  234. end select
  235. ! clear pointers:
  236. call pa_Done( ll2 )
  237. call pa_Done( gg2 )
  238. call pa_Done( sh2 )
  239. ! *** vegetation types
  240. case ( 'tv01', 'tv02', 'tv03', 'tv04', 'tv05', 'tv06', 'tv07', 'tv08', 'tv09', 'tv10', &
  241. 'tv11', 'tv12', 'tv13', 'tv14', 'tv15', 'tv16', 'tv17', 'tv18', 'tv19', 'tv20' )
  242. ! extract number from name
  243. read (paramkey(3:4),'(i2)') iveg
  244. ! low vegetation types:
  245. call mf_ReadRecord_1( mf, 'tvl', tday, t1, t2, nuv, nw, &
  246. gridtype, levi, &
  247. lli, ll, sp_ll, &
  248. ggi, gg, sp_gg, &
  249. shi, sh, lnsp_sh, &
  250. tmi, status )
  251. IF_NOTOK_RETURN(status=1)
  252. ! set elements that match requested vegetation type to 100%, zero elsewhere
  253. select case ( gridtype )
  254. case ( 'gg' )
  255. where ( nint(gg(:,1)) == iveg )
  256. gg(:,1) = 100.0 ! %
  257. elsewhere
  258. gg(:,1) = 0.0
  259. end where
  260. case default
  261. write (gol,'("unsupported gridtype for vegetation fractions :",a)') gridtype; call goErr
  262. TRACEBACK; status=1; return
  263. end select
  264. ! init pointer:
  265. call pa_Init( ll2 )
  266. call pa_Init( gg2 )
  267. call pa_Init( sh2 )
  268. ! high vegetation types:
  269. call mf_ReadRecord_1( mf, 'tvh', tday, t1, t2, nuv, nw, &
  270. gridtype, levi, &
  271. lli, ll2, sp_ll, &
  272. ggi, gg2, sp_gg, &
  273. shi, sh2, lnsp_sh, &
  274. tmi, status )
  275. IF_NOTOK_RETURN(status=1)
  276. ! set elements that match requested vegetation type to 100%:
  277. select case ( gridtype )
  278. case ( 'gg' )
  279. where ( nint(gg2(:,1)) == iveg )
  280. gg(:,1) = 100.0 ! %
  281. end where
  282. case default
  283. write (gol,'("unsupported gridtype for vegetation fractions :",a)') gridtype; call goErr
  284. TRACEBACK; status=1; return
  285. end select
  286. ! clear pointers:
  287. call pa_Done( ll2 )
  288. call pa_Done( gg2 )
  289. call pa_Done( sh2 )
  290. ! *** default
  291. case default
  292. call mf_ReadRecord_1( mf, paramkey, tday, t1, t2, nuv, nw, &
  293. gridtype, levi, &
  294. lli, ll, sp_ll, &
  295. ggi, gg, sp_gg, &
  296. shi, sh, lnsp_sh, &
  297. tmi, status )
  298. IF_NOTOK_RETURN(status=1)
  299. end select
  300. ! unit change ?
  301. unit_change = .true.
  302. select case ( paramkey )
  303. case ( 'lsm' ) ; unit_fac = 100.0 ! 0-1 -> 0-100%
  304. case default ; unit_change = .false.
  305. end select
  306. ! apply ?
  307. if ( unit_change ) then
  308. select case ( gridtype )
  309. case ( 'll' ) ; ll = ll * unit_fac
  310. case ( 'gg' ) ; gg = gg * unit_fac
  311. !case ( 'sh' ) ; sh = sh * unit_fac
  312. case default
  313. write (gol,'("unsupported gridtype for unit change :",a)') gridtype; call goErr
  314. TRACEBACK; status=1; return
  315. end select
  316. end if
  317. ! ok
  318. status = 0
  319. end subroutine mf_ReadRecord
  320. ! ***
  321. subroutine mf_ReadRecord_1( mf, paramkey, tday, t1, t2, nuv, nw, &
  322. gridtype, levi, &
  323. lli, ll, sp_ll, &
  324. ggi, gg, sp_gg, &
  325. shi, sh, lnsp_sh, &
  326. tmi, status )
  327. use parray , only : pa_Init, pa_Done
  328. use GO , only : TDate, TIncrDate, operator(<), operator(-), rTotal, wrtgol
  329. use Grid , only : TLevelInfo
  330. use Grid , only : TllGridInfo, TggGridInfo, TshGridInfo
  331. use tmm_info , only : TMeteoInfo
  332. ! --- in/out -------------------------------
  333. type(TMeteoFile_ecmwf_mars), intent(inout) :: mf
  334. character(len=*), intent(in) :: paramkey
  335. type(TDate), intent(in) :: tday, t1, t2
  336. character(len=1), intent(in) :: nuv
  337. character(len=1), intent(in) :: nw
  338. character(len=2), intent(out) :: gridtype
  339. type(TLevelInfo), intent(out) :: levi
  340. type(TllGridInfo), intent(inout) :: lli
  341. real, pointer :: ll(:,:,:)
  342. real, pointer :: sp_ll(:,:)
  343. type(TggGridInfo), intent(inout) :: ggi
  344. real, pointer :: gg(:,:)
  345. real, pointer :: sp_gg(:)
  346. type(TshGridInfo), intent(inout) :: shi
  347. complex, pointer :: sh(:,:)
  348. complex, pointer :: lnsp_sh(:)
  349. type(TMeteoInfo), intent(out) :: tmi
  350. integer, intent(out) :: status
  351. ! --- const --------------------------------------
  352. character(len=*), parameter :: rname = mname//'/mf_ReadRecord_1'
  353. ! --- local -------------------------------
  354. type(TDate) :: tref
  355. type(TIncrDate) :: tshift
  356. type(TDate) :: trefs, t1s, t2s
  357. real, pointer :: ll1(:,:,:)
  358. real, pointer :: gg1(:,:)
  359. complex, pointer :: sh1(:,:)
  360. real :: dt_sec
  361. ! --- begin ---------------------------------
  362. !write (gol,'("mf_ReadRecord_1: paramkey : ",a)') trim(paramkey); call goPr
  363. !call wrtgol( 'mf_ReadRecord_1: tday : ', tday ); call goPr
  364. !call wrtgol( 'mf_ReadRecord_1: t1 : ', t1 ); call goPr
  365. !call wrtgol( 'mf_ReadRecord_1: t2 : ', t2 ); call goPr
  366. ! accumulated field ?
  367. select case ( paramkey )
  368. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  369. ! accumulated fields
  370. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  371. case ( 'lsp', 'cp', 'sf', 'sshf', 'slhf', &
  372. 'ssr', 'ssrd', 'str', 'strd', &
  373. 'ewss', 'nsss', &
  374. 'UDMF', 'UDDR', 'DDMF', 'DDDR' )
  375. ! get reference time for requested time interval:
  376. call GetGribTime( mf, mf%treskey, tday, t1, t2, status, tref=tref )
  377. IF_NOTOK_RETURN(status=1)
  378. ! should be a time interval ...
  379. if ( .not. (t1 < t2) ) then
  380. write (gol,'("accumulated fields requires time interval:")'); call goErr
  381. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  382. call wrtgol( ' t1 : ', t1 ); call goErr
  383. call wrtgol( ' t2 : ', t2 ); call goErr
  384. TRACEBACK; status=1; return
  385. end if
  386. ! read field accumulated over [tref,t2] :
  387. call wrtgol( ' accum ', tref, ' - ', t2 ); call goPr
  388. call mf_ReadRecord_2( mf, paramkey, tref, t2, t2, nuv, nw, &
  389. gridtype, levi, &
  390. lli, ll, sp_ll, &
  391. ggi, gg, sp_gg, &
  392. shi, sh, lnsp_sh, &
  393. tmi, status )
  394. IF_NOTOK_RETURN(status=1)
  395. ! substract [tref,t1] if necessary:
  396. if ( tref < t1 ) then
  397. ! init pointer:
  398. call pa_Init( ll1 )
  399. call pa_Init( gg1 )
  400. call pa_Init( sh1 )
  401. ! read field accumulated over [tref,t1] :
  402. call wrtgol( ' accum ', tref, ' - ', t1 ); call goPr
  403. call mf_ReadRecord_2( mf, paramkey, tref, t1, t1, nuv, nw, &
  404. gridtype, levi, &
  405. lli, ll1, sp_ll, &
  406. ggi, gg1, sp_gg, &
  407. shi, sh1, lnsp_sh, &
  408. tmi, status )
  409. IF_NOTOK_RETURN(status=1)
  410. ! substract:
  411. select case ( gridtype )
  412. case ( 'll' ) ; ll = ll - ll1
  413. case ( 'gg' ) ; gg = gg - gg1
  414. case ( 'sh' ) ; sh = sh - sh1
  415. case default
  416. write (gol,'("unsupported gridtype for substract :",a)') gridtype; call goErr
  417. TRACEBACK; status=1; return
  418. end select
  419. ! clear pointers:
  420. call pa_Done( ll1 )
  421. call pa_Done( gg1 )
  422. call pa_Done( sh1 )
  423. end if
  424. ! return time averages only:
  425. dt_sec = rTotal( t2 - t1, 'sec' )
  426. select case ( gridtype )
  427. case ( 'll' ) ; ll = ll / dt_sec
  428. case ( 'gg' ) ; gg = gg / dt_sec
  429. case ( 'sh' ) ; sh = sh / dt_sec
  430. case default
  431. write (gol,'("unsupported gridtype for time average :",a)') gridtype; call goErr
  432. TRACEBACK; status=1; return
  433. end select
  434. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  435. ! instantaneous fields
  436. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  437. case default
  438. ! get reference time for requested time interval;
  439. ! eventually shift for analysed fields in case of forecasts:
  440. call GetGribTime( mf, mf%treskey, tday, t1, t2, status, tref=tref, tshift=tshift )
  441. IF_NOTOK_RETURN(status=1)
  442. ! shift times (might be zero):
  443. trefs = tref - tshift
  444. t1s = t1 - tshift
  445. t2s = t2 - tshift
  446. ! just read ..
  447. call mf_ReadRecord_2( mf, paramkey, trefs, t1s, t2s, nuv, nw, &
  448. gridtype, levi, &
  449. lli, ll, sp_ll, &
  450. ggi, gg, sp_gg, &
  451. shi, sh, lnsp_sh, &
  452. tmi, status )
  453. IF_NOTOK_RETURN(status=1)
  454. end select
  455. ! ok
  456. status = 0
  457. end subroutine mf_ReadRecord_1
  458. ! ***
  459. subroutine mf_ReadRecord_2( mf, paramkey, tref, t1, t2, nuv, nw, &
  460. gridtype, levi, &
  461. lli, ll, sp_ll, &
  462. ggi, gg, sp_gg, &
  463. shi, sh, lnsp_sh, &
  464. tmi, status )
  465. use GO , only : TDate, wrtgol, Get, NewDate, operator(>)
  466. use GO , only : goWriteKeyNum
  467. use Grid , only : Init, Done
  468. use Grid , only : TLevelInfo
  469. use Grid , only : TllGridInfo, TggGridInfo, TshGridInfo
  470. use Grid , only : Interpol
  471. use file_grib , only : TGribFile
  472. use file_grib , only : Init, Done, ReadRecord, Get, Check
  473. use file_grib , only : levtype_sfc, levtype_hyb, levtype_land
  474. use file_grib , only : gridtype_ll, gridtype_gg, gridtype_sh
  475. use grib_table, only : GetPid
  476. use PArray , only : pa_Init, pa_Done, pa_SetShape
  477. use tmm_info , only : TMeteoInfo, Init, AddHistory
  478. ! --- in/out -------------------------------
  479. type(TMeteoFile_ecmwf_mars), intent(inout) :: mf
  480. character(len=*), intent(in) :: paramkey
  481. type(TDate), intent(in) :: tref, t1, t2
  482. character(len=1), intent(in) :: nuv
  483. character(len=1), intent(in) :: nw
  484. character(len=2), intent(out) :: gridtype
  485. type(TLevelInfo), intent(out) :: levi
  486. type(TllGridInfo), intent(inout) :: lli
  487. real, pointer :: ll(:,:,:)
  488. real, pointer :: sp_ll(:,:)
  489. type(TggGridInfo), intent(inout) :: ggi
  490. real, pointer :: gg(:,:)
  491. real, pointer :: sp_gg(:)
  492. type(TshGridInfo), intent(inout) :: shi
  493. complex, pointer :: sh(:,:)
  494. complex, pointer :: lnsp_sh(:)
  495. type(TMeteoInfo), intent(out) :: tmi
  496. integer, intent(out) :: status
  497. ! --- const --------------------------------------
  498. character(len=*), parameter :: rname = mname//'/mf_ReadRecord_2'
  499. ! --- local -------------------------------
  500. character(len=16) :: ec_class, ec_type
  501. character(len=16) :: ec_grid, gridN, gridT
  502. character(len=16) :: levs
  503. character(len=16) :: treskey
  504. logical :: constant
  505. type(TGribFile) :: grib
  506. logical :: do_spm
  507. character(len=256) :: spm_fname
  508. type(TGribFile) :: spm_grib
  509. logical :: spm_lnsp
  510. logical :: spm_lnsp2sp
  511. integer :: pid
  512. character(len=7) :: gribcode
  513. character(len=16) :: spm_levs, spm_paramkey, ec_paramkey
  514. type(TDate) :: tfile
  515. integer :: ccyy, mm, dd, hh, mn
  516. type(TDate) :: tc
  517. logical :: exist
  518. logical :: isfirst
  519. logical :: reopened
  520. integer :: pid, tabid
  521. integer :: nlev, glevtype, glevel
  522. integer :: level
  523. integer :: nlev_out, level_out
  524. integer :: ggridtype
  525. real :: lon_first, lon_inc
  526. integer :: lon_n
  527. real :: lat_first, lat_inc
  528. integer :: lat_n
  529. integer :: ggN
  530. integer :: shT
  531. integer :: greftime(5), gtimerange(4)
  532. character(len=64) :: key
  533. integer :: ilat
  534. real, pointer :: pat(:,:)
  535. type(TshGridInfo) :: tmp_shi
  536. complex, pointer :: tmp_sh(:)
  537. ! --- begin ---------------------------------
  538. !write (gol,'("mf_ReadRecord_2: paramkey : ",a)') trim(paramkey); call goPr
  539. !call wrtgol( 'mf_ReadRecord_2: tref : ', tref ); call goPr
  540. !call wrtgol( 'mf_ReadRecord_2: t1 : ', t1 ); call goPr
  541. !call wrtgol( 'mf_ReadRecord_2: t2 : ', t2 ); call goPr
  542. ! no fluxes through boundaries ...
  543. if ( nuv /= 'n' ) then
  544. write (gol,'("unsupported nuv key : ",a)') nuv; call goErr
  545. TRACEBACK; status=1; return
  546. end if
  547. ! limitted support of fluxes ..
  548. !if ( nw /= 'n' ) then
  549. ! write (gol,'("unsupported nw key : ",a)') nw; call goErr
  550. ! TRACEBACK; status=1; return
  551. !end if
  552. ! init pointer arrays:
  553. call pa_Init( pat )
  554. !
  555. ! ~~~ 3d field or 2d stored in 3d array
  556. !
  557. ! grid : T159, N80, etc
  558. call goWriteKeyNum( gridT, 'T', mf%ec_sh )
  559. call goWriteKeyNum( gridN, 'N', mf%ec_gg )
  560. ! defaults
  561. ec_paramkey = paramkey
  562. ec_class = mf%ec_class
  563. ec_type = mf%ec_type
  564. levs = mf%ec_levs
  565. ec_grid = gridN
  566. treskey = mf%treskey
  567. constant = .false.
  568. do_spm = .false.
  569. spm_lnsp = .false.
  570. spm_lnsp2sp = .false.
  571. ! specials
  572. select case ( paramkey )
  573. case ( 'LNSP' )
  574. ec_grid = gridT
  575. case ( 'VO', 'D' )
  576. ec_grid = gridT
  577. do_spm = .true.
  578. spm_lnsp = .true.
  579. case ( 'T', 'W', 'Q', 'CLWC', 'CIWC', 'CC', 'UDMF', 'UDDR', 'DDMF', 'DDDR' )
  580. do_spm = .true.
  581. spm_lnsp2sp = (ec_class == 'e4') .or. (ec_class == 'ei')
  582. case ( 'oro', 'lsm' )
  583. levs = 'sfc'
  584. constant = .true.
  585. treskey = 'const'
  586. ! tmm_convec tries to read oro using the default sourcekey,
  587. ! which probably contains type=fc ; force to use an for oro ...
  588. ec_type = 'an'
  589. ! case ( 'cvl', 'cvh', 'tvl', 'tvh', 'sr', 'albedo', 'lsrh' )
  590. ! ec_type = 'an'
  591. ! levs = 'sfc'
  592. case ( 'ci', 'sst', 'swvl1', 'swvl2', 'swvl3', 'swvl4', '10fg', 'sd', 'lsp', &
  593. 'cp', 'sf', 'sshf', 'slhf', 'blh', 'u10m', 'v10m', 't2m', 'd2m', &
  594. 'ssr', 'ewss', 'nsss', 'sstr' ,'src', 'skt',
  595. 'stl1', 'stl2', 'stl3', 'stl4' )
  596. levs = 'sfc'
  597. case ( 'sp' )
  598. spm_lnsp2sp = (ec_class == 'e4') .or. (ec_class == 'ei')
  599. if ( spm_lnsp2sp ) then
  600. ec_paramkey = 'LNSP'
  601. else
  602. levs = 'sfc'
  603. end if
  604. end select
  605. ! write gribcode
  606. call GetPid( 'ec', ec_paramkey, pid, tabid, status )
  607. IF_NOTOK_RETURN(status=1)
  608. write (gribcode,'(i3,".",i3.3)') pid, tabid
  609. gribcode = adjustl(gribcode)
  610. ! convert input times to file name times:
  611. call GetGribTime( mf, treskey, tref, t1, t2, status, tfile=tfile )
  612. IF_NOTOK_RETURN(status=1)
  613. !call wrtgol( 'mf_ReadRecord_2: tfile : ', tfile ); call goPr
  614. ! extract time values:
  615. call Get( tfile, year=ccyy, month=mm, day=dd, hour=hh, min=mn )
  616. ! create file name:
  617. ! dir/od-fc-20000101-1200-ml60-138-T159.gb
  618. !
  619. ! filename includes date:
  620. write (mf%fname,'(a,"/",a,"-",a,"-",i4.4,2i2.2,"-",2i2.2,"-",a,"-",a,"-",a,".gb")') &
  621. trim(mf%dir), &
  622. trim(ec_class), trim(ec_type), ccyy, mm, dd, hh, mn, &
  623. trim(levs), trim(gribcode), trim(ec_grid)
  624. ! exist ?
  625. inquire( file=mf%fname, exist=exist )
  626. if ( .not. exist ) then
  627. write (gol,'("grib file does not exist:")'); call goErr
  628. write (gol,'(" ",a)') trim(mf%fname); call goErr
  629. TRACEBACK; status=1; return
  630. end if
  631. ! open grib file
  632. call Init( grib, mf%fname, 'r', status )
  633. IF_NOTOK_RETURN(status=1)
  634. ! arrays and grids not defined yet
  635. isfirst = .true.
  636. reopened = .false.
  637. ! loop over records
  638. level = 0
  639. do
  640. !
  641. ! read gribsection in file buffer
  642. !
  643. call ReadRecord( grib, status )
  644. select case ( status )
  645. case ( 0 )
  646. ! no error
  647. case ( 1 )
  648. ! eof
  649. if ( .not. reopened ) then
  650. !write (gol,'("grib read record: re-open ...")'); call goPr
  651. ! close:
  652. call Done( grib, status )
  653. IF_NOTOK_RETURN(status=1)
  654. ! reopen:
  655. call Init( grib, mf%fname, 'r', status )
  656. IF_NOTOK_RETURN(status=1)
  657. reopened = .true.
  658. cycle
  659. else
  660. write (gol,'("reached eof before requested record was found")'); call goErr
  661. write (gol,'(" file : ",a)') trim(mf%fname); call goErr
  662. write (gol,'(" paramkey : ",a)') trim(paramkey); call goErr
  663. call wrtgol( ' tref : ', tref ); call goErr
  664. call wrtgol( ' t1 : ', t1 ); call goErr
  665. call wrtgol( ' t2 : ', t2 ); call goErr
  666. write (gol,'("tips:")'); call goErr
  667. write (gol,'(" o grib file corrupted or zero ?")'); call goErr
  668. write (gol,'(" o if accumulatd field,")'); call goErr
  669. write (gol,'(" check list of accumulated fields in mf_ReadRecord_1")'); call goErr
  670. TRACEBACK; status=1; return
  671. end if
  672. case default
  673. write (gol,'("error from grib ReadRecord; status=",i6)') status; call goErr
  674. TRACEBACK; status=1; return
  675. end select
  676. !
  677. ! checks ...
  678. !
  679. ! get param id for the requested field from grib table:
  680. select case ( paramkey )
  681. case ( 'spm' ) ; call GetPid( 'ec', 'SP' , pid, tabid, status )
  682. case default ; call GetPid( 'ec', ec_paramkey, pid, tabid, status )
  683. end select
  684. IF_NOTOK_RETURN(status=1)
  685. ! check parameter; continue if not ok:
  686. call Check( grib, pid=pid, debug=0, status=status )
  687. if (status/=0) cycle
  688. ! fill times ?
  689. if ( .not. constant ) then
  690. ! extract time fields from grib, store in mf%tref/mf%t1/mf%t2
  691. call SetTime( mf, grib, status )
  692. IF_NOTOK_RETURN(status=1)
  693. ! check time:
  694. call CheckTime( mf, tref, t1, t2, status )
  695. if (status/=0) then
  696. !write (gol,'("grib read record: wrong time; skip ...")'); call goPr
  697. !write (gol,'(" paramkey : ",a)') paramkey; call goPr
  698. !call wrtgol( ' req. tref : ', tref ); call goPr
  699. !call wrtgol( ' t1 : ', t1 ); call goPr
  700. !call wrtgol( ' t2 : ', t2 ); call goPr
  701. !call wrtgol( ' grib tref : ', mf%tref ); call goPr
  702. !call wrtgol( ' t1 : ', mf%t1 ); call goPr
  703. !call wrtgol( ' t2 : ', mf%t2 ); call goPr
  704. !write (gol,'(" grib file : ",a)') trim(mf%fname); call goPr
  705. cycle
  706. end if
  707. end if ! time checking
  708. ! extract level stuff:
  709. call Get( grib, nlev=nlev, levtype=glevtype, level=glevel, status=status )
  710. IF_NOTOK_RETURN(status=1)
  711. ! check level type:
  712. select case ( glevtype )
  713. case ( levtype_sfc, levtype_land )
  714. ! surface field
  715. nlev = 1
  716. glevel = 1
  717. case ( levtype_hyb )
  718. select case ( paramkey )
  719. case ( 'LNSP', 'sp' )
  720. nlev = 1
  721. case default
  722. ! level in 3d field
  723. end select
  724. case default
  725. write (gol,'("found unexpected level type: ")'); call goErr
  726. write (gol,'(" leveltype : ",i3)') glevtype; call goErr
  727. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  728. write (gol,'(" grib file : ",a)') trim(mf%fname); call goErr
  729. TRACEBACK; status=1; return
  730. end select
  731. ! check level:
  732. level = level + 1
  733. if ( glevel /= level ) then
  734. write (gol,'("found unexpected level: ")'); call goErr
  735. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  736. write (gol,'(" req. level : ",i6)') level; call goErr
  737. write (gol,'(" grib level : ",i6)') glevel; call goErr
  738. write (gol,'(" grib nlev : ",i6)') nlev; call goErr
  739. write (gol,'(" grib file : ",a)') trim(mf%fname); call goErr
  740. TRACEBACK; status=1; return
  741. end if
  742. ! number of output levels:
  743. nlev_out = nlev
  744. if ( nw == 'w' ) nlev_out = nlev_out + 1
  745. !
  746. ! define grids and arrays
  747. !
  748. if ( isfirst ) then
  749. !
  750. ! * info
  751. !
  752. ! example of history:
  753. ! model=ecmwf;class=od;type=fc;tref=2000,12,31,12,00; ...
  754. ! trange=001,234,240,001;sh=159;nlev=60
  755. !
  756. call Init( tmi, paramkey, 'unknown', status )
  757. call AddHistory( tmi, 'model==ecmwf', status )
  758. call AddHistory( tmi, 'class=='//trim(mf%ec_class), status )
  759. call AddHistory( tmi, 'type=='//trim(mf%ec_type) , status )
  760. !
  761. call Get( grib, status, reftime=greftime, timerange=gtimerange )
  762. IF_NOTOK_RETURN(status=1)
  763. write (key,'("tref==",i4.4,4(",",i2.2))') greftime
  764. call AddHistory( tmi, trim(key), status )
  765. write (key,'("trange==",i3.3,3(",",i3.3))') gtimerange
  766. call AddHistory( tmi, trim(key), status )
  767. !
  768. write (key,'("nlev==",i3.3)') nlev
  769. call AddHistory( tmi, trim(key), status )
  770. !
  771. ! * define horizontal grid:
  772. !
  773. ! extract grid type:
  774. call Get( grib, status, gridtype=ggridtype )
  775. IF_NOTOK_RETURN(status=1)
  776. ! setup:
  777. select case ( ggridtype )
  778. ! o lat/lon
  779. case ( gridtype_ll )
  780. ! routine returns lat/lon grid:
  781. gridtype = 'll'
  782. ! grib storage is north pole to south pole:
  783. call Get( grib, status, &
  784. lon_first=lon_first, lon_inc=lon_inc, lon_n=lon_n, &
  785. lat_last =lat_first, lat_inc=lat_inc, lat_n=lat_n )
  786. IF_NOTOK_RETURN(status=1)
  787. ! define grid structure:
  788. call Init( lli, lon_first, lon_inc, lon_n, &
  789. lat_first, lat_inc, lat_n, status )
  790. IF_NOTOK_RETURN(status=1)
  791. ! init array to store 2d field from grib file (north-south order):
  792. call pa_SetShape( pat, lon_n, lat_n )
  793. ! allocate output:
  794. call pa_SetShape( ll, lon_n, lat_n, nlev_out )
  795. ll = 0.0
  796. ! add to history:
  797. write (key,'("longrid==",f7.2,",",f6.2,",",i4)') lon_first, lon_inc, lon_n
  798. call AddHistory( tmi, trim(key), status )
  799. write (key,'("latgrid==",f7.2,",",f6.2,",",i4)') lat_first, lat_inc, lat_n
  800. call AddHistory( tmi, trim(key), status )
  801. ! o gaussian grid
  802. case ( gridtype_gg )
  803. ! routine returns gg grid:
  804. gridtype = 'gg'
  805. ! extract grid number:
  806. call Get( grib, status, N=ggN )
  807. IF_NOTOK_RETURN(status=1)
  808. ! define grid structure:
  809. call Init( ggi, ggN, .true., status )
  810. IF_NOTOK_RETURN(status=1)
  811. ! allocate output:
  812. call pa_SetShape( gg, ggi%np, nlev_out )
  813. gg = 0.0
  814. ! add to history:
  815. write (key,'("gg==",i4.4)') ggN
  816. call AddHistory( tmi, trim(key), status )
  817. ! o spectral field:
  818. case ( gridtype_sh )
  819. ! routine returns sh grid:
  820. gridtype = 'sh'
  821. ! extract spectral truncation:
  822. call Get( grib, status, T=shT )
  823. IF_NOTOK_RETURN(status=1)
  824. ! intialize spherical harmonic field info:
  825. call Init( shi, shT, status )
  826. IF_NOTOK_RETURN(status=1)
  827. ! allocate output:
  828. call pa_SetShape( sh, shi%np, nlev_out )
  829. sh = cmplx(0.0,0.0)
  830. ! add to history:
  831. write (key,'("sh==",i4.4)') shT
  832. call AddHistory( tmi, trim(key), status )
  833. case default
  834. write (gol,'("unsupported gridtype for setup : ",i6)') ggridtype; call goErr
  835. TRACEBACK; status=1; return
  836. end select
  837. !
  838. ! * levels
  839. !
  840. select case ( nlev )
  841. case ( 1 )
  842. call Init( levi, nlev, (/0.0,0.0/), (/0.0,0.0/), status )
  843. IF_NOTOK_RETURN(status=1)
  844. case ( 60 )
  845. call Init( levi, 'ec60', status )
  846. IF_NOTOK_RETURN(status=1)
  847. case ( 91 )
  848. call Init( levi, 'ec91', status )
  849. IF_NOTOK_RETURN(status=1)
  850. case default
  851. write (gol,'("do not how to init levi for nlev = ",i6)') nlev; call goErr
  852. TRACEBACK; status=1; return
  853. end select
  854. ! not again ...
  855. isfirst = .false.
  856. end if ! isfirst (grid definition and allocation)
  857. !
  858. ! store
  859. !
  860. ! layers numbered 1..nlev, half levels numberd 1..nlev+1
  861. ! top-down, thus 1 is space and nlev+1 is surface
  862. if ( nw == 'w' ) then
  863. ! store half levels
  864. select case ( paramkey )
  865. ! store in upper half level of a layer:
  866. case ( 'UDMF', 'DDMF' )
  867. level_out = level
  868. !! store in lower half level of a layer:
  869. !case ( 'dummy' )
  870. ! level_out = level + 1
  871. ! to be implemented ...
  872. case default
  873. write (gol,'("do not if data is on upper or lower half level ...")'); call goErr
  874. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  875. TRACEBACK; status=1; return
  876. end select
  877. else
  878. ! store full levels
  879. level_out = level
  880. end if
  881. select case ( ggridtype )
  882. case ( gridtype_ll )
  883. ! read 2d pat from grib; storred from north to south
  884. call Get( grib, status, ll=pat )
  885. IF_NOTOK_RETURN(status=1)
  886. ! store from south to north:
  887. do ilat = 1, lat_n
  888. ll(:,ilat,level_out) = pat(:,lat_n+1-ilat)
  889. end do
  890. case ( gridtype_gg )
  891. ! read 2d pat from grib:
  892. call Get( grib, status, gg=gg(:,level_out) )
  893. IF_NOTOK_RETURN(status=1)
  894. ! convert from lnsp to sp ?
  895. if ( paramkey == 'sp' .and. spm_lnsp2sp ) gg(:,level_out) = exp(gg(:,level_out)) ! Pa
  896. case ( gridtype_sh )
  897. ! read 2d pat from grib:
  898. call Get( grib, status, sh=sh(:,level_out) )
  899. IF_NOTOK_RETURN(status=1)
  900. case default
  901. write (gol,'("unsupported gridtype for 2d pat : ",i6)') gridtype; call goErr
  902. TRACEBACK; status=1; return
  903. end select
  904. ! last record for this field ?
  905. if ( glevel == nlev ) exit
  906. end do ! records
  907. ! close grib file
  908. call Done( grib, status )
  909. IF_NOTOK_RETURN(status=1)
  910. !
  911. ! ~~~ surface pressure
  912. !
  913. if ( do_spm ) then
  914. ! read lnsp and covert to sp ?
  915. if ( spm_lnsp .or. spm_lnsp2sp) then
  916. spm_levs = levs
  917. spm_paramkey = 'LNSP'
  918. else
  919. spm_levs = 'sfc'
  920. spm_paramkey = 'SP'
  921. end if
  922. ! write gribcode:
  923. call GetPid( 'ec', spm_paramkey, pid, tabid, status )
  924. IF_NOTOK_RETURN(status=1)
  925. write (gribcode,'(i3,".128")') pid
  926. gribcode = adjustl(gribcode)
  927. ! create file name:
  928. ! dir/od-fc-2000-01-ml60-T159-T_20000101_fg006up4tr3.gb
  929. write (spm_fname,'(a,"/",a,"-",a,"-",i4.4,2i2.2,"-",2i2.2,"-",a,"-",a,"-",a,".gb")') &
  930. trim(mf%dir), &
  931. trim(ec_class), trim(ec_type), ccyy, mm, dd, hh, mn, &
  932. trim(spm_levs), trim(gribcode), trim(ec_grid)
  933. ! exist ?
  934. inquire( file=spm_fname, exist=exist )
  935. if ( .not. exist ) then
  936. write (gol,'("grib file does not exist:")'); call goErr
  937. write (gol,'(" ",a)') trim(spm_fname); call goErr
  938. TRACEBACK; status=1; return
  939. end if
  940. ! open grib file
  941. call Init( spm_grib, spm_fname, 'r', status )
  942. IF_NOTOK_RETURN(status=1)
  943. ! loop over time records
  944. do
  945. ! read gribsection in file buffer
  946. call ReadRecord( spm_grib, status )
  947. IF_NOTOK_RETURN(status=1)
  948. ! fill times
  949. call SetTime( mf, spm_grib, status )
  950. IF_NOTOK_RETURN(status=1)
  951. ! check time:
  952. call CheckTime( mf, tref, t1, t2, status )
  953. if (status/=0) then
  954. !write (gol,'("grib read record: spm wrong time; skip ...")'); call goPr
  955. cycle
  956. !write (gol,'("found unexpected times in grib file:")'); call goErr
  957. !write (gol,'(" paramkey : ",a)') paramkey; call goErr
  958. !call wrtgol( ' req. t1 : ', t1 ); call goErr
  959. !call wrtgol( ' t2 : ', t2 ); call goErr
  960. !call wrtgol( ' grib t1 : ', mf%t1 ); call goErr
  961. !call wrtgol( ' t2 : ', mf%t2 ); call goErr
  962. !write (gol,'(" grib file : ",a)') trim(spm_fname); call goErr
  963. !TRACEBACK; status=1; return
  964. end if
  965. ! time ok
  966. exit
  967. end do ! time loop
  968. ! set param id:
  969. select case ( ggridtype )
  970. case ( gridtype_ll )
  971. call GetPid( 'ec', 'SP', pid, tabid, status )
  972. IF_NOTOK_RETURN(status=1)
  973. case ( gridtype_gg )
  974. if ( spm_lnsp2sp ) then
  975. call GetPid( 'ec', 'LNSP', pid, tabid, status )
  976. IF_NOTOK_RETURN(status=1)
  977. else
  978. call GetPid( 'ec', 'SP', pid, tabid, status )
  979. IF_NOTOK_RETURN(status=1)
  980. end if
  981. case ( gridtype_sh )
  982. call GetPid( 'ec', 'LNSP', pid, tabid, status )
  983. IF_NOTOK_RETURN(status=1)
  984. case default
  985. write (gol,'("unsupported gridtype for setup sp/lnsp : ",i6)') ggridtype; call goErr
  986. TRACEBACK; status=1; return
  987. end select
  988. ! check parameter:
  989. call Check( spm_grib, pid=pid, debug=1, status=status )
  990. IF_NOTOK_RETURN(status=1)
  991. ! check level:
  992. call Get( spm_grib, levtype=glevtype, level=glevel, status=status )
  993. IF_NOTOK_RETURN(status=1)
  994. select case ( ggridtype )
  995. case ( gridtype_ll )
  996. if ( glevtype /= levtype_sfc ) then
  997. write (gol,'("found unexpected level type ")'); call goErr
  998. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  999. write (gol,'(" sfc level type : ",i6)') levtype_sfc; call goErr
  1000. write (gol,'(" grib level type : ",i6)') glevtype; call goErr
  1001. write (gol,'(" grib file : ",a)') trim(spm_fname); call goErr
  1002. TRACEBACK; status=1; return
  1003. end if
  1004. case ( gridtype_gg )
  1005. if ( spm_lnsp2sp ) then
  1006. if ( (glevtype /= levtype_hyb) .or. (glevel /= 1) ) then
  1007. write (gol,'("found unexpected level type (lnsp for 3d gg)")'); call goErr
  1008. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  1009. write (gol,'(" hyb level type : ",i6)') levtype_hyb; call goErr
  1010. write (gol,'(" grib level type : ",i6)') glevtype; call goErr
  1011. write (gol,'(" grib level : ",i6)') glevel; call goErr
  1012. write (gol,'(" grib file : ",a)') trim(spm_fname); call goErr
  1013. TRACEBACK; status=1; return
  1014. end if
  1015. else
  1016. if ( glevtype /= levtype_sfc ) then
  1017. write (gol,'("found unexpected level type ")'); call goErr
  1018. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  1019. write (gol,'(" sfc level type : ",i6)') levtype_sfc; call goErr
  1020. write (gol,'(" grib level type : ",i6)') glevtype; call goErr
  1021. write (gol,'(" grib file : ",a)') trim(spm_fname); call goErr
  1022. TRACEBACK; status=1; return
  1023. end if
  1024. end if
  1025. case ( gridtype_sh )
  1026. if ( (glevtype /= levtype_hyb) .or. (glevel /= 1) ) then
  1027. write (gol,'("found unexpected level type ")'); call goErr
  1028. write (gol,'(" paramkey : ",a)') paramkey; call goErr
  1029. write (gol,'(" hyb level type : ",i6)') levtype_hyb; call goErr
  1030. write (gol,'(" grib level type : ",i6)') glevtype; call goErr
  1031. write (gol,'(" grib level : ",i6)') glevel; call goErr
  1032. write (gol,'(" grib file : ",a)') trim(spm_fname); call goErr
  1033. TRACEBACK; status=1; return
  1034. end if
  1035. case default
  1036. write (gol,'("unsupported gridtype for sp/lnsp levs : ",i6)') ggridtype; call goErr
  1037. TRACEBACK; status=1; return
  1038. end select
  1039. ! read and store surface pressure field:
  1040. select case ( ggridtype )
  1041. case ( gridtype_ll )
  1042. ! allocate storage
  1043. call pa_SetShape( sp_ll, lon_n, lat_n )
  1044. ! read 2d pat from grib; storred from north to south
  1045. call Get( spm_grib, status, ll=pat )
  1046. IF_NOTOK_RETURN(status=1)
  1047. ! store from south to north:
  1048. do ilat = 1, lat_n
  1049. sp_ll(:,ilat) = pat(:,lat_n+1-ilat)
  1050. end do
  1051. case ( gridtype_gg )
  1052. ! allocate storage
  1053. call pa_SetShape( sp_gg, ggi%np )
  1054. ! read gg field from grib:
  1055. call Get( spm_grib, status, gg=sp_gg )
  1056. IF_NOTOK_RETURN(status=1)
  1057. ! convert from lnsp to sp ?
  1058. if ( spm_lnsp2sp ) sp_gg = exp(sp_gg) ! Pa
  1059. case ( gridtype_sh )
  1060. ! allocate storage
  1061. call pa_SetShape( lnsp_sh, shi%np )
  1062. ! read spectral coeff:
  1063. call Get( spm_grib, status, sh=lnsp_sh )
  1064. IF_NOTOK_RETURN(status=1)
  1065. case default
  1066. write (gol,'("unsupported gridtype for reading sp/lnsp : ",i6)') ggridtype; call goErr
  1067. TRACEBACK; status=1; return
  1068. end select
  1069. ! close grib file
  1070. call Done( spm_grib, status )
  1071. IF_NOTOK_RETURN(status=1)
  1072. end if
  1073. !
  1074. ! ~~~ end
  1075. !
  1076. ! deallocate arrays
  1077. call pa_Done( pat )
  1078. ! ok
  1079. status = 0
  1080. end subroutine mf_ReadRecord_2
  1081. ! ****************************************************************************
  1082. !
  1083. ! In gribfile:
  1084. ! reftime : for example time at which forecast is made
  1085. ! timerange: increment or interval
  1086. !
  1087. ! arguments: ok if:
  1088. !
  1089. ! time1 == time2 time1==time2 == reftime+timerange
  1090. !
  1091. ! time1 == 0 time2 == reftime+timerange
  1092. !
  1093. ! time2 == 0 time1 == reftime
  1094. !
  1095. ! time1 < time2 time1 == refitme, time2 == reftime+timerange
  1096. !
  1097. !
  1098. ! grib [t1-----------t2]
  1099. ! status
  1100. ! o time1/2 record too old 1
  1101. ! o time1/2 ok 0
  1102. ! o time1/2 record too new 2
  1103. !
  1104. !
  1105. ! SetTime( mf, status )
  1106. !
  1107. ! Extracts time values from current grib record,
  1108. ! store in mf%t1, mf%t2
  1109. !
  1110. ! return status:
  1111. ! 0 : ok
  1112. ! other : some error
  1113. !
  1114. ! CheckTime( mf, time1, time2, status )
  1115. !
  1116. ! return status:
  1117. ! 0 : times match
  1118. ! 1 : times do not match, try next record
  1119. ! 2 : current record is newer than requested (reopen ?)
  1120. ! 3 : some error
  1121. !
  1122. ! ***
  1123. !
  1124. ! Return time parameters for grib files:
  1125. ! o tfile : date in filename
  1126. ! o trange : time interval covered by fields in file
  1127. ! o tref : reference time (forecast start?) for tday,[t1,t2]
  1128. !
  1129. ! Called as:
  1130. ! call GetGribTime( treskey, tday, t1, t2, status, tref=tref )
  1131. ! call GetGribTime( treskey, tref, t1, t2, status, tfile=tfile )
  1132. !
  1133. subroutine GetGribTime( mf, treskey, t0, t1, t2, status, tfile, trange, tref, tshift )
  1134. use GO, only : TDate, TIncrDate, Get, Set, wrtgol, NewDate, IncrDate
  1135. use GO, only : operator(+), operator(-), operator(<), rTotal, iTotal
  1136. use GO, only : AnyDate, IsAnyDate
  1137. use GO, only : wrtgol
  1138. ! --- in/out --------------------------------
  1139. type(TMeteoFile_ecmwf_mars), intent(in) :: mf ! for adhoc flags ...
  1140. character(len=*), intent(in) :: treskey
  1141. type(TDate), intent(in) :: t0, t1, t2
  1142. integer, intent(out) :: status
  1143. type(TDate), intent(out), optional :: tfile
  1144. type(TDate), intent(out), optional :: trange(2)
  1145. type(TDate), intent(out), optional :: tref
  1146. type(TIncrDate), intent(out), optional :: tshift
  1147. ! --- const --------------------------------------
  1148. character(len=*), parameter :: rname = mname//'/GetGribTime'
  1149. ! --- local --------------------------------
  1150. integer :: hour0, hour1, dhour2, time6(6)
  1151. integer :: dd, hh, mn, step
  1152. integer :: anhh
  1153. real :: ddr
  1154. ! --- begin --------------------------------
  1155. !! debug ...
  1156. !write (gol,'(" GetGribTime: treskey : ",a)') trim(treskey); call goPr
  1157. !call wrtgol( ' GetGribTime: t0 : ', t0 ); call goPr
  1158. !call wrtgol( ' GetGribTime: t1 : ', t1 ); call goPr
  1159. !call wrtgol( ' GetGribTime: t2 : ', t2 ); call goPr
  1160. ! files opend upon reading, thus no particular time range for which file is valid:
  1161. if ( present(trange) ) then
  1162. trange(1) = AnyDate()
  1163. trange(2) = AnyDate()
  1164. !! debug ...
  1165. !call wrtgol( ' GetGribTime: set trange to : ', trange(1), ' - ', trange(2) ); call goPr
  1166. end if
  1167. ! zero shift by default
  1168. if ( present(tshift) ) then
  1169. tshift = IncrDate(hour=0)
  1170. !! debug ...
  1171. !call wrtgol( ' GetGribTime: set tshift to : ', tshift ); call goPr
  1172. endif
  1173. ! set day shift, start hour, and step
  1174. select case ( treskey )
  1175. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1176. ! constant field
  1177. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1178. case ( 'const' )
  1179. ! only t0 is usefull ...
  1180. if ( present(tfile ) ) tfile = t0
  1181. !if ( present(trange) ) trange = (/t1,t2/) ! any, any
  1182. if ( present(tref ) ) tref = t0 ! dummy ...
  1183. ! take analysed fields always at least 24 hour old,
  1184. ! since these are the only analysed fields available in forecast mode
  1185. if ( present(tshift) ) then
  1186. tshift = IncrDate(hour=24)
  1187. ! FIX for start of ml91 test suite; try if 12 is ok too ...
  1188. !tshift = IncrDate(hour=12)
  1189. end if
  1190. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1191. ! fc, 3 hourly
  1192. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1193. case ( '_fc012up2tr3' )
  1194. ! hour of t0:
  1195. call Get( t0, hour=hour0 )
  1196. ! end hour counted from t0 00:00
  1197. dhour2 = iTotal( t2 - t0, 'hour' )
  1198. !! debug ...
  1199. !write (gol,'(" GetGribTime: slhf4cnv : ",l2)') mf%slhf_for_convec; call goPr
  1200. !write (gol,'(" GetGribTime: hour0 : ",i3)') hour0; call goPr
  1201. !write (gol,'(" GetGribTime: dhour2 : ",i3)') dhour2; call goPr
  1202. ! set forecast start time and step given hour:
  1203. ! ~~ first the unusual stuf ...
  1204. if ( mf%slhf_for_convec .and. (hour0 == 0) .and. (dhour2 == 3) ) then
  1205. ! (00,03] fc day 0 (slhf for convec) take from (-12,03]
  1206. dd = -1 ; hh = 12 ; step = 15
  1207. else if ( mf%slhf_for_convec .and. (hour0 == 12) .and. (dhour2 == 3) ) then
  1208. ! (12,15] fc day 0 (slhf for convec) take from (00,15]
  1209. dd = 0 ; hh = 00 ; step = 15
  1210. ! ~~ now the usual settings ...
  1211. else if ( dhour2 == 0 ) then
  1212. ! 00:00 fc day 0
  1213. dd = -1 ; hh = 12 ; step = 12
  1214. else if ( (dhour2 <= 12) .and. (modulo(dhour2,3) == 0) ) then
  1215. ! (00,12] fc day 0
  1216. dd = 0
  1217. hh = 00
  1218. step = dhour2
  1219. else if ( ( (t0 < NewDate(year=2006,month=03,day=14) ) .and. &
  1220. ( ((dhour2 <= 12+ 72) .and. (modulo(dhour2,3) == 0)) .or. &
  1221. ((dhour2 <= 12+240) .and. (modulo(dhour2,6) == 0)) ) ) &
  1222. .or. ( ( ((dhour2 <= 12+ 96) .and. (modulo(dhour2,3) == 0)) .or. &
  1223. ((dhour2 <= 12+240) .and. (modulo(dhour2,6) == 0)) ) ) ) then
  1224. ! (12,240] fc days 1-10
  1225. dd = 0
  1226. hh = 12
  1227. step = dhour2 - 12
  1228. else
  1229. write (gol,'("unsupported hour :")'); call goErr
  1230. write (gol,'(" hour0 : ",i3)') hour0; call goErr
  1231. write (gol,'(" dhour2 : ",i3)') dhour2; call goErr
  1232. write (gol,'(" treskey : ",a )') treskey; call goErr
  1233. call wrtgol( ' time1 : ', t1 ); call goErr
  1234. write (gol,'(" timslhf_for_convec : ",l2)') mf%slhf_for_convec; call goErr
  1235. TRACEBACK; status=1; return
  1236. end if
  1237. !! fields valid for hh+(00,12] :
  1238. !if ( present(trange) ) then
  1239. ! trange(1) = t2
  1240. ! call Set( trange(1), hour=hh, min=0, sec=0, mili=0 ) ! hh:00
  1241. ! ! trap 00:00, this should be previous day:
  1242. ! call Get( t2, time6=time6 )
  1243. ! if ( all(time6(4:6)==0) ) trange(1) = trange(1) - IncrDate(day=1)
  1244. ! ! complete (00,12]
  1245. ! trange(2) = trange(1) + IncrDate(hour=12) ! 24:00
  1246. ! call Set( trange(1), mili=1 ) ! > 00:00
  1247. !end if
  1248. ! reference time = start of forecast
  1249. if ( present(tref) ) then
  1250. tref = t0 + IncrDate(day=dd,hour=hh)
  1251. !! debug ...
  1252. !write (gol,*) ' GetGribTime: dd, hh, step : ', dd, hh, step; call goPr
  1253. !call wrtgol( ' GetGribTime: set tref to : ', tref ); call goPr
  1254. end if
  1255. ! adhoc: if tfile is requested, probably the 'tref' returned before
  1256. ! is now in input 't0' ...
  1257. if ( present(tfile) ) then
  1258. tfile = t0 ! .. is tref !
  1259. !! debug ...
  1260. !call wrtgol( ' GetGribTime: set tfile to : ', tfile ); call goPr
  1261. end if
  1262. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1263. ! analysis, files for hours 0, 6, 12, and 18
  1264. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1265. case ( '_an0tr6' )
  1266. ! reference time = analysis time
  1267. if ( present(tref) ) then
  1268. tref = t1
  1269. end if
  1270. ! t0 t1,t2 ->
  1271. ! -+-----------------------+-----------------------+--->
  1272. ! 00 06 12 18 00 06 12 18 00
  1273. ! a
  1274. ! a a
  1275. ! a a(-----------] 00 analysis/forecast
  1276. ! a a(--------------> 12 forecast
  1277. !
  1278. ! -24 -24 -24 -24 -48 shift
  1279. ! take analysed fields always at least 24 hour old,
  1280. ! since these are the only analysed fields available in forecast mode,
  1281. ! and to obtain a contineous time line
  1282. ! t0 is always 00:00
  1283. if ( present(tshift) ) then
  1284. ! difference between t1 and t0 00:00 in fraction of days:
  1285. ddr = rTotal( t1 - t0, 'day' )
  1286. ! set time shift in days:
  1287. tshift = IncrDate(day=floor(ddr)+1)
  1288. !! FIX for start of ml91 test suite; try if 12 is ok too ...
  1289. !tshift = IncrDate(day=floor(ddr),hour=12)
  1290. end if
  1291. ! one file for each time:
  1292. if ( present(tfile) ) then
  1293. tfile = t1
  1294. end if
  1295. !! fields in file valid for instant time:
  1296. !if ( present(trange) ) then
  1297. ! trange(1) = t1
  1298. ! trange(2) = t1
  1299. !end if
  1300. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1301. ! fc, 1 hourly emissions
  1302. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1303. case ( '_fc000up24tr1' )
  1304. ! macc emisisons
  1305. ! fields in grib file have time stamps 00:30, 01:30, etc
  1306. ! reference time = start of forecast
  1307. if ( present(tref) ) then
  1308. ! end hour counted from t0 00:00
  1309. dhour2 = iTotal( t2 - t0, 'hour' )
  1310. ! field from current day, 'forecast time' 00:30 etc:
  1311. tref = t0 + IncrDate(day=0,hour=dhour2-1,min=30)
  1312. ! no steps in this timing convention:
  1313. step = 0
  1314. !! debug ...
  1315. !write (gol,*) ' GetGribTime: dd, hh, step : ', dd, hh, step; call goPr
  1316. !call wrtgol( ' GetGribTime: set tref to : ', tref ); call goPr
  1317. end if
  1318. ! adhoc: if tfile is requested, probably the 'tref' returned before
  1319. ! is now in input 't0' ...
  1320. if ( present(tfile) ) then
  1321. tfile = t0 ! .. is tref !
  1322. !! debug ...
  1323. !call wrtgol( ' GetGribTime: set tfile to : ', tfile ); call goPr
  1324. end if
  1325. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1326. ! ???
  1327. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
  1328. case default
  1329. write (gol,'("unsupported time resolution key:")'); call goErr
  1330. write (gol,'(" `",a,"`")') trim(treskey); call goErr
  1331. TRACEBACK; status=1; return
  1332. end select
  1333. ! ok
  1334. status = 0
  1335. end subroutine GetGribTime
  1336. ! ***
  1337. !
  1338. ! Extract time fields of current grib record,
  1339. ! store in mf%tref, mf%t1, mf%t2
  1340. !
  1341. subroutine SetTime( mf, grib, status )
  1342. use GO, only : TDate, NewDate, IncrDate, operator(+), wrtgol
  1343. use file_grib, only : TGribFile, Check, Get
  1344. ! --- const -------------------------------------
  1345. character(len=*), parameter :: rname = mname//'/SetTime'
  1346. ! --- in/out -------------------------------
  1347. type(TMeteoFile_ecmwf_mars), intent(inout) :: mf
  1348. type(TGribFile), intent(inout) :: grib
  1349. integer, intent(out) :: status
  1350. ! --- local -------------------------------
  1351. integer :: reftime(5), timerange(4)
  1352. ! --- begin -------------------------------
  1353. ! extract time fields from grib record:
  1354. call Get( grib, status, reftime=reftime, timerange=timerange )
  1355. IF_NOTOK_RETURN(status=1)
  1356. ! Fill t1 and t2 with the time information; might be equal.
  1357. ! Check time range indicator (WMO code table 5):
  1358. select case ( timerange(4) )
  1359. case ( 0, 1 )
  1360. !
  1361. ! 0 = Forecast product valid for reference time + P1 (P1>0),
  1362. ! or uninitialized analysis product for reference time (P1=0)
  1363. !
  1364. ! 1 = Initialized analysis product for reference time (P1=0).
  1365. !
  1366. ! fill reference time:
  1367. mf%tref = NewDate( time5=reftime )
  1368. ! fill t1 with reftime+timerange;
  1369. ! add P1 in hours; check time unit (WMO code table 4)
  1370. mf%t1 = NewDate( time5=reftime )
  1371. select case ( timerange(1) )
  1372. case ( 1 ) ! hours
  1373. mf%t1 = mf%t1 + IncrDate( hour=timerange(2) )
  1374. case default
  1375. write (gol,'("grib timerange units other than hours not supported yet")'); call goErr
  1376. write (gol,'(" reftime : ",i4,4i3)') reftime; call goErr
  1377. write (gol,'(" timerange : ",4i3)') timerange; call goErr
  1378. write (gol,'(" file : ",a)') trim(grib%fname); call goErr
  1379. TRACEBACK; status=1; return
  1380. end select
  1381. ! instant time:
  1382. mf%t2 = mf%t1
  1383. case ( 2 )
  1384. !
  1385. ! 2 = Product with a valid time ranging between
  1386. ! reference time + P1 and reference time + P2
  1387. !
  1388. ! fill reftime:
  1389. mf%tref = NewDate( time5=reftime )
  1390. ! fill t1 with reftime+P1;
  1391. ! add P1 in hours; check time unit (WMO code table 4)
  1392. mf%t1 = mf%tref
  1393. select case ( timerange(1) )
  1394. case ( 1 ) ! hours
  1395. mf%t1 = mf%t1 + IncrDate( hour=timerange(2) )
  1396. case default
  1397. write (gol,'("grib timerange units other than hours not supported yet")'); call goErr
  1398. write (gol,'(" file : ",a)') trim(grib%fname); call goErr
  1399. TRACEBACK; status=1; return
  1400. end select
  1401. ! fill t2 with reftime+P2;
  1402. ! add P2 in hours; check time unit (WMO code table 4)
  1403. mf%t2 = mf%tref
  1404. select case ( timerange(1) )
  1405. case ( 1 ) ! hours
  1406. mf%t2 = mf%t2 + IncrDate( hour=timerange(3) )
  1407. case default
  1408. write (gol,'("grib timerange units other than hours not supported yet")'); call goErr
  1409. write (gol,'(" file : ",a)') trim(grib%fname); call goErr
  1410. TRACEBACK; status=1; return
  1411. end select
  1412. case default
  1413. write (gol,'("unsupported time range indicator:")'); call goErr
  1414. write (gol,'(" indicator : ",i6)') timerange(4); call goErr
  1415. write (gol,'(" file : ",a)') trim(grib%fname); call goErr
  1416. TRACEBACK; status=1; return
  1417. end select
  1418. ! ok
  1419. status = 0
  1420. end subroutine SetTime
  1421. ! ***
  1422. subroutine CheckTime( mf, tref, t1, t2, status )
  1423. use GO
  1424. use file_grib, only : Check, Get
  1425. ! --- const -------------------------------------
  1426. character(len=*), parameter :: rname = mname//'/CheckTime'
  1427. ! --- in/out -------------------------------
  1428. type(TMeteoFile_ecmwf_mars), intent(in) :: mf
  1429. type(TDate), intent(in) :: tref, t1, t2
  1430. integer, intent(out) :: status
  1431. ! --- local -------------------------------
  1432. integer :: year1, year2
  1433. ! --- begin -------------------------------
  1434. !call wrtgol( 'CheckTime: (', tref, ') ', t1, ' - ', t2 ); call goPr
  1435. !call wrtgol( 'CheckTime: (', mf%tref, ') ', mf%t1, ' - ', mf%t2 ); call goPr
  1436. ! requested year zero ? always ok
  1437. call Get( t1, year=year1 )
  1438. call Get( t2, year=year2 )
  1439. if ( (year1 == 0) .and. (year2 == 0) ) then
  1440. ! requested constant field, always ok
  1441. status = 0; return ! ok
  1442. end if
  1443. if ( year1 == 0 ) then
  1444. ! do not test t1, only t2
  1445. if ( t2 == mf%t2 ) then
  1446. status = 0; return ! ok
  1447. else
  1448. status = 1; return ! not ok, try next
  1449. end if
  1450. else if ( year2 == 0 ) then
  1451. ! do not test t2, only t1
  1452. if ( t1 == mf%t1 ) then
  1453. status = 0; return ! ok
  1454. else
  1455. status = 1; return ! not ok, try next
  1456. end if
  1457. end if
  1458. ! ! interval or instant time
  1459. ! if ( t1 < t2 ) then
  1460. !
  1461. ! !! time interval: [t1,t2] should be inside [t1,t2]
  1462. ! !if ( (t1 >= mf%t1) .and. (t2 <= mf%t2) ) then
  1463. ! ! status = 0; return ! ok
  1464. ! !else
  1465. ! ! status = 1; return ! try next
  1466. ! !end if
  1467. !
  1468. ! ! time interval: [t1,t2] should be equal to [t1,t2]
  1469. ! if ( (t1 == mf%t1) .and. (t2 == mf%t2) ) then
  1470. ! status = 0; return ! ok
  1471. ! else
  1472. ! status = 1; return ! try next
  1473. ! end if
  1474. !
  1475. ! else if ( t1 == t2 ) then
  1476. !
  1477. ! ! instant time: t2 should match t2
  1478. ! if ( t2 == mf%t2 ) then
  1479. ! status = 0; return ! ok
  1480. ! else
  1481. ! status = 1; return ! try next
  1482. ! end if
  1483. !
  1484. ! else
  1485. !
  1486. ! write (gol,'("t1 should not exceed t2:")'); call goErr
  1487. ! call wrtgol( ' t1 : ', t1 ); call goErr
  1488. ! call wrtgol( ' t2 : ', t2 ); call goErr
  1489. ! TRACEBACK; status=3; return
  1490. !
  1491. ! end if
  1492. ! compare all:
  1493. if ( (tref == mf%tref) .and. (t1 == mf%t1) .and. (t2 == mf%t2) ) then
  1494. status = 0; return ! ok
  1495. else
  1496. status = 1; return ! try next
  1497. end if
  1498. ! some error ...
  1499. status = 1
  1500. end subroutine CheckTime
  1501. end module tmm_mf_ecmwf_mars