emission_ss.F90 23 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610
  1. !
  2. #define TRACEBACK write (gol,'("in ",a," (",a,", line",i5,")")') rname, __FILE__, __LINE__; call goErr
  3. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  4. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  5. !
  6. #include "tm5.inc"
  7. !
  8. !-----------------------------------------------------------------------------
  9. ! TM5 !
  10. !-----------------------------------------------------------------------------
  11. !BOP
  12. !
  13. ! !MODULE: EMISSION_SS
  14. !
  15. ! !DESCRIPTION: Perform Sea-salt emissions needed for M7 version.
  16. ! Sea-salt 2 modes are accumulation and coarse soluble.
  17. !\\
  18. !\\
  19. ! !INTERFACE:
  20. !
  21. MODULE EMISSION_SS
  22. !
  23. ! !USES:
  24. !
  25. USE GO, ONLY : gol, goErr, goPr
  26. USE global_types, ONLY : d3_data, emis_data
  27. IMPLICIT NONE
  28. PRIVATE
  29. !
  30. ! !PUBLIC MEMBER FUNCTIONS:
  31. !
  32. PUBLIC :: calc_emission_ss, emission_ss_init
  33. !
  34. ! !PRIVATE DATA MEMBERS:
  35. !
  36. CHARACTER(len=*), PARAMETER :: mname = 'emission_ss'
  37. !
  38. ! !REVISION HISTORY:
  39. ! ? ??? 2005 - Elisabetta Vignati - changed for coupling with M7
  40. !
  41. ! ? ??? 2006 - EV and MK - changed for introducing the sea salt
  42. ! source function developed from Gong (2003) in two modes
  43. ! ? ??? ???? - AJS - switch from default aerocom emissions to Gong
  44. ! parameterisation if 'seasalt_emis_gong' is defined.
  45. ! 1 Sep 2010 - Achim Strunk - deleted procedures
  46. ! - removed with_seasalt-switch
  47. ! - introduced vertical splitting
  48. ! 25 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  49. ! April 2015 - T. van Noije - modified mode prefactors
  50. !
  51. ! !REMARKS:
  52. !
  53. !EOP
  54. !------------------------------------------------------------------------
  55. CONTAINS
  56. !--------------------------------------------------------------------------
  57. ! TM5 !
  58. !--------------------------------------------------------------------------
  59. !BOP
  60. !
  61. ! !IROUTINE: EMISSION_SS_INIT
  62. !
  63. ! !DESCRIPTION:
  64. !\\
  65. !\\
  66. ! !INTERFACE:
  67. !
  68. SUBROUTINE emission_ss_init( status )
  69. !
  70. ! !USES:
  71. !
  72. USE dims, ONLY : iglbsfc
  73. USE meteo, ONLY : Set
  74. USE meteodata, ONLY : wspd_dat, ci_dat, lsmask_dat
  75. USE meteodata, ONLY : sst_dat
  76. !
  77. ! !OUTPUT PARAMETERS:
  78. !
  79. INTEGER, INTENT(out) :: status
  80. !
  81. ! !REVISION HISTORY:
  82. ! 1 Jul 2013 - Ph Le Sager - v0
  83. !
  84. ! !REMARKS:
  85. !
  86. !EOP
  87. !------------------------------------------------------------------------
  88. !BOC
  89. CALL Set( wspd_dat(iglbsfc), status, used=.TRUE. )
  90. CALL Set( ci_dat(iglbsfc), status, used=.TRUE. )
  91. CALL Set( lsmask_dat(iglbsfc), status, used=.TRUE. )
  92. CALL Set( sst_dat(iglbsfc), status, used=.TRUE. )
  93. END SUBROUTINE emission_ss_init
  94. !EOC
  95. !--------------------------------------------------------------------------
  96. ! TM5 !
  97. !--------------------------------------------------------------------------
  98. !BOP
  99. !
  100. ! !IROUTINE: CALC_EMISSION_SS
  101. !
  102. ! !DESCRIPTION: Calculate emitted numbers and mass as function of the ten-meter
  103. ! wind speed as described in Vignati et al. (2010) and Gong
  104. ! (2003). Sea salt is emitted in both the soluble accumulation
  105. ! mode and the soluble coarse mode.
  106. !
  107. ! Ref: Vignati et al. (2010) : Global scale emission and
  108. ! distribution of sea-spray aerosol: Sea-salt and organic
  109. ! enrichment, Atmos. Environ., 44, 670-677,
  110. ! doi:10.1016/j.atmosenv.2009.11.013
  111. !
  112. ! Gong (2003) : A parameterization of sea-salt aerosol source
  113. ! function for sub- and super-micron particles, Global
  114. ! Biogeochem. Cy., 17, 1097, doi:10.1029/2003GB002079
  115. !\\
  116. !\\
  117. ! !INTERFACE:
  118. !
  119. SUBROUTINE calc_emission_ss( status )
  120. !
  121. ! !USES:
  122. !
  123. USE Binas, ONLY : pi, T0
  124. USE dims, ONLY : okdebug, itaur, nsrce, tref
  125. USE datetime, ONLY : tau2date
  126. USE dims, ONLY : nlon360, nlat180, idate, itau, nregions, staggered, dxy11
  127. USE dims, ONLY : sec_month, iglbsfc, im, jm, lm
  128. USE chem_param, ONLY : sigma_lognormal, ss_density, nmodes, mode_acs, mode_cos
  129. USE chem_param, ONLY : iacs_n, icos_n, issacs, isscos
  130. USE toolbox, ONLY : coarsen_emission
  131. USE emission_data, ONLY : emis_mass, emis_number, emis_temp, msg_emis
  132. USE emission_data, ONLY : radius_ssa, radius_ssc
  133. USE partools, ONLY : isRoot
  134. USE tm5_distgrid, ONLY : dgrid, get_distgrid, scatter, gather
  135. USE meteodata, ONLY : wspd_dat, ci_dat, lsmask_dat
  136. USE meteodata, ONLY : sst_dat
  137. USE emission_data, ONLY : vd_class_name_len, emission_vdist_by_sector
  138. !
  139. ! !OUTPUT PARAMETERS:
  140. !
  141. INTEGER, INTENT(out) :: status
  142. !
  143. ! !REVISION HISTORY:
  144. ! 1 Oct 2010 - Achim Strunk - introduced vertical splitting
  145. ! 25 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  146. !
  147. ! !REMARKS:
  148. ! (1) done on 1x1 grid.
  149. ! (2) this routine basicaly is the "declare" routine for sea-salt.
  150. !
  151. !EOP
  152. !------------------------------------------------------------------------
  153. !BOC
  154. CHARACTER(len=*), PARAMETER :: rname = mname//'/calc_emission_ss'
  155. REAL, DIMENSION(:,:), ALLOCATABLE :: number, mass, glb_arr
  156. INTEGER :: region
  157. REAL :: xwind, rg1, rg2, dens, xsea
  158. !>>> TvN
  159. REAL :: norm
  160. REAL, DIMENSION(:,:), ALLOCATABLE :: emis_fac
  161. REAL :: area_frac
  162. ! ratio of prefactors in Eqs. (2) and (1) of Salisbury et al. (GRL, 2014):
  163. !REAL, PARAMETER :: W10_fac = 4.60e-3/3.84e-4
  164. !REAL, PARAMETER :: W10_exp = 2.26
  165. REAL, DIMENSION(:,:), ALLOCATABLE :: temperature
  166. REAL :: tt, t_scale
  167. !<<< TvN
  168. INTEGER, PARAMETER :: add_field = 0
  169. INTEGER, PARAMETER :: amonth = 2
  170. INTEGER :: i, j, i1, i2, j1, j2, i01, i02, j01, j02
  171. TYPE(d3_data) :: emis3d
  172. TYPE(emis_data), DIMENSION(nregions) :: emis_glb
  173. CHARACTER(len=vd_class_name_len) :: splittype
  174. ! --- begin -----------------------------------------
  175. ! vertical splitting is that class
  176. splittype = 'surface'
  177. !===================
  178. ! Work arrays
  179. !===================
  180. CALL GET_DISTGRID( dgrid(iglbsfc), I_STRT=i01, I_STOP=i02, J_STRT=j01, J_STOP=j02 )
  181. ALLOCATE(number(i01:i02,j01:j02))
  182. ALLOCATE(mass (i01:i02,j01:j02))
  183. !>>> TvN
  184. ALLOCATE(emis_fac(i01:i02,j01:j02)) ; emis_fac=0.
  185. ALLOCATE(temperature(i01:i02,j01:j02))
  186. temperature = sst_dat(iglbsfc)%data(:,:,1)
  187. ! convert to celsius
  188. temperature = temperature - T0
  189. !<<< TvN
  190. ! no zoom region if region #1 is at 1x1
  191. IF ( (iglbsfc /= 1) .OR. okdebug) THEN
  192. IF (isRoot) THEN
  193. ALLOCATE(glb_arr(nlon360,nlat180))
  194. DO region = 1, nregions
  195. ALLOCATE(emis_glb(region)%surf(im(region), jm(region)))
  196. END DO
  197. ELSE
  198. ALLOCATE(glb_arr(1,1))
  199. DO region = 1, nregions
  200. ALLOCATE(emis_glb(region)%surf(1,1))
  201. END DO
  202. END IF
  203. END IF
  204. !>>> TvN
  205. ! The parameterization of Gong (2003)
  206. ! gives the particle number flux as a function
  207. ! of the radius and the 10-m wind speed u_10:
  208. ! dF/dr80 = f(u_10) x g(r80).
  209. ! The dependence on wind speed is given by
  210. ! the power law f(u_10) = u_10^3.41,
  211. ! and does not affect the size distribution.
  212. ! r80 is the radius at 80% humidity,
  213. ! which is almost exactly 2.0 times the dry radius
  214. ! (Lewis and Schwartz, Sea Salt Aerosol Production, 2004).
  215. ! Note also that in Eq. (2) of Gong
  216. ! B is defined in terms of log(10,r) not ln(r).
  217. !
  218. ! The number mean radii defined in chem_param.F90,
  219. ! i.e. 0.09 um for the accumulation mode,
  220. ! and 0.794 for the coarse mode,
  221. ! were obtained by fitting an accumulation
  222. ! and coarse mode to the dF/dln(r),
  223. ! with r the dry particle radius
  224. ! (see presentation E. Vignati, 21 December 2005).
  225. ! It was verified that the size distribution of Gong
  226. ! can be reasonable well described using these mode radii.
  227. !
  228. ! It is not totally clear how the corresponding prefactors
  229. ! for the two modes were obtained.
  230. ! One way to calculate these prefactors,
  231. ! is to define two size ranges
  232. ! and require that the numbers of particles emitted
  233. ! in both ranges correspond to Gong
  234. ! This procedure is similar to that used in emission_dust.F90,
  235. ! but for particle number instead of mass.
  236. ! Using ranges r1 and r2 for the dry particle radii
  237. ! with r1 from 0.05 to 0.5 um and r2 from 0.5 to 5 um,
  238. ! the resulting prefactors are 9.62013e-3 and 9.05658e-4
  239. ! for the accumulation and coarse mode, respectively.
  240. ! These numbers are very close to the values
  241. ! of 9e-3 and 9e-4, obtained by E. Vignati.
  242. ! They are also insensitive to the
  243. ! value used for the upper bound of r2.
  244. ! (Using a value of 10 instead of 5 um,
  245. ! the prefactors would be 9.62536e-3 and 8.91184e-3.)
  246. !
  247. DO j=j01,j02
  248. norm = 1.e4 * dxy11(j) * sec_month
  249. DO i=i01,i02
  250. ! sea fraction
  251. xsea=1.-lsmask_dat(iglbsfc)%data(i,j,1)/100.
  252. ! sea salt is emitted only over sea without ice cover
  253. area_frac = xsea * (1.-ci_dat(iglbsfc)%data(i,j,1))
  254. if (area_frac .lt. 1.e-10) cycle
  255. emis_fac(i,j) = norm * area_frac
  256. ! Wind speed dependence following W10 parameterization
  257. ! of Salisbury et al. (JGR, 2013; GRL, 2014),
  258. ! which replaces the dependence of Gong.
  259. ! Salisbury et al. suggest that "at this stage ...
  260. ! use of W10 is preferable to W37".
  261. ! In TM5, W10 gives better results than W37,
  262. ! which leads to high AOD values compared to MODIS C6,
  263. ! at mid- and high latitudes.
  264. ! The same is true for the wind dependence
  265. ! proposed by Albert et al. (ACPD, 2015).
  266. xwind = wspd_dat(iglbsfc)%data(i,j,1)
  267. ! Revert to original wind speed dependence of Monahan et al. (1986)
  268. ! used by Gong (2003):
  269. !emis_fac(i,j) = emis_fac(i,j) * W10_fac * xwind**W10_exp
  270. emis_fac(i,j) = emis_fac(i,j) * xwind**3.41
  271. ENDDO
  272. ENDDO
  273. !<<< TvN
  274. !===================
  275. ! Accumulation mode
  276. !===================
  277. ! emitted numbers
  278. ! ---------------
  279. DO j=j01,j02
  280. DO i=i01,i02
  281. !>>> TvN
  282. ! sea fraction
  283. !xsea=1.-lsmask_dat(iglbsfc)%data(i,j,1)/100.
  284. !xwind=SQRT(u10m_dat(iglbsfc)%data(i,j,1)**2+v10m_dat(iglbsfc)%data (i,j,1)**2)
  285. !number(i,j) = 0.009*xwind**(3.4269)*1e4*dxy11(j)*xsea*(1.-ci_dat(iglbsfc)%data(i,j,1))*sec_month ! #/gridbox/month
  286. number(i,j) = emis_fac(i,j)*9.62013e-3 ! #/gridbox/month
  287. !Include temperature dependence following Salter et al. (ACP, 2015).
  288. !We assume that the emissions in our accumulation mode
  289. !follow the temperature dependence of the smallest mode described by Santer et al.
  290. !In reality our mode is in between the smallest and middle mode of Santer et al.,
  291. !so the temperature dependence might also be.
  292. !The corresponding third-order polynomial for the smallest mode
  293. !decreases with temperature between -1 and about 15 degC,
  294. !remains almost constant between 15 degC and 25 degC,
  295. !and shows a very small decrease between 25 degC and 30 degC.
  296. !The latter dependence seems an artefact of the fitting procedure,
  297. !and is neglected here.
  298. !As for the coarse mode,
  299. !we rescale the polynomial of Salter et al.
  300. !so that it goes through 1 at some reference temperature,
  301. !which can currently be set to either 15 or 20 degC.
  302. !Because of the small difference between 15 and 20 degC,
  303. !it really doesn't matter which reference value
  304. !is chosen for the accumulation mode.
  305. !This temperature dependency is used to produce
  306. !the pre-industrial aerosol climatology v4.0,
  307. !which is used in all CMIP6 versions of EC-Earth3 that use MACv2-SP,
  308. !and will therefore also be used in EC-Earth3-AerChem.
  309. tt=max(-1.,temperature(i,j))
  310. !For a reference temperature of 20 degC, use:
  311. !if (tt .lt. 20.) then
  312. ! t_scale = -8.88108055e-5*tt*tt*tt + 5.64728654e-3*tt*tt &
  313. ! -0.118363619*tt + 1.81884421
  314. ! number(i,j) = number(i,j) * t_scale
  315. !endif
  316. !For a reference temperature of 15 degC, use:
  317. if (tt .lt. 15.) then
  318. t_scale = -8.75593266e-5*tt*tt*tt + 5.56770771e-3*tt*tt &
  319. -0.116695696*tt + 1.79321393
  320. number(i,j) = number(i,j) * t_scale
  321. endif
  322. !Using the above relation we still underestimate CDNC by a factor 2 to 4
  323. !over the Soutern Ocean.
  324. !The enhancement factor from Salter et al. (2015) is only 1.9 at -1 degC,
  325. !while earlier studies measured a stronger enhancement
  326. !in the number of submicron particles as SST decreases
  327. !from about 9 to -1 (or -1.3) degC
  328. !(Salter et al., JGR, 2014, Bowyer et al., JGR, 1990;
  329. !Zabori et al. (ACP, 2012a; 2012b)
  330. !In these studies the enhancement at about -1 degC varies
  331. !from ~3 (Salter et al.; Bowyer et al., 1990),
  332. !~7 (Zabori et al., 2012b), to ~10 (Zabori et al., 2012a).
  333. !As the studies don't provide an expression for the enhancement factor
  334. !as a function of temperature,
  335. !I have approximated it by a quadratic function
  336. !which reaches 1 at t0=9 degC with a zero slope:
  337. ![(F-1)/(t0-tmin)^2]*(t-t0)^2 + 1
  338. !where F is the enhancement factor at tmin=-1 degC.
  339. !Using t0=9 and tmin=-1, the quadratic coefficient
  340. !becomes (F-1)/100.
  341. !I have implemented this expression for F=3, F=5, and F=7.
  342. !With F=5 or F=7, the zonal mean AOD over the Southern Ocean
  343. !is overestimated at high latitudes
  344. !compared to observational estimates from AATSR.
  345. !For v5.0 of the pre-industrial aerosol climatology
  346. !we have therefore used F=3.
  347. !tt=max(-1.,temperature(i,j))
  348. !if (tt .lt. 9.) then
  349. !F=3:
  350. !t_scale = 0.02*(tt-9.)*(tt-9.) + 1.0
  351. !F=5:
  352. !t_scale = 0.04*(tt-9.)*(tt-9.) + 1.0
  353. !F=7:
  354. !t_scale = 0.06*(tt-9.)*(tt-9.) + 1.0
  355. !number(i,j) = number(i,j) * t_scale
  356. !endif
  357. !<<< TvN
  358. END DO
  359. END DO
  360. CALL fill_target_array( number, 'number acc', 'ss_number acc' ) ! fill emis_temp(region)
  361. IF_NOTOK_RETURN(status=1)
  362. ! vertically distribute according to sector
  363. DO region = 1, nregions
  364. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  365. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  366. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  367. emis_number(region,mode_acs)%d4(:,:,:,4) = emis3d%d3 !#/grid/month
  368. DEALLOCATE(emis3d%d3)
  369. END DO
  370. ! emitted mass
  371. ! ------------
  372. dens = ss_density*0.001 ! in g/cm3
  373. rg1 = radius_ssa*1e6 ! in micron
  374. mass(:,:) = number(:,:)*pi*4./3. *(rg1*1e-4)**3 * EXP(9./2.*alog(sigma_lognormal(3))**2)*dens*1e-3 ! kg/gridbox/month
  375. CALL fill_target_array( mass, 'mass acc', 'ss_mass acc' ) ! fill emis_temp(region)
  376. IF_NOTOK_RETURN(status=1)
  377. ! vertically distribute according to sector
  378. DO region = 1, nregions
  379. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  380. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  381. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  382. emis_mass (region,mode_acs)%d4(:,:,:,4) = emis3d%d3 !kg/grid/month
  383. DEALLOCATE(emis3d%d3)
  384. END DO
  385. !===================
  386. ! Coarse mode
  387. !===================
  388. ! emitted numbers
  389. ! ---------------
  390. DO j=j01,j02
  391. DO i=i01,i02
  392. ! >>> TvN
  393. ! sea fraction
  394. !xsea=1.-lsmask_dat(iglbsfc)%data(i,j,1)/100.
  395. !xwind=SQRT(u10m_dat(iglbsfc)%data(i,j,1)**2+v10m_dat(iglbsfc)%data (i,j,1)**2)
  396. !number(i,j) = 0.0009*xwind**(3.4195)*dxy11(j)*1e4*xsea*(1.-ci_dat(iglbsfc)%data(i,j,1))*sec_month !#/cm2/s-->#/gridbox/month
  397. number(i,j) = emis_fac(i,j)*9.05658e-4 !#/cm2/s-->#/gridbox/month
  398. !Include temperature dependence following Salter et al. (ACP, 2015).
  399. !We assume that the emissions in our coarse mode
  400. !follow the temperature dependence of the largest mode described by Santer et al.
  401. !Indeed, the mode radii of these modes are close (0.794 vs. 0.75 micron).
  402. !The corresponding second-order polynomial (valid between -1 and 30 degC)
  403. !shows almost linear increase with temperature.
  404. !The higher order terms are therefore irrelevant.
  405. !We rescale the resulting linear function
  406. !so that it goes through 1 at some reference temperature,
  407. !which can currently be set to either 15 or 20 degC.
  408. !The resulting expression is then used as
  409. !an additional factor multiplying the expression from Gong.
  410. !Gong uses the relation between whitecap coverage and wind speed
  411. !from Monahan and Muircheartaigh (1980).
  412. !It was derived from data in the Atlantic and Pacific Ocean,
  413. !where most measurements were made at sea water temperatures
  414. !between 20 and 30 degC: 46 out of 54 points in the Atlantic data set
  415. !and all 36 points in the Pacific data set.
  416. !Based on this, one would argue that a reference temperature
  417. !in the range 20-25 degC would be most reasonable.
  418. !As the emissions in the coarse mode increase with temperature,
  419. !the smaller the reference temperature the higher the emissions.
  420. !For example, the emissions will increase by 18.76%
  421. !when the reference temperature is reduced from 20 to 15 degC.
  422. !The resulting temperature dependence is somewhat weaker than
  423. !the quasi-linear dependence obtained by Ovadnevaite et al.
  424. ![ACP, 2014; see also Gantt et al., GMD, 2015, Eq. (2)].
  425. !The function proposed by Jaegle et al. (JGR, 2011)
  426. !oscillates around or mostly below
  427. !(reference temperature at 20 degC or 15 degC, resp.)
  428. !our linear function up to about 20 degC,
  429. !but increases more strongly at higher temperatures.
  430. !For CMIP6 a reference temperature of 15 degC is used.
  431. tt=min(30.,max(-1.0,temperature(i,j)))
  432. !For a reference temperature of 20 degC, use:
  433. !t_scale = 0.03159982*tt + 0.36800362
  434. !For a reference temperature to 15 degC, use:
  435. t_scale = 0.03752944*tt + 0.43705846
  436. number(i,j) = number(i,j) * t_scale
  437. ! <<< TvN
  438. ENDDO
  439. ENDDO
  440. CALL fill_target_array( number, 'number coa', 'ss_number coa' ) ! fill emis_temp(region)
  441. IF_NOTOK_RETURN(status=1)
  442. ! vertically distribute according to sector
  443. DO region = 1, nregions
  444. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  445. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  446. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  447. emis_number(region,mode_cos)%d4(:,:,:,4) = emis3d%d3 !#/grid/month
  448. DEALLOCATE(emis3d%d3)
  449. END DO
  450. ! emitted mass
  451. ! ------------
  452. dens = ss_density*0.001 ! im g/cm3
  453. rg2 = radius_ssc*1e6 ! in micron
  454. mass(:,:) = number(:,:)*pi*4./3. *(rg2*1e-4)**3 * EXP(9./2.*alog(sigma_lognormal(4))**2)*dens*1e-3 ! kg/gridbox/month
  455. CALL fill_target_array( mass, 'mass coa', 'ss_mass coa' ) ! fill emis_temp(region)
  456. IF_NOTOK_RETURN(status=1)
  457. ! vertically distribute according to sector
  458. DO region = 1, nregions
  459. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  460. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  461. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  462. emis_mass(region,mode_cos)%d4(:,:,:,4) = emis3d%d3 !kg/gridbox/month
  463. DEALLOCATE(emis3d%d3)
  464. END DO
  465. !===================
  466. ! Done
  467. !===================
  468. DEALLOCATE(number, mass)
  469. DEALLOCATE(emis_fac)
  470. DEALLOCATE(temperature)
  471. DO region = 1, nregions
  472. IF (ASSOCIATED(emis_glb(region)%surf)) DEALLOCATE(emis_glb(region)%surf)
  473. END DO
  474. IF (ALLOCATED(glb_arr)) DEALLOCATE(glb_arr)
  475. CONTAINS
  476. !--------------------------------------------------------------------------
  477. ! TM5 !
  478. !--------------------------------------------------------------------------
  479. !BOP
  480. !
  481. ! !IROUTINE: FILL_TARGET_ARRAY
  482. !
  483. ! !DESCRIPTION: Convenient internal program to fill EMIS_TEMP - same as routine in emission_dust.F90
  484. !\\
  485. !\\
  486. ! !INTERFACE:
  487. !
  488. SUBROUTINE fill_target_array( arr_in, str1, str2 )
  489. !
  490. ! !INPUT PARAMETERS:
  491. !
  492. REAL, INTENT(in) :: arr_in(i01:,j01:)
  493. CHARACTER(len=*), INTENT(in) :: str1, str2
  494. !
  495. ! !REVISION HISTORY:
  496. ! 25 Jun 2012 - P. Le Sager - v0
  497. !
  498. !EOP
  499. !------------------------------------------------------------------------
  500. !BOC
  501. ! Test for optimization: if region #1 is at 1x1, assume no zoom region
  502. ! and skip call to coarsen and mpi comm
  503. IF (iglbsfc == 1) THEN
  504. emis_temp(1)%surf = arr_in
  505. IF (okdebug) THEN
  506. CALL gather(dgrid(iglbsfc), arr_in, glb_arr, 0, status)
  507. IF_NOTOK_RETURN(status=1)
  508. !FIXME - really needed? too much messaging
  509. ! IF (isRoot) THEN
  510. ! CALL msg_emis( amonth, ' ', TRIM(str1), 'SS', 1., SUM(glb_arr) )
  511. ! END IF
  512. END IF
  513. ELSE
  514. DO region = 1, nregions
  515. emis_temp(region)%surf = 0.0
  516. END DO
  517. ! gather on 1x1, coarsen to other grid on root, then scatter
  518. !-----------------------------------------------------------
  519. ! FIXME-PERF : Need a coarsen that works independtly on each proc, regardless of zooming... :(
  520. CALL gather(dgrid(iglbsfc), arr_in, glb_arr, 0, status)
  521. IF_NOTOK_RETURN(status=1)
  522. IF (isRoot) THEN
  523. !FIXME - really needed? too much messaging
  524. ! CALL msg_emis( amonth, ' ', TRIM(str1), 'SS', 1., SUM(glb_arr) )
  525. CALL coarsen_emission(TRIM(str2), nlon360, nlat180, glb_arr, emis_glb, add_field, status)
  526. IF_NOTOK_RETURN(status=1)
  527. END IF
  528. DO region = 1, nregions
  529. CALL scatter(dgrid(region), emis_temp(region)%surf, emis_glb(region)%surf, 0, status)
  530. IF_NOTOK_RETURN(status=1)
  531. ENDDO
  532. ENDIF
  533. status=0
  534. END SUBROUTINE FILL_TARGET_ARRAY
  535. !EOC
  536. END SUBROUTINE CALC_EMISSION_SS
  537. !EOC
  538. END MODULE EMISSION_SS