prism_putget.F90 21 KB

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