meteodata.F90 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826
  1. !### macro's #####################################################
  2. !
  3. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  4. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  5. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  6. !
  7. #include "tm5.inc"
  8. !
  9. !#################################################################
  10. module MeteoData
  11. use GO , only : gol, goErr, goPr
  12. use go , only : TDate
  13. use grid , only : TllGridInfo, TLevelInfo
  14. use tmm , only : TMeteoInfo
  15. use dims , only : nregions_all
  16. implicit none
  17. ! --- const ----------------------------------
  18. ! module name
  19. character(len=*), parameter, private :: mname = 'MeteoData'
  20. ! number of surface types in ECMWF
  21. integer, parameter :: nveg = 20
  22. ! --- types -----------------------------------
  23. ! storage for single meteo field:
  24. type TMeteoData
  25. ! in use ?
  26. logical :: used
  27. ! changed ?
  28. logical :: changed
  29. ! description:
  30. character(len=16) :: name ! field name
  31. character(len=16) :: unit ! kg, K, ...
  32. ! time interpolation:
  33. character(len=10) :: tinterp ! const6, interp3, ...
  34. ! shapes:
  35. integer :: is(2), js(2), ls(2) ! shape w/o halo
  36. integer :: halo
  37. ! target data:
  38. real, pointer :: data(:,:,:)
  39. type(TDate) :: tr(2) ! timerange
  40. type(TMeteoInfo) :: tmi ! history info
  41. ! primary data:
  42. real, pointer :: data1(:,:,:)
  43. logical :: filled1
  44. type(TDate) :: tr1(2) ! timerange
  45. type(TMeteoInfo) :: tmi1 ! history info
  46. ! secondary data:
  47. real, pointer :: data2(:,:,:)
  48. logical :: filled2
  49. type(TDate) :: tr2(2) ! timerange
  50. type(TMeteoInfo) :: tmi2 ! history info
  51. ! input:
  52. character(len=256) :: sourcekey
  53. ! output
  54. logical :: putout
  55. character(len=256) :: destkey
  56. end type TMeteoData
  57. ! --- interfaces -----------------------------------
  58. interface Init
  59. module procedure mdat_Init
  60. end interface
  61. interface Done
  62. module procedure mdat_Done
  63. end interface
  64. interface Set
  65. module procedure mdat_Set
  66. end interface
  67. interface Check
  68. module procedure mdat_Check
  69. end interface
  70. interface Alloc
  71. module procedure mdat_Alloc
  72. end interface
  73. interface TimeInterpolation
  74. module procedure mdat_TimeInterpolation
  75. end interface
  76. interface SetData
  77. module procedure mdat_SetData
  78. end interface
  79. ! --- var ---------------------------------------------
  80. ! horizontal grid definitions
  81. type(TllGridInfo), save :: lli(0:nregions_all) ! local (MPI tile) Lat-Lon grid info
  82. type(TllGridInfo), save :: global_lli(0:nregions_all) ! global Lat-Lon grid info
  83. ! zonal grid definitions:
  84. type(TllGridInfo), save :: lli_z(0:nregions_all) ! global zonal grid info
  85. ! vertical grid definition:
  86. type(TLevelInfo) :: levi
  87. type(TLevelInfo) :: levi_ec
  88. ! meteo fields:
  89. type(TMeteoData), save, target :: sp1_dat(1:nregions_all)
  90. type(TMeteoData), save, target :: sp2_dat(1:nregions_all)
  91. type(TMeteoData), save, target :: sp_dat(1:nregions_all)
  92. type(TMeteoData), save, target :: spm_dat(1:nregions_all)
  93. type(TMeteoData), save, target :: tsp_dat(1:nregions_all)
  94. type(TMeteoData), save, target :: phlb_dat(1:nregions_all)
  95. type(TMeteoData), save, target :: m_dat(1:nregions_all)
  96. type(TMeteoData), save, target :: mfu_dat(1:nregions_all)
  97. type(TMeteoData), save, target :: mfv_dat(1:nregions_all)
  98. type(TMeteoData), save, target :: mfw_dat(1:nregions_all)
  99. type(TMeteoData), save, target :: pu_dat(1:nregions_all)
  100. type(TMeteoData), save, target :: pv_dat(1:nregions_all)
  101. type(TMeteoData), save, target :: pw_dat(1:nregions_all)
  102. type(TMeteoData), save, target :: temper_dat(1:nregions_all)
  103. type(TMeteoData), save, target :: humid_dat(1:nregions_all)
  104. type(TMeteoData), save, target :: gph_dat(1:nregions_all)
  105. type(TMeteoData), save, target :: omega_dat(1:nregions_all)
  106. type(TMeteoData), save, target :: lwc_dat(1:nregions_all)
  107. type(TMeteoData), save, target :: iwc_dat(1:nregions_all)
  108. type(TMeteoData), save, target :: cc_dat(1:nregions_all)
  109. type(TMeteoData), save, target :: cco_dat(1:nregions_all)
  110. type(TMeteoData), save, target :: ccu_dat(1:nregions_all)
  111. type(TMeteoData), save, target :: entu_dat(1:nregions_all)
  112. type(TMeteoData), save, target :: entd_dat(1:nregions_all)
  113. type(TMeteoData), save, target :: detu_dat(1:nregions_all)
  114. type(TMeteoData), save, target :: detd_dat(1:nregions_all)
  115. type(TMeteoData), save, target :: kzz_dat(1:nregions_all)
  116. type(TMeteoData), save, target :: oro_dat(1:nregions_all)
  117. type(TMeteoData), save, target :: lsmask_dat(1:nregions_all)
  118. type(TMeteoData), save, target :: albedo_dat(1:nregions_all)
  119. type(TMeteoData), save, target :: sr_ecm_dat(1:nregions_all)
  120. type(TMeteoData), save, target :: sr_ols_dat(1:nregions_all)
  121. type(TMeteoData), save, target :: ci_dat(1:nregions_all)
  122. type(TMeteoData), save, target :: sst_dat(1:nregions_all)
  123. type(TMeteoData), save, target :: u10m_dat(1:nregions_all)
  124. type(TMeteoData), save, target :: v10m_dat(1:nregions_all)
  125. type(TMeteoData), save, target :: wspd_dat(1:nregions_all)
  126. type(TMeteoData), save, target :: sshf_dat(1:nregions_all)
  127. type(TMeteoData), save, target :: slhf_dat(1:nregions_all)
  128. type(TMeteoData), save, target :: ewss_dat(1:nregions_all)
  129. type(TMeteoData), save, target :: nsss_dat(1:nregions_all)
  130. type(TMeteoData), save, target :: lsp_dat(1:nregions_all)
  131. type(TMeteoData), save, target :: cp_dat(1:nregions_all)
  132. type(TMeteoData), save, target :: sd_dat(1:nregions_all)
  133. type(TMeteoData), save, target :: swvl1_dat(1:nregions_all)
  134. type(TMeteoData), save, target :: stl1_dat(1:nregions_all)
  135. type(TMeteoData), save, target :: src_dat(1:nregions_all)
  136. type(TMeteoData), save, target :: d2m_dat(1:nregions_all)
  137. type(TMeteoData), save, target :: t2m_dat(1:nregions_all)
  138. type(TMeteoData), save, target :: ssr_dat(1:nregions_all)
  139. type(TMeteoData), save, target :: ssrd_dat(1:nregions_all)
  140. type(TMeteoData), save, target :: str_dat(1:nregions_all)
  141. type(TMeteoData), save, target :: strd_dat(1:nregions_all)
  142. type(TMeteoData), save, target :: skt_dat(1:nregions_all)
  143. type(TMeteoData), save, target :: tv_dat(1:nregions_all,nveg)
  144. type(TMeteoData), save, target :: cvl_dat(1:nregions_all)
  145. type(TMeteoData), save, target :: cvh_dat(1:nregions_all)
  146. type(TMeteoData), save, target :: blh_dat(1:nregions_all)
  147. type(TMeteoData), save, target :: sf_dat(1:nregions_all)
  148. type(TMeteoData), save, target :: g10m_dat(1:nregions_all)
  149. contains
  150. ! ==========================================================
  151. subroutine mdat_Init( md, name, unit, tinterp, is, js, halo, ls, &
  152. sourcekey, putout, destkey, status )
  153. ! --- in/out -----------------------------------
  154. type(TMeteoData), intent(out) :: md
  155. character(len=*), intent(in) :: name, unit
  156. character(len=*), intent(in) :: tinterp
  157. integer, intent(in) :: is(2), js(2)
  158. integer, intent(in) :: halo
  159. integer, intent(in) :: ls(2)
  160. character(len=*), intent(in) :: sourcekey
  161. logical, intent(in) :: putout
  162. character(len=*), intent(in) :: destkey
  163. integer, intent(out) :: status
  164. ! --- const -------------------------------
  165. character(len=*), parameter :: rname = mname//'/mdat_Init'
  166. ! --- begin --------------------------------
  167. ! not in use:
  168. md%used = .false.
  169. ! not changed yet
  170. md%changed = .false.
  171. ! store description:
  172. md%name = name
  173. md%unit = unit
  174. ! store time info:
  175. md%tinterp = tinterp
  176. ! store data shape:
  177. md%is = is
  178. md%js = js
  179. md%halo = halo
  180. md%ls = ls
  181. ! no data allocated yet:
  182. nullify( md%data )
  183. ! no primary data allocated yet:
  184. nullify( md%data1 )
  185. md%filled1 = .false.
  186. ! no secondary data allocated yet:
  187. nullify( md%data2 )
  188. md%filled2 = .false.
  189. ! store input info:
  190. md%sourcekey = sourcekey
  191. ! store output info:
  192. md%putout = putout
  193. md%destkey = destkey
  194. ! ok
  195. status = 0
  196. end subroutine mdat_Init
  197. ! ***
  198. subroutine mdat_Done( md, status )
  199. ! --- in/out -----------------------------------
  200. type(TMeteoData), intent(inout) :: md
  201. integer, intent(out) :: status
  202. ! --- const -------------------------------
  203. character(len=*), parameter :: rname = mname//'/mdat_Done'
  204. ! --- begin --------------------------------
  205. ! deallocate target data if neccesary:
  206. if ( associated(md%data) ) then
  207. ! target data points to data1 ?
  208. if ( associated(md%data,md%data1) ) then
  209. ! data points to data1 :
  210. nullify( md%data )
  211. else
  212. ! data is allocated; clear:
  213. deallocate( md%data )
  214. end if
  215. end if
  216. ! deallocate primary and seondary data:
  217. if ( associated(md%data1) ) deallocate( md%data1 )
  218. if ( associated(md%data2) ) deallocate( md%data2 )
  219. ! for safety ...
  220. md%used = .false.
  221. md%name = 'none'
  222. md%unit = 'none'
  223. ! ok
  224. status = 0
  225. end subroutine mdat_Done
  226. ! ***
  227. subroutine mdat_Set( md, status, used )
  228. ! --- in/out -----------------------------------
  229. type(TMeteoData), intent(inout) :: md
  230. integer, intent(out) :: status
  231. logical, intent(in), optional :: used
  232. ! --- const -------------------------------
  233. character(len=*), parameter :: rname = mname//'/mdat_Set'
  234. ! --- begin --------------------------------
  235. ! code to quickly find which met field is used in your configuration (useful when developping ecearth)
  236. ! if ( present(used) ) then
  237. ! if (used) then
  238. ! write(gol,*)"USED MET FIELD:",md%name; call goPr
  239. ! endif
  240. ! endif
  241. if ( present(used) ) md%used = used
  242. ! ok
  243. status = 0
  244. end subroutine mdat_Set
  245. ! ***
  246. subroutine mdat_Check( md, status )
  247. ! --- in/out -----------------------------------
  248. type(TMeteoData), intent(inout) :: md
  249. integer, intent(out) :: status
  250. ! --- const -------------------------------
  251. character(len=*), parameter :: rname = mname//'/mdat_Check'
  252. ! --- begin --------------------------------
  253. if ( .not. md%used ) then
  254. write (gol,'("meteo data `",a,"` not in use ...")') trim(md%name); call goErr
  255. TRACEBACK; status=1; return
  256. end if
  257. ! ok
  258. status = 0
  259. end subroutine mdat_Check
  260. ! ***
  261. subroutine mdat_Alloc( md, status )
  262. ! --- in/out -----------------------------------
  263. type(TMeteoData), intent(inout) :: md
  264. integer, intent(out) :: status
  265. ! --- const -------------------------------
  266. character(len=*), parameter :: rname = mname//'/mdat_Alloc'
  267. ! --- begin --------------------------------
  268. ! allocate if field is in use:
  269. if ( md%used ) then
  270. ! allocate target, primary, and/or secondary data array;
  271. ! set all to zero to avoid fpe in halo cells during array operations:
  272. select case ( md%tinterp )
  273. ! computed data is stored in target data array:
  274. case ( 'computed' )
  275. ! allocate target data:
  276. allocate( md%data( md%is(1)-md%halo:md%is(2)+md%halo, &
  277. md%js(1)-md%halo:md%js(2)+md%halo, &
  278. md%ls(1) :md%ls(2) ) )
  279. md%data = 0.0
  280. ! constant data is storred in primary data;
  281. ! target data points to data1 :
  282. case ( 'const', 'month', 'const6', 'const3', 'cpl6', 'cpl3', 'cpl2', 'cpl1' )
  283. ! allocate primary data:
  284. allocate( md%data1( md%is(1)-md%halo:md%is(2)+md%halo, &
  285. md%js(1)-md%halo:md%js(2)+md%halo, &
  286. md%ls(1) :md%ls(2) ) )
  287. md%data1 = 0.0
  288. ! point target data:
  289. md%data => md%data1
  290. case ( 'interp6_3', 'interp6', 'interp3', 'interp2', 'interp1', &
  291. 'aver1', 'aver3', 'aver6', 'aver24', 'aver24_3' )
  292. ! allocate target data:
  293. allocate( md%data( md%is(1)-md%halo:md%is(2)+md%halo, &
  294. md%js(1)-md%halo:md%js(2)+md%halo, &
  295. md%ls(1) :md%ls(2) ) )
  296. md%data = 0.0
  297. ! allocate primary data:
  298. allocate( md%data1( md%is(1)-md%halo:md%is(2)+md%halo, &
  299. md%js(1)-md%halo:md%js(2)+md%halo, &
  300. md%ls(1) :md%ls(2) ) )
  301. md%data1 = 0.0
  302. ! allocate secondary data:
  303. allocate( md%data2( md%is(1)-md%halo:md%is(2)+md%halo, &
  304. md%js(1)-md%halo:md%js(2)+md%halo, &
  305. md%ls(1) :md%ls(2) ) )
  306. md%data2 = 0.0
  307. case default
  308. write (gol,'("unsupported time interpolation:")'); call goErr
  309. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  310. write (gol,'(" md%name : ",a)') trim(md%name); call goErr
  311. TRACEBACK; status=1; return
  312. end select
  313. end if
  314. ! ok
  315. status = 0
  316. end subroutine mdat_Alloc
  317. ! ***
  318. subroutine mdat_SetData( md, mdin, status )
  319. ! Copy MDIN to MD : only %data and %tr
  320. use GO, only : TDate
  321. ! --- in/out -----------------------------------
  322. type(TMeteoData), intent(inout) :: md
  323. type(TMeteoData), intent(in) :: mdin
  324. integer, intent(out) :: status
  325. ! --- const -------------------------------
  326. character(len=*), parameter :: rname = mname//'/mdat_SetData'
  327. ! --- begin --------------------------------
  328. ! skip rest ?
  329. if ( .not. md%used ) then
  330. write (gol,'("WARNING - destination meteo data not in use ...")'); call goPr
  331. write (*,'("WARNING in ",a)') rname; status=1; return
  332. status=0; return
  333. end if
  334. ! check shapes
  335. if ( any(md%is /=mdin%is ) .or. &
  336. any(md%js /=mdin%js ) .or. &
  337. md%halo/=mdin%halo .or. &
  338. any(md%ls /=mdin%ls ) ) then
  339. write (gol,'("destination and source shapes should be the same:")'); call goErr
  340. write (gol,'(" is : ",2i4," , ",2i4)') md%is , mdin%is; call goErr
  341. write (gol,'(" js : ",2i4," , ",2i4)') md%js , mdin%js; call goErr
  342. write (gol,'(" halo : ", i8," , ", i8)') md%halo, mdin%halo; call goErr
  343. write (gol,'(" ls : ",2i4," , ",2i4)') md%ls , mdin%ls; call goErr
  344. TRACEBACK; status=1; return
  345. end if
  346. ! check source data:
  347. if ( .not. associated(mdin%data) ) then
  348. write (gol,'("source data not allocated ...")'); call goErr
  349. TRACEBACK; status=1; return
  350. end if
  351. !if ( .not. mdin%filled ) then
  352. ! write (gol,'("source data not filled")'); call goErr
  353. ! TRACEBACK; status=1; return
  354. !end if
  355. ! check target data:
  356. if ( md%tinterp /= 'computed' ) then
  357. write (gol,'("destination data has wrong tinterp:")'); call goErr
  358. write (gol,'(" expected : ",a)') 'computed'; call goErr
  359. write (gol,'(" found : ",a)') trim(md%tinterp); call goErr
  360. TRACEBACK; status=1; return
  361. end if
  362. if ( .not. associated(md%data) ) then
  363. write (gol,'("destination data not allocated ...")'); call goErr
  364. TRACEBACK; status=1; return
  365. end if
  366. ! check shapes
  367. if ( any( shape(md%data) /= shape(mdin%data) ) ) then
  368. write (gol,'("shapes are not the same:")'); call goErr
  369. write (gol,'(" md : ",3i5)') shape( md%data); call goErr
  370. write (gol,'(" mdin : ",3i5)') shape(mdin%data); call goErr
  371. TRACEBACK; status=1; return
  372. end if
  373. ! copy data:
  374. md%data = mdin%data
  375. md%tr = mdin%tr
  376. ! ok
  377. status = 0
  378. end subroutine mdat_SetData
  379. ! ***
  380. subroutine mdat_TimeInterpolation( md, tr, status )
  381. use go , only : TDate, NewDate, IncrDate, Get
  382. use go , only : wrtgol, InterpolFractions, rTotal
  383. use go , only : operator(/=), operator(<), operator(<=)
  384. use go , only : operator(+), operator(-), operator(/)
  385. use tmm, only : SetHistory, AddHistory
  386. ! --- in/out ----------------------------------
  387. type(TMeteoData), intent(inout) :: md
  388. type(TDate), intent(in) :: tr(2)
  389. integer, intent(out) :: status
  390. ! --- const ----------------------------------
  391. character(len=*), parameter :: rname = mname//'/mdat_TimeInterpolation'
  392. ! --- local ----------------------------------
  393. integer :: dth, baseh
  394. integer :: year, month, day, hour, minu
  395. type(TDate) :: tmid, tc(2)
  396. real :: alfa1, alfa2
  397. ! --- begin -----------------------------
  398. ! not used ? error
  399. if ( .not. md%used ) then
  400. write (gol,'("meteo data `",a,"` not used")') trim(md%name); call goErr
  401. TRACEBACK; status=1; return
  402. end if
  403. ! different actions based on time interpolation type:
  404. select case ( md%tinterp )
  405. !
  406. ! constant data
  407. !
  408. case ( 'const' )
  409. ! md%data points to md%data1, so nothing to be done
  410. !
  411. ! data constant during interval
  412. !
  413. case ( 'month' )
  414. ! check time: tr should be in md%tr1
  415. if ( (tr(2) < md%tr1(1)) .or. (md%tr1(2) < tr(1)) ) then
  416. write (gol,'("model data does not include requested interval:")'); call goErr
  417. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  418. call wrtgol( ' md%tr1 : ', md%tr1(1), ' - ', md%tr1(2) ); call goErr
  419. call wrtgol( ' tr : ', tr(1), ' - ', tr(2) ); call goErr
  420. TRACEBACK; status=1; return
  421. end if
  422. ! md%data points to md%data1, so no changes
  423. case ( 'const6', 'const3' )
  424. select case ( md%tinterp )
  425. case ( 'const3' ) ; baseh = 0 ; dth = 3
  426. case ( 'const6' ) ; baseh = 0 ; dth = 6
  427. end select
  428. ! extract time values for begin of current interval:
  429. call Get( tr(1), year, month, day, hour, minu )
  430. ! round hour to 00/06/12/18 or 00/03/06/09/12/15/18/21 or 09
  431. !WP! changed line to NOT use the nint function, instead just use int+0.5 which is the same
  432. !WP! for positive numbers, see http://h21007.www2.hp.com/portal/download/files/unprot/Fortran/docs/lrm/lrm0299.htm
  433. !WP! For negative numbers, the lines below need to be added, this can happen only if baseh!=0
  434. !WP! This fixes a problem on the NOAA machines that cause unpredictable crashes when using nint()
  435. !WP!
  436. !dummy=real(hour+minu/60.0-baseh)/real(dth)
  437. !if(dummy>0.0) hour=dth*int(dummy+0.5)+baseh
  438. !if(dummy<0.0) hour=dth*int(dummy-0.5)+baseh
  439. !endif
  440. !WP!
  441. !WP!
  442. ! ported from cy2, 24 Jan 2008, ARJ
  443. hour = dth * int(real(hour+minu/60.0-baseh)/real(dth)+0.5) + baseh
  444. ! set mid of 3 or 6 hour interval:
  445. tmid = NewDate( year, month, day, hour )
  446. ! interval with constant field
  447. tc(1) = tmid - IncrDate(hour=dth)/2
  448. tc(2) = tmid + IncrDate(hour=dth)/2
  449. ! check interval
  450. if ( (tr(1) < tc(1)) .or. (tc(2) < tr(2)) ) then
  451. write (gol,'("time intervals do not match:")'); call goErr
  452. call wrtgol( ' requested : ', tr(1), ' - ', tr(2) ); call goErr
  453. call wrtgol( ' mdin valid for : ', tc(1), ' - ', tc(2) ) ; call goErr
  454. write (gol,'(" mdin%tinterp : ",a)') trim(md%tinterp); call goErr
  455. TRACEBACK; status=1; return
  456. end if
  457. ! md%data points to md%data1, so no changes
  458. !
  459. ! coupling: get field for t, use it for [t,t+dt]
  460. !
  461. case ( 'cpl6', 'cpl3', 'cpl2', 'cpl1' )
  462. select case ( md%tinterp )
  463. case ( 'cpl6' ) ; baseh = 0 ; dth = 6
  464. case ( 'cpl3' ) ; baseh = 0 ; dth = 3
  465. case ( 'cpl2' ) ; baseh = 0 ; dth = 2
  466. case ( 'cpl1' ) ; baseh = 0 ; dth = 1
  467. end select
  468. ! extract time values for begin of current interval:
  469. call Get( tr(1), year, month, day, hour, minu )
  470. ! round hour to 00/06/12/18 or 00/03/06/09/12/15/18/21
  471. hour = dth * floor(real(hour+minu/60.0-baseh)/real(dth)) + baseh
  472. ! interval with constant field
  473. tc(1) = NewDate( year, month, day, hour )
  474. tc(2) = tc(1) + IncrDate(hour=dth)
  475. ! check interval
  476. if ( (tr(1) < tc(1)) .or. (tc(2) < tr(2)) ) then
  477. write (gol,'("time intervals do not match:")'); call goErr
  478. call wrtgol( ' requested : ', tr(1), ' - ', tr(2) ); call goErr
  479. call wrtgol( ' mdin valid for : ', tc(1), ' - ', tc(2) ) ; call goErr
  480. write (gol,'(" mdin%tinterp : ",a)') trim(md%tinterp); call goErr
  481. TRACEBACK; status=1; return
  482. end if
  483. ! md%data points to md%data1, so no changes
  484. !
  485. ! linear interpolation between instant times
  486. !
  487. case ( 'interp6', 'interp6_3', 'interp3', 'interp2', 'interp1' )
  488. ! not filled ? error
  489. if ( (.not. md%filled1) .or. (.not. md%filled2) ) then
  490. write (gol,'("meteo data not filled:")'); call goErr
  491. write (gol,'(" name : ",a)') trim(md%name); call goErr
  492. write (gol,'(" filled : ",2l2)') md%filled1, md%filled2; call goErr
  493. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  494. TRACEBACK; status=1; return
  495. end if
  496. ! interpolation between instant times, not between intervals ...
  497. if ( (md%tr1(1) /= md%tr1(2)) .or. (md%tr2(1) /= md%tr2(2)) ) then
  498. write (gol,'("time interpolation not between intervals:")'); call goErr
  499. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  500. call wrtgol( ' tr1 : ', md%tr1(1), ' - ', md%tr1(2) ); call goErr
  501. call wrtgol( ' tr2 : ', md%tr2(1), ' - ', md%tr2(2) ); call goErr
  502. TRACEBACK; status=1; return
  503. end if
  504. ! interpolate to mid of interval:
  505. tmid = tr(1) + (tr(2)-tr(1))/2
  506. ! deterimine weights to data and data2 :
  507. call InterpolFractions( tmid, md%tr1(1), md%tr2(1), alfa1, alfa2, status )
  508. if (status/=0) then; TRACEBACK; return; end if
  509. !$OMP PARALLEL &
  510. !$OMP default (none) &
  511. !$OMP shared (alfa1, alfa2, md)
  512. md%data = alfa1 * md%data1 + alfa2 * md%data2
  513. !$OMP END PARALLEL
  514. ! data is changed ...
  515. md%changed = .true.
  516. !
  517. ! fractions of time average fields:
  518. ! data1 : [tr1(1),tr1(2)]
  519. ! data2 : [tr1(1),tr1(2)]
  520. ! tr : [tr(1),tr(2)]
  521. !
  522. case ( 'aver1', 'aver3', 'aver6', 'aver24', 'aver24_3' )
  523. ! primary data not filled ? error
  524. if ( .not. md%filled1 ) then
  525. write (gol,'("meteo data1 not filled:")'); call goErr
  526. write (gol,'(" name : ",a)') trim(md%name); call goErr
  527. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  528. TRACEBACK; status=1; return
  529. end if
  530. ! tr earlier than tr1 ? error ...
  531. if ( tr(1) < md%tr1(1) ) then
  532. write (gol,'("requested time interval earlier than data:")'); call goErr
  533. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  534. call wrtgol( ' md%tr1 : ', md%tr1(1), ' - ', md%tr1(2) ); call goErr
  535. call wrtgol( ' tr : ', tr(1), ' - ', tr(2) ); call goErr
  536. TRACEBACK; status=1; return
  537. end if
  538. ! tr complete in tr1 ? simple ...
  539. if ( tr(2) <= md%tr1(2) ) then
  540. ! just copy ...
  541. md%data = md%data1
  542. ! data is changed ...
  543. md%changed = .true.
  544. else
  545. ! fractions of data1 and data2
  546. ! secondary data not filled ? error
  547. if ( .not. md%filled2 ) then
  548. write (gol,'("meteo data2 not filled:")'); call goErr
  549. write (gol,'(" name : ",a)') trim(md%name); call goErr
  550. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  551. TRACEBACK; status=1; return
  552. end if
  553. ! time ranges for data1 and data2 should be connected:
  554. if ( md%tr1(2) /= md%tr2(1) ) then
  555. write (gol,'("time intervals not connected:")'); call goErr
  556. call wrtgol( ' md%tr1 : ', md%tr1(1), ' - ', md%tr1(2) ); call goErr
  557. call wrtgol( ' md%tr2 : ', md%tr2(1), ' - ', md%tr2(2) ); call goErr
  558. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  559. TRACEBACK; status=1; return
  560. end if
  561. ! check requested time range:
  562. if ( (tr(1) < md%tr1(1)) .or. (md%tr2(2) < tr(2)) ) then
  563. write (gol,'("requested time interval not covered by data :")'); call goErr
  564. call wrtgol( ' md%tr1 : ', md%tr1(1), ' - ', md%tr1(2) ); call goErr
  565. call wrtgol( ' md%tr2 : ', md%tr2(1), ' - ', md%tr2(2) ); call goErr
  566. call wrtgol( ' tr : ', tr(1), ' - ', tr(2) ); call goErr
  567. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  568. TRACEBACK; status=1; return
  569. end if
  570. ! first fraction:
  571. if ( tr(1) < md%tr1(2) ) then
  572. if ( tr(2) < md%tr1(2) ) then
  573. ! tr complete in tr1 ...
  574. alfa1 = 1.0
  575. else
  576. ! fraction of tr inside tr1 :
  577. alfa1 = rTotal(md%tr1(2)-tr(1),'sec') / rTotal(tr(2)-tr(1),'sec')
  578. end if
  579. else
  580. ! tr later than tr1, thus not covered:
  581. alfa1 = 0.0
  582. end if
  583. ! second fraction:
  584. if ( md%tr2(1) < tr(2) ) then
  585. if ( md%tr2(1) < tr(1) ) then
  586. ! tr complete in tr2 ...
  587. alfa2 = 1.0
  588. else
  589. ! fraction of tr inside tr2 :
  590. alfa2 = rTotal(tr(2)-md%tr2(1),'sec') / rTotal(tr(2)-tr(1),'sec')
  591. end if
  592. else
  593. ! tr before tr2, thus not covered:
  594. alfa2 = 0.0
  595. end if
  596. ! replace data array:
  597. md%data = alfa1 * md%data1 + alfa2 * md%data2
  598. ! data is changed ...
  599. md%changed = .true.
  600. end if
  601. !
  602. ! unknown ...
  603. !
  604. case default
  605. write (gol,'("unsupported time interpolation:")'); call goErr
  606. write (gol,'(" md%tinterp : ",a)') trim(md%tinterp); call goErr
  607. write (gol,'(" md%name : ",a)') trim(md%name); call goErr
  608. TRACEBACK; status=1; return
  609. end select
  610. ! store new time:
  611. md%tr = tr
  612. ! copy history:
  613. call SetHistory( md%tmi, md%tmi1, status )
  614. ! ok
  615. status = 0
  616. end subroutine mdat_TimeInterpolation
  617. end module MeteoData