tm5_restart.F90 44 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141114211431144114511461147114811491150115111521153115411551156115711581159116011611162116311641165116611671168116911701171117211731174117511761177117811791180118111821183118411851186118711881189119011911192119311941195119611971198119912001201120212031204120512061207120812091210121112121213121412151216121712181219122012211222122312241225122612271228122912301231123212331234123512361237123812391240
  1. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  2. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  3. #define IF_NOTOK_MDF(action) if (status/=0) then; TRACEBACK; action; if (isRoot) call MDF_CLose(fid,status); status=1; return; end if
  4. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  5. #include "tm5.inc"
  6. !----------------------------------------------------------------------------
  7. ! TM5 !
  8. !----------------------------------------------------------------------------
  9. !BOP
  10. !
  11. ! !MODULE: RESTART
  12. !
  13. ! !DESCRIPTION: Write and read restart files.
  14. !\\
  15. !\\
  16. ! !INTERFACE:
  17. !
  18. MODULE RESTART
  19. !
  20. ! !USES:
  21. !
  22. use GO , only : gol, goPr, goErr
  23. use dims , only : nregions
  24. implicit none
  25. private
  26. !
  27. ! !PUBLIC MEMBER FUNCTIONS:
  28. !
  29. public :: Restart_Init ! read restart keys in rc file
  30. public :: Restart_Done ! nothing yet
  31. public :: Restart_Save ! wrapper around Restart_Write
  32. public :: Restart_Write ! write a restart file
  33. public :: Restart_Read ! read a restart file
  34. public :: rs_write ! model must write restart
  35. !
  36. ! !PRIVATE DATA MEMBERS:
  37. !
  38. character(len=*), parameter :: mname = 'Restart'
  39. character(len=256) :: rs_write_dir
  40. logical :: rs_write
  41. logical :: rs_write_extra
  42. integer :: rs_write_extra_dhour, rs_write_extra_hour
  43. integer :: fid ! file id for IF_NOTOK_MDF macro
  44. !
  45. ! !REVISION HISTORY:
  46. ! 8 Apr 2011 - P. Le Sager - Close MDF file if error occurs. This is
  47. ! needed for mpi_abort not to hang. See TM5_MPI_Abort in
  48. ! partools, and remarks below. Made IF_NOTOK_MDF macro for
  49. ! that purpose.
  50. ! 28 Apr 2011 - P. Le Sager - Read method : handle restart file with extra
  51. ! tracers.
  52. ! 10 Oct 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  53. !
  54. ! !REMARKS:
  55. ! (1) when an error occurs when accessing MDF files, you should first close
  56. ! the file before returning. The IF_NOTOK_MDF macro takes care of that.
  57. ! The only thing you need is to call it like that :
  58. !
  59. ! IF_NOTOK_MDF(fid=xxxx)
  60. !
  61. ! where you replace xxxx with the integer id (file handler) of the file
  62. ! you are accessing. Note that this does not solve all problems (but
  63. ! probably most of them): it is still possible that MDF_Close hangs...
  64. !
  65. !EOP
  66. !------------------------------------------------------------------------
  67. CONTAINS
  68. !--------------------------------------------------------------------------
  69. ! TM5 !
  70. !--------------------------------------------------------------------------
  71. !BOP
  72. !
  73. ! !IROUTINE: RESTART_INIT
  74. !
  75. ! !DESCRIPTION: read settings from rcfile
  76. !\\
  77. !\\
  78. ! !INTERFACE:
  79. !
  80. SUBROUTINE RESTART_INIT( status )
  81. !
  82. ! !USES:
  83. !
  84. use GO , only : TrcFile, Init, Done, ReadRc
  85. use global_data, only : rcfile
  86. use global_data, only : outdir
  87. use meteodata , only : lli
  88. !
  89. ! !OUTPUT PARAMETERS:
  90. !
  91. integer, intent(out) :: status
  92. !
  93. ! !REVISION HISTORY:
  94. !
  95. !EOP
  96. !------------------------------------------------------------------------
  97. !BOC
  98. character(len=*), parameter :: rname = 'Restart_Init'
  99. type(TrcFile) :: rcF
  100. ! ---- begin
  101. call Init( rcF, rcfile, status )
  102. IF_NOTOK_RETURN(status=1)
  103. ! write restart files at all ?
  104. call ReadRc( rcF, 'restart.write', rs_write, status, default=.false. )
  105. IF_ERROR_RETURN(status=1)
  106. ! further settings ...
  107. if ( rs_write ) then
  108. ! output directory:
  109. call ReadRc( rcF, 'restart.write.dir', rs_write_dir, status, default=outdir )
  110. IF_ERROR_RETURN(status=1)
  111. ! extra restart files ?
  112. call ReadRc( rcF, 'restart.write.extra', rs_write_extra, status, default=.false. )
  113. IF_ERROR_RETURN(status=1)
  114. if ( rs_write_extra ) then
  115. call ReadRc( rcF, 'restart.write.extra.hour', rs_write_extra_hour, status, default=0 )
  116. IF_ERROR_RETURN(status=1)
  117. call ReadRc( rcF, 'restart.write.extra.dhour', rs_write_extra_dhour, status, default=24 )
  118. IF_ERROR_RETURN(status=1)
  119. end if
  120. end if ! write restart files
  121. call Done( rcF, status )
  122. IF_NOTOK_RETURN(status=1)
  123. status = 0
  124. END SUBROUTINE RESTART_INIT
  125. !EOC
  126. !--------------------------------------------------------------------------
  127. ! TM5 !
  128. !--------------------------------------------------------------------------
  129. !BOP
  130. !
  131. ! !IROUTINE: RESTART_DONE
  132. !
  133. ! !DESCRIPTION:
  134. !\\
  135. !\\
  136. ! !INTERFACE:
  137. !
  138. SUBROUTINE RESTART_DONE( status )
  139. !
  140. ! !OUTPUT PARAMETERS:
  141. !
  142. integer, intent(out) :: status
  143. !
  144. ! !REVISION HISTORY:
  145. !
  146. !EOP
  147. !------------------------------------------------------------------------
  148. !BOC
  149. character(len=*), parameter :: rname = 'Restart_Done'
  150. ! --- begin --------------------------------
  151. ! nothing to be done ...
  152. ! ok
  153. status = 0
  154. END SUBROUTINE RESTART_DONE
  155. !EOC
  156. !--------------------------------------------------------------------------
  157. ! TM5 !
  158. !--------------------------------------------------------------------------
  159. !BOP
  160. !
  161. ! !IROUTINE: RESTART_SAVE
  162. !
  163. ! !DESCRIPTION:
  164. !\\
  165. !\\
  166. ! !INTERFACE:
  167. !
  168. SUBROUTINE RESTART_SAVE( status, extra, isfirst )
  169. !
  170. ! !USES:
  171. !
  172. use dims, only : idate
  173. !
  174. ! !OUTPUT PARAMETERS:
  175. !
  176. integer, intent(out) :: status
  177. !
  178. ! !INPUT PARAMETERS:
  179. !
  180. logical, intent(in), optional :: extra
  181. logical, intent(in), optional :: isfirst
  182. !
  183. ! !REVISION HISTORY:
  184. !
  185. !EOP
  186. !------------------------------------------------------------------------
  187. !BOC
  188. character(len=*), parameter :: rname = 'Restart_Save'
  189. logical :: is_extra
  190. real :: t1, t2
  191. ! --- begin --------------------------------
  192. ! options ...
  193. is_extra = .false.
  194. if ( present(extra) ) is_extra = extra
  195. ! write restart files at all ?
  196. if ( rs_write ) then
  197. ! end or extra ?
  198. if ( is_extra ) then
  199. ! save extra restart files ?
  200. if ( rs_write_extra ) then
  201. ! every hour+n*dhour only :
  202. if ( modulo( idate(4) - rs_write_extra_hour, rs_write_extra_dhour ) == 0 .and. &
  203. all( idate(5:6) == 0 ) ) then
  204. ! write restart file for this time:
  205. call Restart_Write( status, isfirst=isfirst )
  206. IF_NOTOK_RETURN(status=1)
  207. end if ! for this hour
  208. end if ! extra restart files ?
  209. else
  210. ! write restart file :
  211. call cpu_time(t1)
  212. call Restart_Write( status, isfirst=isfirst )
  213. IF_NOTOK_RETURN(status=1)
  214. call cpu_time(t2)
  215. write (gol,*) " time to write restart [s]: ", t2-t1 ; call goPr
  216. end if ! not extra
  217. end if ! write at all
  218. ! ok
  219. status = 0
  220. END SUBROUTINE RESTART_SAVE
  221. !EOC
  222. !--------------------------------------------------------------------------
  223. ! TM5 !
  224. !--------------------------------------------------------------------------
  225. !BOP
  226. !
  227. ! !IROUTINE: RESTART_FILENAME
  228. !
  229. ! !DESCRIPTION: Build restart filename from inputs.
  230. !\\
  231. !\\
  232. ! !INTERFACE:
  233. !
  234. SUBROUTINE RESTART_FILENAME( region, fname, status, key, dir, isfirst )
  235. !
  236. ! !USES:
  237. !
  238. use dims , only : idate
  239. use global_data, only : outdir
  240. use meteodata , only : lli
  241. !
  242. ! !INPUT PARAMETERS:
  243. !
  244. integer, intent(in) :: region
  245. logical, intent(in), optional :: isfirst
  246. character(len=*), intent(in), optional :: dir
  247. character(len=*), intent(in), optional :: key
  248. !
  249. ! !OUTPUT PARAMETERS:
  250. !
  251. character(len=*), intent(out) :: fname
  252. integer, intent(out) :: status
  253. !
  254. ! !REVISION HISTORY:
  255. !
  256. !EOP
  257. !------------------------------------------------------------------------
  258. !BOC
  259. character(len=*), parameter :: rname = 'Restart_FileName'
  260. character(len=256) :: adir
  261. character(len=32) :: akey
  262. ! --- begin --------------------------------
  263. ! destination directory:
  264. adir = trim(outdir)
  265. if ( present(dir) ) adir = trim(dir)
  266. ! extra key, for example '_x' to denote that
  267. ! a restart file was dumped after process 'x':
  268. akey = ''
  269. if ( present(key) ) akey = trim(key)
  270. ! if this is the initial time, add an extra key to avoid
  271. ! that the restart file for this hour from the previous
  272. ! run is overwritten:
  273. if ( present(isfirst) ) then
  274. if ( isfirst ) akey = trim(akey)//'_initial'
  275. end if
  276. ! write filename:
  277. write (fname,'(a,"/TM5_restart_",i4.4,2i2.2,"_",2i2.2,"_",a,a,".nc")') &
  278. trim(adir), idate(1:5), trim(lli(region)%name), trim(akey)
  279. ! ok
  280. status = 0
  281. END SUBROUTINE RESTART_FILENAME
  282. !EOC
  283. !--------------------------------------------------------------------------
  284. ! TM5 !
  285. !--------------------------------------------------------------------------
  286. !BOP
  287. !
  288. ! !IROUTINE: RESTART_WRITE
  289. !
  290. ! !DESCRIPTION: write restart
  291. !\\
  292. !\\
  293. ! !INTERFACE:
  294. !
  295. SUBROUTINE RESTART_WRITE( status, key, region, isfirst )
  296. !
  297. ! !USES:
  298. !
  299. use GO , only : Get
  300. use dims , only : nregions, at, bt
  301. use chem_param , only : ntracet, ntrace_chem, ntrace, names
  302. use partools , only : isRoot
  303. use tm5_distgrid, only : dgrid, Get_DistGrid, gather
  304. use global_data , only : mass_dat, chem_dat
  305. use meteodata , only : global_lli, levi
  306. use meteodata , only : sp_dat, phlb_dat, m_dat
  307. use MDF , only : MDF_Create, MDF_EndDef, MDF_Close
  308. use MDF , only : MDF_Def_Dim, MDF_Def_Var
  309. use MDF , only : MDF_Put_Att, MDF_Put_Var
  310. use MDF , only : MDF_REPLACE, MDF_NETCDF4
  311. use MDF , only : MDF_FLOAT, MDF_DOUBLE, MDF_CHAR
  312. !
  313. ! !OUTPUT PARAMETERS:
  314. !
  315. integer, intent(out) :: status
  316. !
  317. ! !INPUT PARAMETERS:
  318. !
  319. character(len=*), intent(in), optional :: key
  320. integer, intent(in), optional :: region
  321. logical, intent(in), optional :: isfirst
  322. !
  323. ! !REVISION HISTORY:
  324. ! 8 Apr 2011 - P. Le Sager - use IF_NOTOK_MDF macro
  325. ! 10 Oct 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  326. !
  327. ! !REMARKS:
  328. !
  329. !EOP
  330. !------------------------------------------------------------------------
  331. !BOC
  332. character(len=*), parameter :: rname = 'Restart_Write'
  333. integer :: imr, jmr, lmr, n
  334. character(len=256) :: fname
  335. integer :: ftype
  336. integer :: ncid
  337. integer :: dimid_lon, dimid_lat, dimid_lev, dimid_hlev
  338. integer :: dimid_lon_sfc, dimid_lat_sfc
  339. integer :: dimid_trace, dimid_trace_transp, dimid_trace_chem
  340. integer :: dimid_name
  341. integer :: varid, varid_at, varid_bt
  342. integer :: varid_sp, varid_ph, varid_m
  343. integer :: varid_names, varid_rm
  344. #ifdef slopes
  345. integer :: varid_rxm, varid_rym, varid_rzm
  346. #endif
  347. integer :: varid_rmc
  348. integer :: rtype
  349. real, allocatable :: arr4d(:,:,:,:), arr3d(:,:,:)
  350. ! --- begin --------------------------------
  351. write (gol,'("write restart file(s) ...")'); call goPr
  352. ! loop over regions:
  353. REG: do n = 1, nregions
  354. ! only selected region ?
  355. if ( present(region) ) then
  356. if ( n /= region ) cycle
  357. end if
  358. ! entire region grid size
  359. imr = global_lli(n)%nlon
  360. jmr = global_lli(n)%nlat
  361. lmr = levi%nlev
  362. ! allocate 3D and 4D global arrays for gathering data
  363. if (isRoot) then
  364. allocate( arr4d(imr,jmr,lmr,ntracet) )
  365. allocate( arr3d(imr,jmr,lmr+1) )
  366. else
  367. allocate( arr4d(1,1,1,1) )
  368. allocate( arr3d(1,1,1) )
  369. endif
  370. ! name of restart file
  371. call Restart_FileName( n, fname, status, key=key, dir=rs_write_dir, isfirst=isfirst )
  372. IF_NOTOK_RETURN(status=1)
  373. write (gol,'(" destination : ",a)') trim(fname); call goPr
  374. if (isRoot) then
  375. !------------------
  376. ! OPEN NETCDF FILE
  377. !------------------
  378. ! overwrite existing files (clobber)
  379. call MDF_Create( trim(fname), MDF_NETCDF4, MDF_REPLACE, ncid, status )
  380. IF_NOTOK_RETURN(status=1)
  381. !------------------
  382. ! DEFINE DIMENSIONS
  383. !------------------
  384. call MDF_Def_Dim( ncid, 'lon', imr, dimid_lon, status )
  385. IF_NOTOK_MDF(fid=ncid)
  386. call MDF_Def_Dim( ncid, 'lat', jmr, dimid_lat, status )
  387. IF_NOTOK_MDF(fid=ncid)
  388. call MDF_Def_Dim( ncid, 'lev', lmr, dimid_lev, status )
  389. IF_NOTOK_MDF(fid=ncid)
  390. call MDF_Def_Dim( ncid, 'hlev', lmr+1, dimid_hlev, status )
  391. IF_NOTOK_MDF(fid=ncid)
  392. call MDF_Def_Dim( ncid, 'trace_transp', ntracet, dimid_trace_transp, status )
  393. IF_NOTOK_MDF(fid=ncid)
  394. if ( ntrace_chem > 0 ) then
  395. call MDF_Def_Dim( ncid, 'trace_chem', ntrace_chem, dimid_trace_chem, status )
  396. IF_NOTOK_MDF(fid=ncid)
  397. else
  398. dimid_trace_chem = -1
  399. end if
  400. call MDF_Def_Dim( ncid, 'trace', ntrace, dimid_trace, status )
  401. IF_NOTOK_MDF(fid=ncid)
  402. call MDF_Def_Dim( ncid, 'name', len(names(1)), dimid_name, status )
  403. IF_NOTOK_MDF(fid=ncid)
  404. !------------------
  405. ! DEFINE VARIABLES
  406. !------------------
  407. select case ( kind(m_dat(n)%data) )
  408. case ( 4 ) ; rtype = MDF_FLOAT
  409. case ( 8 ) ; rtype = MDF_DOUBLE
  410. case default
  411. write (gol,'("unsupported real kind : ",i6)') kind(m_dat(n)%data)
  412. TRACEBACK; status=1; return
  413. end select
  414. ! surface pressure
  415. call MDF_Def_Var( ncid, 'sp', rtype, (/dimid_lon,dimid_lat/), varid, status )
  416. IF_NOTOK_MDF(fid=ncid)
  417. call MDF_Put_Att( ncid, varid, 'long_name', 'surface pressure', status )
  418. IF_NOTOK_MDF(fid=ncid)
  419. call MDF_Put_Att( ncid, varid, 'unit', 'Pa', status )
  420. IF_NOTOK_MDF(fid=ncid)
  421. varid_sp = varid
  422. ! at, bt coefficients for hybrid grid
  423. call MDF_Def_Var( ncid, 'at', rtype, (/dimid_hlev/), varid, status )
  424. IF_NOTOK_MDF(fid=ncid)
  425. call MDF_Put_Att( ncid, varid, 'long_name', 'hybrid grid a_t coefficient', status )
  426. IF_NOTOK_MDF(fid=ncid)
  427. varid_at = varid
  428. call MDF_Def_Var( ncid, 'bt', rtype, (/dimid_hlev/), varid, status )
  429. IF_NOTOK_MDF(fid=ncid)
  430. call MDF_Put_Att( ncid, varid, 'long_name', 'hybrid grid b_t coefficient', status )
  431. IF_NOTOK_MDF(fid=ncid)
  432. varid_bt = varid
  433. ! half level pressure
  434. call MDF_Def_Var( ncid, 'ph', rtype, &
  435. (/dimid_lon,dimid_lat,dimid_hlev/), varid, status )
  436. IF_NOTOK_MDF(fid=ncid)
  437. call MDF_Put_Att( ncid, varid, 'long_name', 'half level pressure', status )
  438. IF_NOTOK_MDF(fid=ncid)
  439. call MDF_Put_Att( ncid, varid, 'unit', 'Pa', status )
  440. IF_NOTOK_MDF(fid=ncid)
  441. varid_ph = varid
  442. ! air mass
  443. call MDF_Def_Var( ncid, 'm', rtype, &
  444. (/dimid_lon,dimid_lat,dimid_lev/), varid, status )
  445. IF_NOTOK_MDF(fid=ncid)
  446. call MDF_Put_Att( ncid, varid, 'long_name', 'air mass', status )
  447. IF_NOTOK_MDF(fid=ncid)
  448. call MDF_Put_Att( ncid, varid, 'unit', 'kg', status )
  449. IF_NOTOK_MDF(fid=ncid)
  450. varid_m = varid
  451. !! accumulated surface fluxes
  452. !!
  453. !call MDF_Def_Var( ncid, 'slhf', rtype, (/dimid_lon_sfc,dimid_lat_sfc/), varid, status )
  454. !IF_NOTOK_MDF(fid=ncid)
  455. !call MDF_Put_Att( ncid, varid, 'long_name', 'surface latent heat flux', status )
  456. !IF_NOTOK_MDF(fid=ncid)
  457. !call MDF_Put_Att( ncid, varid, 'unit', 'W/m2', status )
  458. !IF_NOTOK_MDF(fid=ncid)
  459. !varid_slhf = varid
  460. !!
  461. !call MDF_Def_Var( ncid, 'sshf', rtype, (/dimid_lon_sfc,dimid_lat_sfc/), varid, status )
  462. !IF_NOTOK_MDF(fid=ncid)
  463. !call MDF_Put_Att( ncid, varid, 'long_name', 'surface sensible heat flux', status )
  464. !IF_NOTOK_MDF(fid=ncid)
  465. !call MDF_Put_Att( ncid, varid, 'unit', 'W/m2', status )
  466. !IF_NOTOK_MDF(fid=ncid)
  467. !varid_sshf = varid
  468. ! tracer names
  469. call MDF_Def_Var( ncid, 'names', MDF_CHAR, (/dimid_name,dimid_trace/), varid, status )
  470. IF_NOTOK_MDF(fid=ncid)
  471. call MDF_Put_Att( ncid, varid, 'long_name', 'tracer names', status )
  472. IF_NOTOK_MDF(fid=ncid)
  473. varid_names = varid
  474. ! tracer mass
  475. call MDF_Def_Var( ncid, 'rm', rtype, &
  476. (/dimid_lon,dimid_lat,dimid_lev,dimid_trace_transp/), varid, status )
  477. IF_NOTOK_MDF(fid=ncid)
  478. call MDF_Put_Att( ncid, varid, 'long_name', 'transported tracer mass', status )
  479. IF_NOTOK_MDF(fid=ncid)
  480. call MDF_Put_Att( ncid, varid, 'unit', 'kg', status )
  481. IF_NOTOK_MDF(fid=ncid)
  482. varid_rm = varid
  483. ! tracer mass slopes:
  484. #ifdef slopes
  485. call MDF_Def_Var( ncid, 'rxm', rtype, &
  486. (/dimid_lon,dimid_lat,dimid_lev,dimid_trace_transp/), varid, status )
  487. IF_NOTOK_MDF(fid=ncid)
  488. call MDF_Put_Att( ncid, varid, 'long_name', 'tracer mass slope in x direction', status )
  489. IF_NOTOK_MDF(fid=ncid)
  490. call MDF_Put_Att( ncid, varid, 'unit', 'kg/(half cell)', status )
  491. IF_NOTOK_MDF(fid=ncid)
  492. varid_rxm = varid
  493. call MDF_Def_Var( ncid, 'rym', rtype, &
  494. (/dimid_lon,dimid_lat,dimid_lev,dimid_trace_transp/), varid, status )
  495. IF_NOTOK_MDF(fid=ncid)
  496. call MDF_Put_Att( ncid, varid, 'long_name', 'tracer mass slope in y direction', status )
  497. IF_NOTOK_MDF(fid=ncid)
  498. call MDF_Put_Att( ncid, varid, 'unit', 'kg/(half cell)', status )
  499. IF_NOTOK_MDF(fid=ncid)
  500. varid_rym = varid
  501. call MDF_Def_Var( ncid, 'rzm', rtype, &
  502. (/dimid_lon,dimid_lat,dimid_lev,dimid_trace_transp/), varid, status )
  503. IF_NOTOK_MDF(fid=ncid)
  504. call MDF_Put_Att( ncid, varid, 'long_name', 'tracer mass slope in z direction', status )
  505. IF_NOTOK_MDF(fid=ncid)
  506. call MDF_Put_Att( ncid, varid, 'unit', 'kg/(half cell)', status )
  507. IF_NOTOK_MDF(fid=ncid)
  508. varid_rzm = varid
  509. #endif
  510. ! non-transported tracers:
  511. if ( ntrace_chem > 0 ) then
  512. call MDF_Def_Var( ncid, 'rmc', rtype, &
  513. (/dimid_lon,dimid_lat,dimid_lev,dimid_trace_chem/), varid, status )
  514. IF_NOTOK_MDF(fid=ncid)
  515. call MDF_Put_Att( ncid, varid, 'long_name', 'non-transported tracer mass', status )
  516. IF_NOTOK_MDF(fid=ncid)
  517. call MDF_Put_Att( ncid, varid, 'unit', 'kg', status )
  518. IF_NOTOK_MDF(fid=ncid)
  519. varid_rmc = varid
  520. end if
  521. !------------------
  522. ! END DEFINITION MODE
  523. !------------------
  524. call MDF_EndDef( ncid, status )
  525. IF_NOTOK_MDF(fid=ncid)
  526. endif
  527. !------------------
  528. ! WRITE VARIABLES
  529. !------------------
  530. ! surface pressure
  531. call gather( dgrid(n), sp_dat(n)%data, arr3d(:,:,1:1), sp_dat(n)%halo, status)
  532. IF_NOTOK_RETURN(status=1)
  533. if (isRoot) call MDF_Put_Var( ncid, varid_sp, arr3d(:,:,1), status )
  534. IF_NOTOK_MDF(fid=ncid)
  535. ! half level pressure
  536. call gather( dgrid(n), phlb_dat(n)%data, arr3d, phlb_dat(n)%halo, status)
  537. IF_NOTOK_RETURN(status=1)
  538. if (isRoot) call MDF_Put_Var( ncid, varid_ph, arr3d, status)
  539. IF_NOTOK_MDF(fid=ncid)
  540. ! at, bt coefficients
  541. if (isRoot) then
  542. call MDF_Put_Var( ncid, varid_at, at(1:lmr+1), status )
  543. IF_NOTOK_MDF(fid=ncid)
  544. call MDF_Put_Var( ncid, varid_bt, bt(1:lmr+1), status )
  545. IF_NOTOK_MDF(fid=ncid)
  546. end if
  547. ! air mass
  548. call gather( dgrid(n), m_dat(n)%data, arr4d(:,:,:,1), m_dat(n)%halo, status)
  549. IF_NOTOK_RETURN(status=1)
  550. if (isRoot) call MDF_Put_Var( ncid, varid_m, arr4d(:,:,:,1), status)
  551. IF_NOTOK_MDF(fid=ncid)
  552. !! surface latent heat flux; global surface field !
  553. !call MDF_Put_Var( ncid, varid_slhf, slhf_dat(iglbsfc)%data(1:n360,1:n180,1), status )
  554. !IF_NOTOK_MDF(fid=ncid)
  555. !
  556. !! surface sensible heat flux; global surface field !
  557. !call MDF_Put_Var( ncid, varid_sshf, sshf_dat(iglbsfc)%data(1:n360,1:n180,1), status )
  558. !IF_NOTOK_MDF(fid=ncid)
  559. ! tracer names
  560. if (isRoot) call MDF_Put_Var( ncid, varid_names, names, status )
  561. IF_NOTOK_MDF(fid=ncid)
  562. ! write transported tracers
  563. call gather( dgrid(n), mass_dat(n)%rm, arr4d, mass_dat(n)%halo, status)
  564. IF_NOTOK_RETURN(status=1)
  565. if (isRoot) call MDF_Put_Var( ncid, varid_rm, arr4d, status)
  566. IF_NOTOK_MDF(fid=ncid)
  567. #ifdef slopes
  568. call gather( dgrid(n), mass_dat(n)%rxm, arr4d, mass_dat(n)%halo, status)
  569. IF_NOTOK_RETURN(status=1)
  570. if (isRoot) call MDF_Put_Var( ncid, varid_rxm, arr4d, status)
  571. IF_NOTOK_MDF(fid=ncid)
  572. call gather( dgrid(n), mass_dat(n)%rym, arr4d, mass_dat(n)%halo, status)
  573. IF_NOTOK_RETURN(status=1)
  574. if (isRoot) call MDF_Put_Var( ncid, varid_rym, arr4d, status)
  575. IF_NOTOK_MDF(fid=ncid)
  576. call gather( dgrid(n), mass_dat(n)%rzm, arr4d, mass_dat(n)%halo, status)
  577. IF_NOTOK_RETURN(status=1)
  578. if (isRoot) call MDF_Put_Var( ncid, varid_rzm, arr4d, status)
  579. IF_NOTOK_MDF(fid=ncid)
  580. #endif
  581. ! write transported tracers
  582. call gather( dgrid(n), mass_dat(n)%rm, arr4d, mass_dat(n)%halo, status)
  583. IF_NOTOK_RETURN(status=1)
  584. if (isRoot) call MDF_Put_Var( ncid, varid_rm, arr4d, status)
  585. IF_NOTOK_MDF(fid=ncid)
  586. #ifdef slopes
  587. call gather( dgrid(n), mass_dat(n)%rxm, arr4d, mass_dat(n)%halo, status)
  588. IF_NOTOK_RETURN(status=1)
  589. if (isRoot) call MDF_Put_Var( ncid, varid_rxm, arr4d, status)
  590. IF_NOTOK_MDF(fid=ncid)
  591. call gather( dgrid(n), mass_dat(n)%rym, arr4d, mass_dat(n)%halo, status)
  592. IF_NOTOK_RETURN(status=1)
  593. if (isRoot) call MDF_Put_Var( ncid, varid_rym, arr4d, status)
  594. IF_NOTOK_MDF(fid=ncid)
  595. call gather( dgrid(n), mass_dat(n)%rzm, arr4d, mass_dat(n)%halo, status)
  596. IF_NOTOK_RETURN(status=1)
  597. if (isRoot) call MDF_Put_Var( ncid, varid_rzm, arr4d, status)
  598. IF_NOTOK_MDF(fid=ncid)
  599. #endif
  600. ! write non-transported tracers
  601. if (ntrace_chem > 0) then
  602. call gather( dgrid(n), chem_dat(n)%rm, arr4d(:,:,:,1:ntrace_chem), chem_dat(n)%halo, status)
  603. IF_NOTOK_RETURN(status=1)
  604. if (isRoot) call MDF_Put_Var( ncid, varid_rmc, arr4d(:,:,:,1:ntrace_chem), status)
  605. IF_NOTOK_MDF(fid=ncid)
  606. end if
  607. ! Done
  608. if (isRoot) call MDF_Close( ncid, status )
  609. IF_NOTOK_RETURN(status=1)
  610. deallocate(arr4d, arr3d)
  611. end do REG
  612. status = 0
  613. END SUBROUTINE RESTART_WRITE
  614. !EOC
  615. !--------------------------------------------------------------------------
  616. ! TM5 !
  617. !--------------------------------------------------------------------------
  618. !BOP
  619. !
  620. ! !IROUTINE: RESTART_READ
  621. !
  622. ! !DESCRIPTION: Read restart file. Case of istart=33 (can read any of the
  623. ! available variables) or 32 (can read only tracer mass).
  624. !\\
  625. !\\
  626. ! !INTERFACE:
  627. !
  628. SUBROUTINE RESTART_READ( status, region, &
  629. surface_pressure, pressure, air_mass, surface_fluxes, &
  630. tracer_mass, tendencies, megan_history, nox_pulsing )
  631. !
  632. ! !USES:
  633. !
  634. use GO, only : TrcFile, Init, Done, ReadRc
  635. use GO, only : goMatchValue
  636. use dims, only : nregions, im, jm, istart, idate, idatei
  637. use grid, only : TllGridInfo, TLevelInfo, Init, Done, Fill3D
  638. use chem_param, only : ntracet, ntrace_chem, ntrace
  639. use chem_param, only : names, tracer_name_len
  640. use partools, only : isRoot, par_broadcast
  641. use tm5_distgrid, only : dgrid, gather, scatter
  642. use global_data, only : rcfile, mass_dat, chem_dat
  643. use meteodata, only : levi, global_lli, sp_dat, phlb_dat, m_dat
  644. !use meteodata, only : slhf_dat, sshf_dat
  645. use MDF, only : MDF_Open, MDF_Close, MDF_Inquire_Dimension
  646. use MDF, only : MDF_Inq_VarID, MDF_Inquire_Variable, MDF_Inq_DimID
  647. use MDf, only : MDF_Var_Par_Access, MDF_INDEPENDENT, MDF_COLLECTIVE
  648. use MDF, only : MDF_Get_Att, MDF_Get_Var
  649. use MDF, only : MDF_READ, MDF_NETCDF4
  650. !
  651. ! !OUTPUT PARAMETERS:
  652. !
  653. integer, intent(out) :: status
  654. !
  655. ! !INPUT PARAMETERS:
  656. !
  657. integer, intent(in), optional :: region
  658. logical, intent(in), optional :: surface_pressure, pressure, air_mass, surface_fluxes
  659. logical, intent(in), optional :: tracer_mass, tendencies, megan_history, nox_pulsing
  660. !
  661. ! !REVISION HISTORY:
  662. ! 8 Apr 2011 - P. Le Sager - use IF_NOTOK_MDF macro
  663. ! 28 Apr 2011 - P. Le Sager - Check on tracer availability in restart file.
  664. ! - Allows for more tracers in restart file than needed
  665. ! 10 May 2011 - P. Le Sager - Added deallocate statement to work with zoom regions
  666. ! 10 Oct 2011 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  667. !
  668. ! !REMARKS:
  669. ! - logically, if we need to remap, then meteo is not read from restart
  670. ! (but from met field and used for remapping): in other words, only
  671. ! tracer mass is read, and istart should be 32.
  672. !
  673. !EOP
  674. !------------------------------------------------------------------------
  675. !BOC
  676. character(len=*), parameter :: rname = mname//'/Restart_Read'
  677. character(len=tracer_name_len), allocatable :: values_names(:)
  678. character(len=256) :: rs_read_dir, fname
  679. type(TrcFile) :: rcF
  680. logical :: exist
  681. logical :: do_sp, do_ph, do_m, do_sflux, do_rm, do_megan, do_pulse
  682. integer :: imr, jmr, lmr, imr_restart, jmr_restart, lmr_restart
  683. integer :: n, region1, region2
  684. integer :: ncid
  685. integer :: varid_sp, varid_ph, varid_m, varid_rm, varid_rmc, varid_names
  686. !integer :: varid_slhf, varid_sshf
  687. integer :: itr, itr_file
  688. integer :: ntracet_restart, dimid
  689. integer :: shp(2)
  690. #ifdef slopes
  691. integer :: varid_rxm, varid_rym, varid_rzm
  692. #endif
  693. ! global work arrays to read data
  694. real, allocatable :: tmp3d(:,:,:)
  695. real, allocatable :: rmt(:,:,:,:),rms(:,:,:,:), rmx(:,:,:,:),rmy(:,:,:,:), rmz(:,:,:,:)
  696. ! for remapping:
  697. logical :: need_vremap, need_hremap, need_remap
  698. integer :: varid_at, varid_bt
  699. real :: dx, dy
  700. real, allocatable :: sp_gbl(:,:,:)
  701. real, allocatable :: at_restart(:), bt_restart(:)
  702. real, allocatable :: src_glb(:,:,:)
  703. type(TllGridInfo) :: lli_restart
  704. type(TLevelInfo) :: levi_restart
  705. ! --- begin --------------------------------
  706. if ( istart /= 33 .and. istart /= 32 ) then
  707. write (gol,'(" skip read restart; istart not 33 or 32 but ",i2)') istart; call goPr
  708. status=0; return
  709. endif
  710. if ( any( idate /= idatei ) ) then
  711. write (gol,'(" skip read restart; idate not idatei but ",i4,5i2.2)') idate; call goPr
  712. status=0; return
  713. endif
  714. ! input directory:
  715. call Init( rcF, rcfile, status )
  716. IF_NOTOK_RETURN(status=1)
  717. call ReadRc( rcF, 'restart.read.dir', rs_read_dir, status )
  718. IF_NOTOK_RETURN(status=1)
  719. call Done( rcF, status )
  720. IF_NOTOK_RETURN(status=1)
  721. ! region range:
  722. if ( present(region) ) then
  723. region1 = region
  724. region2 = region
  725. else
  726. region1 = 1
  727. region2 = nregions
  728. end if
  729. ! data sets:
  730. do_rm = .false. ; if ( present(tracer_mass ) ) do_rm = tracer_mass
  731. do_m = .false. ; if ( present(air_mass ).and.(istart==33) ) do_m = air_mass
  732. do_sp = .false. ; if ( present(surface_pressure ).and.(istart==33) ) do_sp = surface_pressure
  733. do_ph = .false. ; if ( present(pressure ).and.(istart==33) ) do_ph = pressure
  734. do_sflux = .false. ; if ( present(surface_fluxes ).and.(istart==33) ) do_sflux = surface_fluxes
  735. do_megan = .false. ; if ( present(megan_history ).and.(istart==33) ) do_megan = megan_history
  736. do_pulse = .false. ; if ( present(nox_pulsing ).and.(istart==33) ) do_pulse = nox_pulsing
  737. ! sorry ..
  738. if ( do_sflux ) then
  739. write (gol,'("no surface fluxes in restart files until somebody")') ; call goErr
  740. write (gol,'("has a good idea on what should be storred:")') ; call goErr
  741. write (gol,'(" o global surface field (1x1 ?)")') ; call goErr
  742. write (gol,'(" o zoom regions")') ; call goErr
  743. write (gol,'(" o both")') ; call goErr
  744. TRACEBACK; status=1; return
  745. end if
  746. ! do we need anything?
  747. if(.not.(do_rm.or.do_m.or.do_sp.or.do_ph.or.do_sflux.or.do_megan.or.do_pulse))then
  748. status=0; return
  749. endif
  750. REG: do n = region1, region2
  751. imr = global_lli(n)%nlon
  752. jmr = global_lli(n)%nlat
  753. lmr = levi%nlev
  754. ! name of restart file
  755. call Restart_FileName( n, fname, status, dir=trim(rs_read_dir) )
  756. IF_NOTOK_RETURN(status=1)
  757. write (gol,'(" read restart file: ",a)') trim(fname); call goPr
  758. inquire( file=fname, exist=exist )
  759. if ( .not. exist ) then
  760. write (gol,'("restart file not found : ",a)') trim(fname); call goErr
  761. TRACEBACK; status=1; return
  762. end if
  763. ! ** open netcdf file
  764. if (isRoot) then
  765. call MDF_Open( trim(fname), MDF_NETCDF4, MDF_READ, ncid, status )
  766. IF_NOTOK_RETURN(status=1)
  767. ! ** check for dimension compatibility
  768. call MDF_Inq_DimID( ncid, 'lev', dimid, status )
  769. IF_NOTOK_MDF(fid=ncid)
  770. call MDF_Inquire_Dimension( ncid, dimid, status, length=lmr_restart )
  771. IF_NOTOK_MDF(fid=ncid)
  772. call MDF_Inq_DimID( ncid, 'lat', dimid, status )
  773. IF_NOTOK_MDF(fid=ncid)
  774. call MDF_Inquire_Dimension( ncid, dimid, status, length=jmr_restart )
  775. IF_NOTOK_MDF(fid=ncid)
  776. call MDF_Inq_DimID( ncid, 'lon', dimid, status )
  777. IF_NOTOK_MDF(fid=ncid)
  778. call MDF_Inquire_Dimension( ncid, dimid, status, length=imr_restart )
  779. IF_NOTOK_MDF(fid=ncid)
  780. need_vremap = (lmr /= lmr_restart)
  781. need_hremap = (jmr /= jmr_restart) .or. (imr /= imr_restart)
  782. need_remap = need_hremap .or. need_vremap
  783. endif
  784. call par_broadcast( need_remap, status)
  785. IF_NOTOK_RETURN(status=1)
  786. if ((istart /= 32).and.need_remap) then
  787. status=1
  788. write(gol,*)' you must use istart=32 for using a restart file at different resolution'
  789. call goErr
  790. TRACEBACK; return
  791. endif
  792. ! work arrays
  793. if (isRoot) then
  794. allocate( rmt(imr,jmr,lmr,ntracet) )
  795. allocate( rmx(imr,jmr,lmr,ntracet) )
  796. allocate( rmy(imr,jmr,lmr,ntracet) )
  797. allocate( rmz(imr,jmr,lmr,ntracet) )
  798. if ( ntrace_chem > 0 ) allocate( rms(imr,jmr,lmr,ntracet+1:ntracet+ntrace_chem) )
  799. allocate( tmp3d(imr,jmr,lmr+1 ) )
  800. else
  801. allocate( rmt(1,1,1,1) )
  802. allocate( rmx(1,1,1,1) )
  803. allocate( rmy(1,1,1,1) )
  804. allocate( rmz(1,1,1,1) )
  805. if ( ntrace_chem > 0 ) allocate( rms(1,1,1,1) )
  806. allocate( tmp3d(1,1,1) )
  807. endif
  808. ! prepare for remap
  809. if (need_remap .and. do_rm) then
  810. write (gol,'(" remap tracer from restart file")') ; call goPr
  811. if (isRoot) then
  812. allocate( sp_gbl(imr,jmr,1) )
  813. allocate( src_glb(imr_restart,jmr_restart,lmr_restart))
  814. else
  815. allocate(sp_gbl(1,1,1))
  816. allocate(src_glb(1,1,1))
  817. endif
  818. call gather( dgrid(n), sp_dat(n)%data, sp_gbl, sp_dat(n)%halo, status)
  819. IF_NOTOK_RETURN(status=1)
  820. ! init to 0 in case of data not found in file
  821. rmt=0.
  822. rms=0.
  823. ! init lli_restart/levi_restart
  824. allocate(at_restart(lmr_restart+1))
  825. allocate(bt_restart(lmr_restart+1))
  826. !
  827. call MDF_Inq_VarID( ncid, 'at', varid_at, status )
  828. IF_NOTOK_MDF(fid=ncid)
  829. !
  830. call MDF_Get_Var( ncid, varid_at, at_restart(1:(lmr_restart+1)), status )
  831. IF_NOTOK_MDF(fid=ncid)
  832. !
  833. call MDF_Inq_VarID( ncid, 'bt', varid_bt, status )
  834. IF_NOTOK_MDF(fid=ncid)
  835. !
  836. call MDF_Get_Var( ncid, varid_bt, bt_restart(1:(lmr_restart+1)), status )
  837. IF_NOTOK_MDF(fid=ncid)
  838. !
  839. call Init( levi_restart, lmr_restart, at_restart, bt_restart, status )
  840. IF_NOTOK_RETURN(status=1)
  841. !
  842. deallocate(at_restart,bt_restart)
  843. !
  844. dx=360./imr_restart
  845. dy=180./jmr_restart
  846. call Init( lli_restart, -180.+0.5*dx, dx, imr_restart, &
  847. -90.+0.5*dy, dy, jmr_restart, status )
  848. IF_NOTOK_RETURN(status=1)
  849. endif
  850. ! ** get variables id
  851. if (isRoot) then
  852. ! surface pressure
  853. if ( do_sp ) call MDF_Inq_VarID( ncid, 'sp', varid_sp, status )
  854. IF_NOTOK_MDF(fid=ncid)
  855. ! half level pressure
  856. if ( do_ph ) call MDF_Inq_VarID( ncid, 'ph', varid_ph, status )
  857. IF_NOTOK_MDF(fid=ncid)
  858. ! air mass
  859. if ( do_m ) call MDF_Inq_VarID( ncid, 'm', varid_m, status )
  860. IF_NOTOK_MDF(fid=ncid)
  861. !! surface fluxes
  862. !if ( do_sflux ) then
  863. !end if
  864. ! tracer mass
  865. if ( do_rm ) then
  866. call MDF_Inq_VarID( ncid, 'names', varid_names, status )
  867. if ( status /= 0 ) then
  868. write (gol,'("could not find variable `names` in restart file;")'); call goErr
  869. write (gol,'(" using an old restart file to initialize the model ?")'); call goErr
  870. status=1
  871. end if
  872. IF_NOTOK_MDF(fid=ncid)
  873. ! get dimension of "names"
  874. call MDF_Inquire_Variable( ncid, varid_names, status, shp=shp )
  875. IF_NOTOK_MDF(fid=ncid)
  876. ! get number of transported tracer in restart file
  877. call MDF_Inq_DimID( ncid, 'trace_transp', dimid, status )
  878. IF_NOTOK_MDF(fid=ncid)
  879. call MDF_Inquire_Dimension( ncid, dimid, status, length=ntracet_restart )
  880. IF_NOTOK_MDF(fid=ncid)
  881. ! tracers mass id
  882. call MDF_Inq_VarID( ncid, 'rm', varid_rm, status )
  883. IF_NOTOK_MDF(fid=ncid)
  884. #ifdef slopes
  885. call MDF_Inq_VarID( ncid, 'rxm', varid_rxm, status )
  886. IF_NOTOK_MDF(fid=ncid)
  887. call MDF_Inq_VarID( ncid, 'rym', varid_rym, status )
  888. IF_NOTOK_MDF(fid=ncid)
  889. call MDF_Inq_VarID( ncid, 'rzm', varid_rzm, status )
  890. IF_NOTOK_MDF(fid=ncid)
  891. #endif
  892. ! read non-transported tracers if any
  893. if ( ntrace_chem > 0 ) then
  894. call MDF_Inq_VarID( ncid, 'rmc', varid_rmc, status )
  895. IF_NOTOK_MDF(fid=ncid)
  896. end if
  897. end if
  898. end if
  899. ! *** READ VARIABLES ***
  900. if ( do_sp ) then
  901. write (gol,'(" restore surface pressure ...")'); call goPr
  902. if (isRoot) call MDF_Get_Var( ncid, varid_sp, tmp3d(:,:,1), status )
  903. IF_NOTOK_MDF(fid=ncid)
  904. call scatter( dgrid(n), sp_dat(n)%data, tmp3d(:,:,1:1), sp_dat(n)%halo, status)
  905. IF_NOTOK_RETURN(status=1)
  906. end if
  907. if ( do_ph ) then
  908. write (gol,'(" restore half level pressure ...")'); call goPr
  909. if (isRoot) call MDF_Get_Var( ncid, varid_ph, tmp3d, status)
  910. IF_NOTOK_MDF(fid=ncid)
  911. call scatter( dgrid(n), phlb_dat(n)%data, tmp3d, phlb_dat(n)%halo, status)
  912. IF_NOTOK_RETURN(status=1)
  913. end if
  914. if ( do_m ) then
  915. write (gol,'(" restore air mass ...")'); call goPr
  916. if (isRoot) call MDF_Get_Var( ncid, varid_m, tmp3d(:,:,1:lmr), status)
  917. IF_NOTOK_MDF(fid=ncid)
  918. call scatter( dgrid(n), m_dat(n)%data, tmp3d(:,:,1:lmr), m_dat(n)%halo, status)
  919. IF_NOTOK_RETURN(status=1)
  920. end if
  921. ! tracer mass
  922. READRM: if ( do_rm ) then
  923. write (gol,'(" restore tracer mass ...")'); call goPr
  924. ! read list with tracer names in rcfile
  925. allocate( values_names(shp(2)) )
  926. call MDF_Get_Var( ncid, varid_names, values_names, status )
  927. IF_NOTOK_MDF(fid=ncid)
  928. ! loop over all model tracers
  929. do itr = 1, ntrace
  930. if (isRoot) then
  931. ! search in list:
  932. call goMatchValue( names(itr), values_names, itr_file, status )
  933. if ( status < 0 ) then
  934. write(gol,'("*WARNING* Requested tracer `",a,"` not FOUND in restart file!")') trim(names(itr))
  935. if (istart /= 32) then
  936. call goErr
  937. IF_NOTOK_MDF(fid=ncid)
  938. else
  939. status=0
  940. call goPr
  941. if ( itr <= ntracet ) then
  942. rmt(:,:,:,itr) = 1.e-30
  943. write(gol,'("*WARNING* Requested TRANSPORTED tracer `",a,"` has been SET to a default value of 1.e30")') trim(names(itr))
  944. else
  945. rms(:,:,:,itr) = 1.e-30
  946. write(gol,'("*WARNING* Requested SHORT-LIVED tracer `",a,"` has been SET to a default value of 1.e30")') trim(names(itr))
  947. endif
  948. call goPr
  949. endif
  950. else
  951. ! transported or short lived tracer ?
  952. if ( itr <= ntracet ) then
  953. if ( itr_file > ntracet_restart ) then
  954. write (gol,'("tracer `",a,"` is transported but seems to be not-transported in restart file")') trim(names(itr)); call goErr
  955. status=1
  956. IF_NOTOK_MDF(fid=ncid)
  957. end if
  958. if (need_remap) then
  959. call MDF_Get_Var( ncid, varid_rm, src_glb, status, start=(/1,1,1,itr_file/))
  960. IF_NOTOK_MDF(fid=ncid)
  961. call Fill3D( global_lli(n), levi, 'n', sp_gbl(:,:,1), rmt(:,:,:,itr), &
  962. lli_restart, levi_restart, src_glb, 'sum', status )
  963. IF_NOTOK_RETURN(status=1)
  964. else
  965. call MDF_Get_Var( ncid, varid_rm, rmt(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  966. IF_NOTOK_MDF(fid=ncid)
  967. endif
  968. #ifdef slopes
  969. ! read slopes
  970. if (.not. need_remap) then
  971. if (isRoot) call MDF_Get_Var( ncid, varid_rxm, rmx(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  972. IF_NOTOK_MDF(fid=ncid)
  973. if (isRoot) call MDF_Get_Var( ncid, varid_rym, rmy(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  974. IF_NOTOK_MDF(fid=ncid)
  975. if (isRoot) call MDF_Get_Var( ncid, varid_rzm, rmz(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  976. IF_NOTOK_MDF(fid=ncid)
  977. ! Scale methane concentration slopes by a factor specified in the rc file
  978. if ( (factor_ch4 /= 1.) .and. (itr == ich4) ) then
  979. mass_dat(n)%rxm(:,:,:,itr)= mass_dat(n)%rxm(:,:,:,itr) * factor_ch4
  980. mass_dat(n)%rym(:,:,:,itr)= mass_dat(n)%rym(:,:,:,itr) * factor_ch4
  981. mass_dat(n)%rzm(:,:,:,itr)= mass_dat(n)%rzm(:,:,:,itr) * factor_ch4
  982. endif
  983. endif
  984. #endif
  985. else ! short lived tracer:
  986. if ( itr_file <= ntracet_restart ) then
  987. write (gol,'("tracer `",a,"` is not-transported but seems to be transported in restart file")') trim(names(itr)); call goErr
  988. status=1
  989. IF_NOTOK_MDF(fid=ncid)
  990. end if
  991. itr_file = itr_file - ntracet_restart
  992. if (need_remap) then
  993. call MDF_Get_Var( ncid, varid_rmc, src_glb, status, start=(/1,1,1,itr_file/) )
  994. IF_NOTOK_MDF(fid=ncid)
  995. call Fill3D( global_lli(n), levi, 'n', sp_gbl(:,:,1), rms(:,:,:,itr), &
  996. lli_restart, levi_restart, src_glb, 'sum', status )
  997. IF_NOTOK_RETURN(status=1)
  998. else
  999. call MDF_Get_Var( ncid, varid_rmc, rms(:,:,:,itr), status, start=(/1,1,1,itr_file/) )
  1000. IF_NOTOK_MDF(fid=ncid)
  1001. endif
  1002. end if ! transported or short-lived
  1003. endif ! in the file
  1004. endif ! root
  1005. end do ! tracers
  1006. ! distribute
  1007. call scatter( dgrid(n), mass_dat(n)%rm, rmt, mass_dat(n)%halo, status)
  1008. IF_NOTOK_RETURN(status=1)
  1009. if ( ntrace_chem > 0 ) then
  1010. call scatter( dgrid(n), chem_dat(n)%rm, rms, chem_dat(n)%halo, status)
  1011. IF_NOTOK_RETURN(status=1)
  1012. endif
  1013. #ifdef slopes
  1014. if (.not. need_remap) then
  1015. call scatter( dgrid(n), mass_dat(n)%rxm, rmx, mass_dat(n)%halo, status)
  1016. IF_NOTOK_RETURN(status=1)
  1017. call scatter( dgrid(n), mass_dat(n)%rym, rmy, mass_dat(n)%halo, status)
  1018. IF_NOTOK_RETURN(status=1)
  1019. call scatter( dgrid(n), mass_dat(n)%rzm, rmz, mass_dat(n)%halo, status)
  1020. IF_NOTOK_RETURN(status=1)
  1021. else
  1022. ! Ensure that slopes are initialized to "unset" values of 0.0. Wouter says that
  1023. ! we could remap levels for rxm et al., but 0.0 will also work. The noise
  1024. ! induced from remapping the rm array is almost certainly bigger than any issues
  1025. ! from having this "default=0.0" slopes information. -ARJ 1 Jan 12
  1026. mass_dat(n)%rxm = 0.0
  1027. mass_dat(n)%rym = 0.0
  1028. mass_dat(n)%rzm = 0.0
  1029. endif
  1030. #endif
  1031. ! free mem for next region
  1032. deallocate( values_names)
  1033. if (need_remap) deallocate(sp_gbl,src_glb)
  1034. ENDIF READRM
  1035. ! clean
  1036. deallocate(rmt)
  1037. if ( ntrace_chem > 0 ) deallocate(rms)
  1038. #ifdef slopes
  1039. deallocate(rmx, rmy, rmz)
  1040. #endif
  1041. deallocate( tmp3d )
  1042. if (isRoot) call MDF_Close( ncid, status )
  1043. IF_NOTOK_RETURN(status=1)
  1044. ENDDO REG
  1045. status = 0
  1046. END SUBROUTINE RESTART_READ
  1047. !EOC
  1048. END MODULE RESTART