tm5_restart.F90 45 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277
  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. ! - If we need to remap, then meteo is not read from restart.
  670. ! Airmass is still read but only to convert tracer masses to mixing ratios.
  671. ! 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(:,:,:), airmass(:,:,:), run_airmass(:,:,:)
  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. allocate( airmass(imr_restart,jmr_restart,lmr_restart ) )
  801. if (istart==32) allocate( run_airmass(imr,jmr,lmr) )
  802. else
  803. allocate( rmt(1,1,1,1) )
  804. allocate( rmx(1,1,1,1) )
  805. allocate( rmy(1,1,1,1) )
  806. allocate( rmz(1,1,1,1) )
  807. if ( ntrace_chem > 0 ) allocate( rms(1,1,1,1) )
  808. allocate( airmass(1,1,1) )
  809. if (istart==32) allocate( run_airmass(1,1,1) )
  810. allocate( tmp3d(1,1,1) )
  811. endif
  812. if (istart==32) then
  813. CALL GATHER( dgrid(n), m_dat(n)%data, run_airmass, m_dat(n)%halo, status )
  814. IF_NOTOK_RETURN(status=1)
  815. endif
  816. ! prepare for remap
  817. if (need_remap .and. do_rm) then
  818. write (gol,'(" remap tracer from restart file")') ; call goPr
  819. if (isRoot) then
  820. allocate( sp_gbl(imr,jmr,1) )
  821. allocate( src_glb(imr_restart,jmr_restart,lmr_restart))
  822. else
  823. allocate(sp_gbl(1,1,1))
  824. allocate(src_glb(1,1,1))
  825. endif
  826. call gather( dgrid(n), sp_dat(n)%data, sp_gbl, sp_dat(n)%halo, status)
  827. IF_NOTOK_RETURN(status=1)
  828. ! init to 0 in case of data not found in file
  829. rmt=0.
  830. rms=0.
  831. ! init lli_restart/levi_restart
  832. if (isRoot) then
  833. allocate(at_restart(lmr_restart+1))
  834. allocate(bt_restart(lmr_restart+1))
  835. !
  836. call MDF_Inq_VarID( ncid, 'at', varid_at, status )
  837. IF_NOTOK_MDF(fid=ncid)
  838. !
  839. call MDF_Get_Var( ncid, varid_at, at_restart(1:(lmr_restart+1)), status )
  840. IF_NOTOK_MDF(fid=ncid)
  841. !
  842. call MDF_Inq_VarID( ncid, 'bt', varid_bt, status )
  843. IF_NOTOK_MDF(fid=ncid)
  844. !
  845. call MDF_Get_Var( ncid, varid_bt, bt_restart(1:(lmr_restart+1)), status )
  846. IF_NOTOK_MDF(fid=ncid)
  847. !
  848. call Init( levi_restart, lmr_restart, at_restart, bt_restart, status )
  849. IF_NOTOK_RETURN(status=1)
  850. !
  851. deallocate(at_restart,bt_restart)
  852. !
  853. dx=360./imr_restart
  854. dy=180./jmr_restart
  855. call Init( lli_restart, -180.+0.5*dx, dx, imr_restart, &
  856. -90.+0.5*dy, dy, jmr_restart, status )
  857. IF_NOTOK_RETURN(status=1)
  858. endif
  859. endif
  860. ! ** get variables id
  861. if (isRoot) then
  862. ! surface pressure
  863. if ( do_sp ) call MDF_Inq_VarID( ncid, 'sp', varid_sp, status )
  864. IF_NOTOK_MDF(fid=ncid)
  865. ! half level pressure
  866. if ( do_ph ) call MDF_Inq_VarID( ncid, 'ph', varid_ph, status )
  867. IF_NOTOK_MDF(fid=ncid)
  868. ! air mass
  869. call MDF_Inq_VarID( ncid, 'm', varid_m, status )
  870. IF_NOTOK_MDF(fid=ncid)
  871. !! surface fluxes
  872. !if ( do_sflux ) then
  873. !end if
  874. ! tracer mass
  875. if ( do_rm ) then
  876. call MDF_Inq_VarID( ncid, 'names', varid_names, status )
  877. if ( status /= 0 ) then
  878. write (gol,'("could not find variable `names` in restart file;")'); call goErr
  879. write (gol,'(" using an old restart file to initialize the model ?")'); call goErr
  880. status=1
  881. end if
  882. IF_NOTOK_MDF(fid=ncid)
  883. ! get dimension of "names"
  884. call MDF_Inquire_Variable( ncid, varid_names, status, shp=shp )
  885. IF_NOTOK_MDF(fid=ncid)
  886. ! get number of transported tracer in restart file
  887. call MDF_Inq_DimID( ncid, 'trace_transp', dimid, status )
  888. IF_NOTOK_MDF(fid=ncid)
  889. call MDF_Inquire_Dimension( ncid, dimid, status, length=ntracet_restart )
  890. IF_NOTOK_MDF(fid=ncid)
  891. ! tracers mass id
  892. call MDF_Inq_VarID( ncid, 'rm', varid_rm, status )
  893. IF_NOTOK_MDF(fid=ncid)
  894. #ifdef slopes
  895. call MDF_Inq_VarID( ncid, 'rxm', varid_rxm, status )
  896. IF_NOTOK_MDF(fid=ncid)
  897. call MDF_Inq_VarID( ncid, 'rym', varid_rym, status )
  898. IF_NOTOK_MDF(fid=ncid)
  899. call MDF_Inq_VarID( ncid, 'rzm', varid_rzm, status )
  900. IF_NOTOK_MDF(fid=ncid)
  901. #endif
  902. ! read non-transported tracers if any
  903. if ( ntrace_chem > 0 ) then
  904. call MDF_Inq_VarID( ncid, 'rmc', varid_rmc, status )
  905. IF_NOTOK_MDF(fid=ncid)
  906. end if
  907. end if
  908. end if
  909. ! *** READ VARIABLES ***
  910. if ( do_sp ) then
  911. write (gol,'(" restore surface pressure ...")'); call goPr
  912. if (isRoot) call MDF_Get_Var( ncid, varid_sp, tmp3d(:,:,1), status )
  913. IF_NOTOK_MDF(fid=ncid)
  914. call scatter( dgrid(n), sp_dat(n)%data, tmp3d(:,:,1:1), sp_dat(n)%halo, status)
  915. IF_NOTOK_RETURN(status=1)
  916. end if
  917. if ( do_ph ) then
  918. write (gol,'(" restore half level pressure ...")'); call goPr
  919. if (isRoot) call MDF_Get_Var( ncid, varid_ph, tmp3d, status)
  920. IF_NOTOK_MDF(fid=ncid)
  921. call scatter( dgrid(n), phlb_dat(n)%data, tmp3d, phlb_dat(n)%halo, status)
  922. IF_NOTOK_RETURN(status=1)
  923. end if
  924. ! get air mass in all cases
  925. if (isRoot) call MDF_Get_Var( ncid, varid_m, airmass, status)
  926. IF_NOTOK_MDF(fid=ncid)
  927. if ( do_m ) then
  928. write (gol,'(" restore air mass ...")'); call goPr
  929. call scatter( dgrid(n), m_dat(n)%data, airmass, m_dat(n)%halo, status)
  930. IF_NOTOK_RETURN(status=1)
  931. end if
  932. ! tracer mass
  933. READRM: if ( do_rm ) then
  934. write (gol,'(" restore tracer mass ...")'); call goPr
  935. ! read list with tracer names in rcfile
  936. allocate( values_names(shp(2)) )
  937. if (isRoot) call MDF_Get_Var( ncid, varid_names, values_names, status )
  938. IF_NOTOK_MDF(fid=ncid)
  939. ! loop over all model tracers
  940. do itr = 1, ntrace
  941. if (isRoot) then
  942. ! search in list:
  943. call goMatchValue( names(itr), values_names, itr_file, status )
  944. if ( status < 0 ) then
  945. write(gol,'("*WARNING* Requested tracer `",a,"` not FOUND in restart file!")') trim(names(itr))
  946. if (istart /= 32) then
  947. call goErr
  948. IF_NOTOK_MDF(fid=ncid)
  949. else
  950. status=0
  951. call goPr
  952. if ( itr <= ntracet ) then
  953. rmt(:,:,:,itr) = 1.e-30
  954. write(gol,'("*WARNING* Requested TRANSPORTED tracer `",a,"` has been SET to a default value of 1.e-30")') trim(names(itr))
  955. else
  956. rms(:,:,:,itr) = 1.e-30
  957. write(gol,'("*WARNING* Requested SHORT-LIVED tracer `",a,"` has been SET to a default value of 1.e-30")') trim(names(itr))
  958. endif
  959. call goPr
  960. endif
  961. else
  962. ! transported or short lived tracer ?
  963. if ( itr <= ntracet ) then
  964. if ( itr_file > ntracet_restart ) then
  965. write (gol,'("tracer `",a,"` is transported but seems to be not-transported in restart file")') trim(names(itr)); call goErr
  966. status=1
  967. IF_NOTOK_MDF(fid=ncid)
  968. end if
  969. if (need_remap) then
  970. call MDF_Get_Var( ncid, varid_rm, src_glb, status, start=(/1,1,1,itr_file/))
  971. IF_NOTOK_MDF(fid=ncid)
  972. src_glb = src_glb / airmass
  973. call Fill3D( global_lli(n), levi, 'n', sp_gbl(:,:,1), rmt(:,:,:,itr), &
  974. lli_restart, levi_restart, src_glb, 'mass-aver', status )
  975. IF_NOTOK_RETURN(status=1)
  976. rmt(:,:,:,itr) = rmt(:,:,:,itr) * run_airmass
  977. else
  978. call MDF_Get_Var( ncid, varid_rm, rmt(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  979. IF_NOTOK_MDF(fid=ncid)
  980. if (istart==32) then
  981. rmt(:,:,:,itr) = rmt(:,:,:,itr) * run_airmass / airmass
  982. endif
  983. endif
  984. #ifdef slopes
  985. ! read slopes
  986. if ((.not. need_remap) .and. (istart==33)) then
  987. if (isRoot) call MDF_Get_Var( ncid, varid_rxm, rmx(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  988. IF_NOTOK_MDF(fid=ncid)
  989. if (isRoot) call MDF_Get_Var( ncid, varid_rym, rmy(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  990. IF_NOTOK_MDF(fid=ncid)
  991. if (isRoot) call MDF_Get_Var( ncid, varid_rzm, rmz(:,:,:,itr), status, start=(/1,1,1,itr_file/))
  992. IF_NOTOK_MDF(fid=ncid)
  993. endif
  994. #endif
  995. else ! short lived tracer:
  996. if ( itr_file <= ntracet_restart ) then
  997. write (gol,'("tracer `",a,"` is not-transported but seems to be transported in restart file")') trim(names(itr)); call goErr
  998. status=1
  999. IF_NOTOK_MDF(fid=ncid)
  1000. end if
  1001. itr_file = itr_file - ntracet_restart
  1002. if (need_remap) then
  1003. call MDF_Get_Var( ncid, varid_rmc, src_glb, status, start=(/1,1,1,itr_file/) )
  1004. IF_NOTOK_MDF(fid=ncid)
  1005. src_glb = src_glb / airmass
  1006. call Fill3D( global_lli(n), levi, 'n', sp_gbl(:,:,1), rms(:,:,:,itr), &
  1007. lli_restart, levi_restart, src_glb, 'mass-aver', status )
  1008. IF_NOTOK_RETURN(status=1)
  1009. rms(:,:,:,itr) = rms(:,:,:,itr) * run_airmass
  1010. else
  1011. call MDF_Get_Var( ncid, varid_rmc, rms(:,:,:,itr), status, start=(/1,1,1,itr_file/) )
  1012. IF_NOTOK_MDF(fid=ncid)
  1013. if (istart==32) then
  1014. rms(:,:,:,itr) = rms(:,:,:,itr) * run_airmass / airmass
  1015. endif
  1016. endif
  1017. end if ! transported or short-lived
  1018. endif ! in the file
  1019. endif ! root
  1020. end do ! tracers
  1021. ! distribute
  1022. call scatter( dgrid(n), mass_dat(n)%rm, rmt, mass_dat(n)%halo, status)
  1023. IF_NOTOK_RETURN(status=1)
  1024. if ( ntrace_chem > 0 ) then
  1025. call scatter( dgrid(n), chem_dat(n)%rm, rms, chem_dat(n)%halo, status)
  1026. IF_NOTOK_RETURN(status=1)
  1027. endif
  1028. #ifdef slopes
  1029. if ((.not. need_remap).and.(istart==33)) then
  1030. call scatter( dgrid(n), mass_dat(n)%rxm, rmx, mass_dat(n)%halo, status)
  1031. IF_NOTOK_RETURN(status=1)
  1032. call scatter( dgrid(n), mass_dat(n)%rym, rmy, mass_dat(n)%halo, status)
  1033. IF_NOTOK_RETURN(status=1)
  1034. call scatter( dgrid(n), mass_dat(n)%rzm, rmz, mass_dat(n)%halo, status)
  1035. IF_NOTOK_RETURN(status=1)
  1036. else
  1037. ! Ensure that slopes are initialized to "unset" values of 0.0. Wouter says that
  1038. ! we could remap levels for rxm et al., but 0.0 will also work. The noise
  1039. ! induced from remapping the rm array is almost certainly bigger than any issues
  1040. ! from having this "default=0.0" slopes information. -ARJ 1 Jan 12
  1041. mass_dat(n)%rxm = 0.0
  1042. mass_dat(n)%rym = 0.0
  1043. mass_dat(n)%rzm = 0.0
  1044. endif
  1045. #endif
  1046. ! free mem for next region
  1047. deallocate( values_names)
  1048. if (need_remap) then
  1049. deallocate(sp_gbl,src_glb)
  1050. if (isRoot) then
  1051. call Done( levi_restart, status )
  1052. IF_NOTOK_RETURN(status=1)
  1053. call Done( lli_restart, status )
  1054. IF_NOTOK_RETURN(status=1)
  1055. endif
  1056. endif
  1057. ENDIF READRM
  1058. ! clean
  1059. deallocate(rmt)
  1060. if ( ntrace_chem > 0 ) deallocate(rms)
  1061. #ifdef slopes
  1062. deallocate(rmx, rmy, rmz)
  1063. #endif
  1064. deallocate( tmp3d )
  1065. deallocate( airmass)
  1066. if (istart==32) deallocate(run_airmass)
  1067. if (isRoot) call MDF_Close( ncid, status )
  1068. IF_NOTOK_RETURN(status=1)
  1069. ENDDO REG
  1070. status = 0
  1071. END SUBROUTINE RESTART_READ
  1072. !EOC
  1073. END MODULE RESTART