prism_putget.F90 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652
  1. !
  2. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  3. !
  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. #define PRISM_ERR call prism_error(status,gol); call goErr
  8. #define IF_PRISM_NOTOK_RETURN(action) if (status/=0) then; PRISM_ERR; TRACEBACK; action; return; end if
  9. !
  10. #include "tm5.inc"
  11. !
  12. !-----------------------------------------------------------------------------
  13. ! TM5 !
  14. !-----------------------------------------------------------------------------
  15. !BOP
  16. !
  17. ! !MODULE: PRISM_PUTGET
  18. !
  19. ! !DESCRIPTION: "Get" and "Put" PRISM
  20. !\\
  21. !\\
  22. ! !INTERFACE:
  23. !
  24. MODULE PRISM_PUTGET
  25. !
  26. ! !USES:
  27. !
  28. USE MOD_OASIS
  29. use GO, only : gol, goPr, goErr
  30. use GO, only : TDate, TIncrDate, IncrDate, Pretty, operator(+)
  31. use TM5_PRISM, only : comp_id, coupled_to_lpj, coupled_to_pis
  32. use TM5_Prism, only : ifs_cpl_freq, lpj_cpl_freq, pis_cpl_freq, SetPrismTime, CplVar, ncplvar
  33. use dims, only : im, jm, lm, lme, ndyn_max
  34. use partools, only : isroot, par_broadcast, myid
  35. use tm5_distgrid, only : dgrid, Get_DistGrid, GATHER, SCATTER, Dist_Arr_Stat
  36. use tracer_data, only : mass_dat
  37. use chem_param, only : xmc, xmco2, ico2
  38. use binas, only : xmair
  39. !
  40. implicit none
  41. private
  42. !
  43. ! !PUBLIC MEMBER FUNCTIONS:
  44. !
  45. public :: TM5_Prism_Puts, TM5_Prism_Gets
  46. !
  47. ! !PRIVATE DATA MEMBERS:
  48. !
  49. character(len=*), parameter :: mname = 'prism_putget'
  50. character(len=256) :: error_message
  51. !
  52. ! !REVISION HISTORY:
  53. ! 14 Dec 2010 - Ph. Le Sager - added vertical regridding in "put".
  54. ! 8 Oct 2013 - Ph. Le Sager - dummy CO2 exchange with LPJ-Guess
  55. !
  56. ! !REMARKS:
  57. !
  58. !EOP
  59. !------------------------------------------------------------------------
  60. CONTAINS
  61. !--------------------------------------------------------------------------
  62. ! TM5 !
  63. !--------------------------------------------------------------------------
  64. !BOP
  65. !
  66. ! !IROUTINE: TM5_PRISM_PUTS
  67. !
  68. ! !DESCRIPTION:
  69. !\\
  70. !\\
  71. ! !INTERFACE:
  72. !
  73. SUBROUTINE TM5_PRISM_PUTS( t, status )
  74. !
  75. ! !USES:
  76. !
  77. use meteodata, only : m_dat
  78. #ifdef with_ecearth_optics
  79. use ecearth_optics, only : optics_dat
  80. use meteodata, only : gph_dat
  81. #endif
  82. ! for vertical regridding
  83. use Grid, only : FillLevelsParents
  84. use meteodata, only : levi, sp_dat
  85. #ifdef with_m7
  86. use chem_param, only : isoanus, isoaais, isoaacs, isoacos, isoaaii
  87. #endif
  88. !
  89. ! !INPUT PARAMETERS:
  90. !
  91. type(TDate), intent(in) :: t
  92. !
  93. ! !OUTPUT PARAMETERS:
  94. !
  95. integer, intent(out) :: status
  96. !
  97. ! !REVISION HISTORY:
  98. ! 14 Dec 2010 - P. Le Sager - added vertical regridding of fields that
  99. ! are passed back to IFS
  100. ! 18 Sep 2013 - P. Le Sager - TM5-MP
  101. !EOP
  102. !------------------------------------------------------------------------
  103. !BOC
  104. character(len=*), parameter :: rname = mname//'/TM5_Prism_Puts'
  105. integer :: ivar, ilev, prism_t, region, imr, jmr, lmr
  106. integer :: i1, i2, j1, j2, k
  107. integer :: soa_itr
  108. real, allocatable, target :: mmr_tm5(:,:,:) ! local field on TM5 levels
  109. real, allocatable, target :: mmr_ifs(:,:,:) ! local field on IFS levels
  110. real, pointer :: field3D(:,:,:) ! field being send
  111. ! real(ip_realwp_p), allocatable :: field3D_ip(:,:,:) ! global data to send, with oasis precision
  112. logical :: subset ! TM5 levels = subset of IFS levels ?
  113. character(len=10) :: interp
  114. real, dimension(:,:), pointer :: sp
  115. real, dimension(:,:,:), pointer :: mass, masssoa, airm, gph
  116. type(TDate) :: lag_date
  117. type(TIncrDate) :: deltat
  118. real :: fscaleCO2
  119. fscaleCO2 = xmair / xmco2
  120. !-----------------------------------------------------------------
  121. ! Check if fields are effectively sent or not, to avoid extra work
  122. ! !! THIS ASSUME THAT COUPLED FREQUENCY FOR TM5-LPJG IS MULTIPLE OF THE ONE FOR TM5-IFS
  123. !-----------------------------------------------------------------
  124. ! by design, the lag defined in the namcouple is exactly one dynamic time step, and positive
  125. deltat = IncrDate(sec=ndyn_max)
  126. lag_date=t+deltat
  127. call SetPrismTime( prism_t, lag_date, status )
  128. IF_NOTOK_RETURN(status=1)
  129. ! avoid unneeded work (only if coupled to IFS only)
  130. if ((.not. coupled_to_lpj).and.(.not. coupled_to_pis)) then
  131. if (modulo( prism_t, ifs_cpl_freq*3600) /= 0) return
  132. end if
  133. ! Convert from TM5 time structure to OASIS time structure
  134. call SetPrismTime( prism_t, t, status )
  135. IF_ERROR_RETURN(status=1)
  136. !DBG if ( isroot ) then
  137. !DBG write (gol,'(" prism time : ",i6)') prism_t; call goPr
  138. !DBG end if
  139. ! Region dims
  140. region = 1
  141. imr = im(region) ; jmr = jm(region) ; lmr = lm(region)
  142. subset = lmr /= lme ! check that if nb of levels exchanges ifs_nlev(=lme) read in tm5_prism is same as nb of TM5 levels
  143. #ifdef with_ecearth_optics
  144. gph => gph_dat(region)%data
  145. #endif
  146. !------------------------------------------------------------------------
  147. ! COUPLED VARIABLES
  148. !------------------------------------------------------------------------
  149. ! storage (could go into init)
  150. !------------------
  151. call Get_DistGrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  152. allocate( mmr_tm5(i1:i2, j1:j2, lmr) ) ! local tracer on source level (TM)
  153. if (subset) allocate( mmr_ifs(i1:i2, j1:j2, lme) ) ! local tracer on target level (IFS)
  154. #ifdef parallel_cplng
  155. field3D => mmr_tm5
  156. #else
  157. if (isroot) then
  158. allocate( field3D(imr,jmr,lme) ) ! global gathered field
  159. else
  160. allocate( field3D(1,1,1) )
  161. endif
  162. #endif
  163. ! Send to oasis3
  164. !----------------------------------------------------
  165. VAR: do ivar = 1, ncplvar
  166. if ( CplVar(ivar)%intent /= 'out' ) cycle
  167. select case ( CplVar(ivar)%name )
  168. #ifdef with_ecearth_optics
  169. case ( 'AOD_01', 'AOD_02', 'AOD_03', 'AOD_04', 'AOD_05', 'AOD_06', 'AOD_07', &
  170. 'AOD_08', 'AOD_09', 'AOD_10', 'AOD_11', 'AOD_12', 'AOD_13', 'AOD_14', &
  171. 'SSA_01', 'SSA_02', 'SSA_03', 'SSA_04', 'SSA_05', 'SSA_06', 'SSA_07', &
  172. 'SSA_08', 'SSA_09', 'SSA_10', 'SSA_11', 'SSA_12', 'SSA_13', 'SSA_14', &
  173. 'ASF_01', 'ASF_02', 'ASF_03', 'ASF_04', 'ASF_05', 'ASF_06', 'ASF_07', &
  174. 'ASF_08', 'ASF_09', 'ASF_10', 'ASF_11', 'ASF_12', 'ASF_13', 'ASF_14' )
  175. select case ( CplVar(ivar)%name )
  176. case ( 'AOD_01' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,1)
  177. case ( 'SSA_01' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,1)
  178. case ( 'ASF_01' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,1)
  179. case ( 'AOD_02' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,2)
  180. case ( 'SSA_02' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,2)
  181. case ( 'ASF_02' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,2)
  182. case ( 'AOD_03' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,3)
  183. case ( 'SSA_03' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,3)
  184. case ( 'ASF_03' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,3)
  185. case ( 'AOD_04' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,4)
  186. case ( 'SSA_04' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,4)
  187. case ( 'ASF_04' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,4)
  188. case ( 'AOD_05' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,5)
  189. case ( 'SSA_05' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,5)
  190. case ( 'ASF_05' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,5)
  191. case ( 'AOD_06' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,6)
  192. case ( 'SSA_06' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,6)
  193. case ( 'ASF_06' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,6)
  194. case ( 'AOD_07' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,7)
  195. case ( 'SSA_07' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,7)
  196. case ( 'ASF_07' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,7)
  197. case ( 'AOD_08' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,8)
  198. case ( 'SSA_08' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,8)
  199. case ( 'ASF_08' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,8)
  200. case ( 'AOD_09' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,9)
  201. case ( 'SSA_09' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,9)
  202. case ( 'ASF_09' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,9)
  203. case ( 'AOD_10' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,10)
  204. case ( 'SSA_10' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,10)
  205. case ( 'ASF_10' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,10)
  206. case ( 'AOD_11' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,11)
  207. case ( 'SSA_11' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,11)
  208. case ( 'ASF_11' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,11)
  209. case ( 'AOD_12' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,12)
  210. case ( 'SSA_12' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,12)
  211. case ( 'ASF_12' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,12)
  212. case ( 'AOD_13' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,13)
  213. case ( 'SSA_13' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,13)
  214. case ( 'ASF_13' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,13)
  215. case ( 'AOD_14' ) ; mmr_tm5 = optics_dat(region)%Ext(:,:,:,14)
  216. case ( 'SSA_14' ) ; mmr_tm5 = optics_dat(region)%a (:,:,:,14)
  217. case ( 'ASF_14' ) ; mmr_tm5 = optics_dat(region)%g (:,:,:,14)
  218. case default
  219. write (gol,'("unsupported optical variable name: ",a)') trim(CplVar(ivar)%name); call goErr
  220. TRACEBACK; status=1; return
  221. end select
  222. if (subset) then ! vertical remapping
  223. ! From email exchange b/w KNMI experts:
  224. ! How to regrid the optical properties (extinction coefficient
  225. ! (Ext) [1/m], single-scattering albedo (SSA) [-], and asymmetry
  226. ! factor (AF) [-])?
  227. !
  228. ! For a temporal averaging, let us take the monthly mean AOD, and
  229. ! calculate the monthly means of the other optical parameters SSA
  230. ! and ASYM in the proper way, by weighting with the amount of
  231. ! scattering and absorption. The recipe is:
  232. !
  233. ! monthly mean AOD = <AOD> = 1/N Sum (AOD(i)) (i=1,...,N; i = day, N = number of days in the month)
  234. ! monthly mean SSA = <SSA> = 1/N Sum (SSA(i) x AOD(i)) / <AOD>
  235. ! monthly mean ASYM = <ASYM> = 1/N Sum (ASYM(i) x SSA(i) x AOD(i)) / (<SSA> x <AOD>)
  236. !
  237. ! BUT HERE WE REGRID VERTICALLY, NOT IN TIME:
  238. ! The SSA and AF values in the IFS sublayers should be the same
  239. ! as in the corresponding merged TM5 layer. Basically they are
  240. ! intensive variables w/r/t pressure levels.
  241. !
  242. ! For Ext we can do as follows:
  243. !
  244. ! - convert to AOD = Ext x layer height
  245. ! - assume that AOD is distributed according to the mass
  246. ! distribution in the sublayers. This is consistent with the
  247. ! assumption of constant (?) aerosol mixing ratios.
  248. ! - convert back to Ext - could be done in IFS?
  249. !
  250. select case ( CplVar(ivar)%name(1:3) )
  251. case('AOD')
  252. ! convert to AOD
  253. !PLS if new_units: case('AOD', 'SSA', 'ASF' )
  254. !PLS if new_units: ! all optical properties are now proportional to extinction (1/m)
  255. do ilev=1,lmr
  256. mmr_tm5(:,:,ilev) = mmr_tm5(:,:,ilev) * (gph(i1:i2, j1:j2,ilev+1) - gph(i1:i2, j1:j2,ilev))
  257. end do
  258. interp='sum'
  259. case('SSA','ASF') ! remove this line if new_units
  260. interp='mass-aver' ! remove this line if new_units
  261. case default
  262. TRACEBACK; status=1; return
  263. end select
  264. sp => sp_dat(region)%data(i1:i2,j1:j2,1)
  265. call FillLevelsParents(levi,'n',mmr_tm5, trim(interp), mmr_ifs, status, sp)
  266. IF_NOTOK_RETURN(status=1)
  267. ! ! convert AOD back to Ext (dh_ifs unknown)
  268. ! if (trim(interp) == 'sum') then
  269. ! mmr_ifs = mmr_ifs * dh_ifs
  270. ! endif
  271. #ifdef parallel_cplng
  272. field3D => mmr_ifs
  273. #else
  274. call GATHER( dgrid(region), mmr_ifs, field3D, 0, status )
  275. IF_NOTOK_RETURN(status=1)
  276. #endif
  277. else
  278. !PLS If using new_units in ./ecearth_optics.F90 then this should be:
  279. !PLS ! convert extinctions to optical depths
  280. !PLS if (CplVar(ivar)%name(1:3) == 'AOD' .or. &
  281. !PLS CplVar(ivar)%name(1:3) == 'SSA' .or. &
  282. !PLS CplVar(ivar)%name(1:3) == 'ASF') then
  283. !PLS and IFS code needs changes too!
  284. ! convert Ext to AOD
  285. if (CplVar(ivar)%name(1:3) == 'AOD') then
  286. do ilev=1,lmr
  287. mmr_tm5(:,:,ilev) = mmr_tm5(:,:,ilev) * (gph(i1:i2, j1:j2,ilev+1) - gph(i1:i2, j1:j2,ilev))
  288. end do
  289. endif
  290. #ifdef parallel_cplng
  291. field3D => mmr_tm5
  292. #else
  293. call GATHER( dgrid(region), mmr_tm5, field3D, 0, status )
  294. IF_NOTOK_RETURN(status=1)
  295. #endif
  296. endif
  297. #endif
  298. case ('ppmCO2') ! CO2 concentration (ppmv) to LPJ-Guess/PISCES
  299. mass => mass_dat(region)%rm(i1:i2,j1:j2,:,CplVar(ivar)%itr)
  300. airm => m_dat(region)%data(i1:i2,j1:j2,:)
  301. mmr_tm5(:,:,1) = mass(:,:,1)*fscaleCO2*1e6 / airm(:,:,1) ! Convert from mass to ppmv
  302. #ifdef parallel_cplng
  303. field3D => mmr_tm5
  304. #else
  305. call GATHER( dgrid(region), mmr_tm5(:,:,1), field3D(:,:,1), 0, status )
  306. IF_NOTOK_RETURN(status=1)
  307. #endif
  308. case default ! TRACER MASS
  309. mass => mass_dat(region)%rm(i1:i2,j1:j2,:,CplVar(ivar)%itr)
  310. airm => m_dat(region)%data(i1:i2,j1:j2,:)
  311. sp => sp_dat(region)%data(i1:i2,j1:j2,1)
  312. select case ( CplVar(ivar)%name )
  313. #ifdef with_m7
  314. case ('OC1','OC2','OC3','OC4','OC5')
  315. ! Select itracer for SOA compound to be added to OCn field to be send to IFS. This
  316. ! can be done because SOA and POA are assumed to have the same mass density,
  317. ! hygroscopicity, and LW absorption coefficient. The nucleation mode can be ignored,
  318. ! because it is neglected in the cloud activation scheme as well as in the mass-based
  319. ! calculation of LW absorption.
  320. select case (CplVar(ivar)%name)
  321. case ('OC1')
  322. soa_itr=isoanus
  323. case ('OC2')
  324. soa_itr=isoaais
  325. case ('OC3')
  326. soa_itr=isoaacs
  327. case ('OC4')
  328. soa_itr=isoacos
  329. case ('OC5')
  330. soa_itr=isoaaii
  331. end select
  332. masssoa => mass_dat(region)%rm(i1:i2,j1:j2,:,soa_itr)
  333. mmr_tm5 = (mass + masssoa) / airm ! Convert from mass to mass-mixing-ratio
  334. #endif
  335. case default
  336. mmr_tm5 = mass / airm ! Convert from mass to mass-mixing-ratio
  337. end select
  338. if (subset) then ! vertical remapping
  339. call FillLevelsParents(levi,'n',mmr_tm5,'mass-aver', mmr_ifs, status, sp)
  340. IF_NOTOK_RETURN(status=1)
  341. #ifdef parallel_cplng
  342. field3D => mmr_ifs
  343. #else
  344. call GATHER( dgrid(region), mmr_ifs, field3D, 0, status )
  345. IF_NOTOK_RETURN(status=1)
  346. #endif
  347. else
  348. #ifdef parallel_cplng
  349. field3D => mmr_tm5
  350. #else
  351. call GATHER( dgrid(region), mmr_tm5, field3D, 0, status )
  352. IF_NOTOK_RETURN(status=1)
  353. #endif
  354. endif
  355. end select
  356. ! ---- SEND ----
  357. #ifndef parallel_cplng
  358. if ( isroot ) then
  359. #endif
  360. do ilev = 1, CplVar(ivar)%nlev
  361. ! reverse layer order if needed, but skip surface fields
  362. if ( CplVar(ivar)%nlev /= 1 ) then
  363. k=CplVar(ivar)%nlev+1-ilev
  364. if (subset) k=ilev+lme-CplVar(ivar)%nlev
  365. else
  366. k = ilev ! i.e. 1
  367. endif
  368. call oasis_put( CplVar(ivar)%var_id(ilev), prism_t, field3D(:,:,k), status )
  369. !DBG if(isroot .and.(ilev == 1).and.(ivar==53).and.(status/=OASIS_Ok)) then
  370. !DBG write (gol,*)" prism sent : ", prism_t, status; call goPr
  371. !DBG endif
  372. select case ( status )
  373. case ( OASIS_Sent, OASIS_LocTrans, OASIS_ToRest, OASIS_Output, &
  374. OASIS_SentOut, OASIS_ToRestOut, OASIS_WaitGroup, OASIS_Ok )
  375. continue
  376. case default
  377. TRACEBACK
  378. write (error_message,'("from OASIS_PUT : ",i6)') status
  379. call oasis_abort( comp_id, rname, error_message )
  380. end select
  381. end do
  382. #ifndef parallel_cplng
  383. end if
  384. #endif
  385. end do VAR
  386. ! Done
  387. deallocate( mmr_tm5 )
  388. #ifndef parallel_cplng
  389. deallocate( field3D )
  390. #endif
  391. if (subset) deallocate(mmr_ifs)
  392. status = 0
  393. END SUBROUTINE TM5_PRISM_PUTS
  394. !EOC
  395. !--------------------------------------------------------------------------
  396. ! TM5 !
  397. !--------------------------------------------------------------------------
  398. !BOP
  399. !
  400. ! !IROUTINE: TM5_Prism_Gets
  401. !
  402. ! !DESCRIPTION: To receive tracer fields from other models; Meteo from IFS
  403. ! is received in tmm_mf_prism.F90.
  404. !\\
  405. !\\
  406. ! !INTERFACE:
  407. !
  408. SUBROUTINE TM5_PRISM_GETS( t, status )
  409. !
  410. ! !USES:
  411. !
  412. use grid, only : TllGridInfo, AreaOper, done
  413. !
  414. ! !INPUT PARAMETERS:
  415. !
  416. type(TDate), intent(in) :: t
  417. !
  418. ! !OUTPUT PARAMETERS:
  419. !
  420. integer, intent(out) :: status
  421. !
  422. ! !REMARKS:
  423. !
  424. !EOP
  425. !------------------------------------------------------------------------
  426. !BOC
  427. character(len=*), parameter :: rname = mname//'/TM5_Prism_Gets'
  428. integer :: ivar, prism_t, region, imr, jmr, istat
  429. integer :: i1, i2, j1, j2
  430. integer :: ii,jj,kk
  431. type(TllGridInfo) :: lli
  432. real, allocatable :: mmr_tm5(:,:) ! local field
  433. real, allocatable :: tempo_fld(:,:) ! global field to receive fields
  434. real, allocatable :: oasis_fld(:,:) ! global field to accumulate received field
  435. logical :: has_data
  436. ! --- begin -----------------------------------
  437. if ((.not. coupled_to_lpj) .and. (.not. coupled_to_pis)) then
  438. status=0
  439. return
  440. end if
  441. !DBG write (gol,'("get C fluxes from LPJ-Guess and/or PISCES")'); call goPr
  442. !
  443. ! if ( (modulo(t%hour,ifs_cpl_freq)/=0) .or. any((/t%min,t%sec,t%mili/)/=0) ) then
  444. ! write (gol,'(" skip; only every ",i2," hour ...")') ifs_cpl_freq; call goPr
  445. ! status=0; return
  446. ! end if
  447. ! Storage
  448. !----------------------------------------------------
  449. region = 1
  450. call Get_DistGrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2, lli=lli )
  451. imr = im(region) ; jmr = jm(region)
  452. #ifdef parallel_cplng
  453. allocate( oasis_fld(i1:i2, j1:j2) )
  454. allocate( tempo_fld(i1:i2, j1:j2) )
  455. #else
  456. allocate( mmr_tm5(i1:i2, j1:j2) )
  457. if (isroot) then
  458. allocate( oasis_fld(imr,jmr) )
  459. allocate( tempo_fld(imr,jmr) )
  460. else
  461. allocate( oasis_fld(1,1) )
  462. allocate( tempo_fld(1,1) )
  463. endif
  464. #endif
  465. ! Prism time
  466. !----------------------------------------------------
  467. call SetPrismTime( prism_t, t, status )
  468. IF_ERROR_RETURN(status=1)
  469. !DBG if ( isroot ) then
  470. !DBG write (gol,'(" prism time : ",i6)') prism_t; call goPr
  471. !DBG end if
  472. ! Receive variable from LPJ-Guess and/or PISCES
  473. !----------------------------------------------------
  474. has_data = .false.
  475. #ifndef parallel_cplng
  476. IF (isroot) THEN
  477. #endif
  478. oasis_fld = 0.0
  479. VARIA: DO ivar = 1, ncplvar
  480. SELECT CASE ( CplVar(ivar)%cpl_name )
  481. CASE ( 'land_c_flux_nat', 'land_c_flux_ant', 'land_c_npp', 'oce_c_flux' )
  482. TEMPO_FLD = 0.0
  483. CALL OASIS_GET( CplVar(ivar)%var_id(1), prism_t, tempo_fld, status )
  484. SELECT CASE ( status )
  485. CASE (OASIS_Recvd, OASIS_FromRest, OASIS_Input, OASIS_RecvOut, OASIS_FromRestOut )
  486. ! Convert C-fluxes into kg(CO2)/m2
  487. SELECT CASE ( CplVar(ivar)%cpl_name )
  488. CASE ( 'land_c_flux_nat' ) ! LPJG: kg(C)/m2/day
  489. TEMPO_FLD = TEMPO_FLD * xmco2/xmc
  490. CASE ( 'land_c_flux_ant' ) ! LPJG: kg(C)/m2/day
  491. TEMPO_FLD = TEMPO_FLD * xmco2/xmc
  492. CASE ( 'land_c_npp' ) ! LPJG: kg(C)/m2/day
  493. TEMPO_FLD = TEMPO_FLD * xmco2/xmc
  494. CASE ( 'oce_c_flux' ) ! PISCES: molC/m2/s
  495. TEMPO_FLD = TEMPO_FLD * pis_cpl_freq * 3600 * xmco2 / 1000
  496. END SELECT
  497. OASIS_FLD = OASIS_FLD + TEMPO_FLD
  498. HAS_DATA = .TRUE.
  499. CASE ( OASIS_OK )
  500. CONTINUE
  501. CASE DEFAULT
  502. TRACEBACK
  503. WRITE (error_message,'("PRISM_PUTGET: Error in OASIS_GET:",i6)') status
  504. CALL OASIS_ABORT( comp_id, rname, error_message )
  505. END SELECT
  506. END SELECT
  507. END DO VARIA
  508. #ifndef parallel_cplng
  509. END IF
  510. call par_broadcast(has_data, status)
  511. IF_NOTOK_RETURN(status=1)
  512. #endif
  513. ! add data if any
  514. if ( has_data ) then
  515. ! DBG EXAMPLE - this will print the min/max/mean of the oasis_fld array
  516. !call Dist_Arr_Stat(dgrid(region), 'oasisfld', oasis_fld, 0, status)
  517. !IF_NOTOK_RETURN(status=1)
  518. call AreaOper( lli, oasis_fld, '*', 'm2', status) ! into kg(CO2)/box
  519. IF_NOTOK_RETURN(status=1)
  520. #ifndef parallel_cplng
  521. CALL SCATTER( dgrid(region), mmr_tm5, oasis_fld, 0, status )
  522. IF_NOTOK_RETURN(status=1)
  523. mass_dat(region)%rm(i1:i2,j1:j2,1,ico2) = mass_dat(region)%rm(i1:i2,j1:j2,1,ico2) + mmr_tm5
  524. #else
  525. ! spread out flux over 2 lowest levels to avoid negative mass
  526. mass_dat(region)%rm(i1:i2,j1:j2,1,ico2) = mass_dat(region)%rm(i1:i2,j1:j2,1,ico2) + oasis_fld/2.0
  527. mass_dat(region)%rm(i1:i2,j1:j2,2,ico2) = mass_dat(region)%rm(i1:i2,j1:j2,2,ico2) + oasis_fld/2.0
  528. #endif
  529. ! safety check for negative co2 concentrations, see issue #706
  530. if ( minval( mass_dat(region)%rm(i1:i2,j1:j2,1:2,ico2) ) .lt. 0 ) then
  531. do kk=1,2
  532. do jj=j1,j2
  533. do ii=i1,i2
  534. if ( ( mass_dat(region)%rm(ii,jj,kk,ico2) ) .lt. 0 ) then
  535. write (gol,'("negative co2 mass after OASIS_GET prism_t= ",i20," ii= ",i6," jj= ",i6," kk= ",i6,")")') prism_t,ii,jj,kk; call goErr
  536. write (gol,'("mass b4= ",e15.6," mass af= ",e15.6," oasisfld= ",e15.6,")")') (mass_dat(region)%rm(ii,jj,kk,ico2)-oasis_fld(ii,jj)/2.0),(mass_dat(region)%rm(ii,jj,kk,ico2)),oasis_fld(ii,jj)/2.0; call goErr
  537. end if
  538. end do
  539. end do
  540. end do
  541. status = 1
  542. write (gol,'("negative co2 mass after OASIS_GET")'); call goErr
  543. call par_broadcast(status, istat, myid)
  544. IF_NOTOK_RETURN(status=1)
  545. endif
  546. endif
  547. ! done
  548. call done(lli, status)
  549. IF_NOTOK_RETURN(status=1)
  550. DEALLOCATE( oasis_fld, tempo_fld )
  551. if(allocated(mmr_tm5)) deallocate(mmr_tm5)
  552. status = 0
  553. END SUBROUTINE TM5_PRISM_GETS
  554. !EOC
  555. END MODULE PRISM_PUTGET