emission_ss.F90 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484
  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 : u10m_dat, v10m_dat, ci_dat, lsmask_dat
  75. !
  76. ! !OUTPUT PARAMETERS:
  77. !
  78. INTEGER, INTENT(out) :: status
  79. !
  80. ! !REVISION HISTORY:
  81. ! 1 Jul 2013 - Ph Le Sager - v0
  82. !
  83. ! !REMARKS:
  84. !
  85. !EOP
  86. !------------------------------------------------------------------------
  87. !BOC
  88. CALL Set( u10m_dat(iglbsfc), status, used=.TRUE. )
  89. CALL Set( v10m_dat(iglbsfc), status, used=.TRUE. )
  90. CALL Set( ci_dat(iglbsfc), status, used=.TRUE. )
  91. CALL Set( lsmask_dat(iglbsfc), status, used=.TRUE. )
  92. END SUBROUTINE emission_ss_init
  93. !EOC
  94. !--------------------------------------------------------------------------
  95. ! TM5 !
  96. !--------------------------------------------------------------------------
  97. !BOP
  98. !
  99. ! !IROUTINE: CALC_EMISSION_SS
  100. !
  101. ! !DESCRIPTION: Calculate emitted numbers and mass as function of the ten-meter
  102. ! wind speed as described in Vignati et al. (2010) and Gong
  103. ! (2003). Sea salt is emitted in both the soluble accumulation
  104. ! mode and the soluble coarse mode.
  105. !
  106. ! Ref: Vignati et al. (2010) : Global scale emission and
  107. ! distribution of sea-spray aerosol: Sea-salt and organic
  108. ! enrichment, Atmos. Environ., 44, 670-677,
  109. ! doi:10.1016/j.atmosenv.2009.11.013
  110. !
  111. ! Gong (2003) : A parameterization of sea-salt aerosol source
  112. ! function for sub- and super-micron particles, Global
  113. ! Biogeochem. Cy., 17, 1097, doi:10.1029/2003GB002079
  114. !\\
  115. !\\
  116. ! !INTERFACE:
  117. !
  118. SUBROUTINE calc_emission_ss( status )
  119. !
  120. ! !USES:
  121. !
  122. USE Binas, ONLY : pi
  123. USE dims, ONLY : okdebug, itaur, nsrce, tref
  124. USE datetime, ONLY : tau2date
  125. USE dims, ONLY : nlon360, nlat180, idate, itau, nregions, staggered, dxy11
  126. USE dims, ONLY : sec_month, iglbsfc, im, jm, lm
  127. USE chem_param, ONLY : sigma_lognormal, ss_density, nmodes, mode_acs, mode_cos
  128. USE chem_param, ONLY : iacs_n, icos_n, issacs, isscos, radius_ssa, radius_ssc
  129. USE toolbox, ONLY : coarsen_emission
  130. USE emission_data, ONLY : emis_mass, emis_number, emis_temp, msg_emis
  131. USE partools, ONLY : isRoot
  132. USE tm5_distgrid, ONLY : dgrid, get_distgrid, scatter, gather
  133. USE meteodata, ONLY : u10m_dat, v10m_dat, ci_dat, lsmask_dat
  134. USE emission_data, ONLY : vd_class_name_len, emission_vdist_by_sector
  135. !
  136. ! !OUTPUT PARAMETERS:
  137. !
  138. INTEGER, INTENT(out) :: status
  139. !
  140. ! !REVISION HISTORY:
  141. ! 1 Oct 2010 - Achim Strunk - introduced vertical splitting
  142. ! 25 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  143. !
  144. ! !REMARKS:
  145. ! (1) done on 1x1 grid.
  146. ! (2) this routine basicaly is the "declare" routine for sea-salt.
  147. !
  148. !EOP
  149. !------------------------------------------------------------------------
  150. !BOC
  151. CHARACTER(len=*), PARAMETER :: rname = mname//'/calc_emission_ss'
  152. REAL, DIMENSION(:,:), ALLOCATABLE :: number, mass, glb_arr
  153. INTEGER :: region
  154. REAL :: xwind, rg1, rg2, dens, xsea
  155. !>>> TvN
  156. REAL :: norm
  157. REAL, DIMENSION(:,:), ALLOCATABLE :: emis_fac
  158. REAL :: area_frac
  159. ! ratio of prefactors in Eqs. (2) and (1) of Salisbury et al. (GRL, 2014):
  160. REAL, PARAMETER :: W10_fac = 4.60e-3/3.84e-4
  161. REAL, PARAMETER :: W10_exp = 2.26
  162. !<<< TvN
  163. INTEGER, PARAMETER :: add_field = 0
  164. INTEGER, PARAMETER :: amonth = 2
  165. INTEGER :: i, j, i1, i2, j1, j2, i01, i02, j01, j02
  166. TYPE(d3_data) :: emis3d
  167. TYPE(emis_data), DIMENSION(nregions) :: emis_glb
  168. CHARACTER(len=vd_class_name_len) :: splittype
  169. ! --- begin -----------------------------------------
  170. ! vertical splitting is that class
  171. splittype = 'surface'
  172. !===================
  173. ! Work arrays
  174. !===================
  175. CALL GET_DISTGRID( dgrid(iglbsfc), I_STRT=i01, I_STOP=i02, J_STRT=j01, J_STOP=j02 )
  176. ALLOCATE(number(i01:i02,j01:j02))
  177. ALLOCATE(mass (i01:i02,j01:j02))
  178. !>>> TvN
  179. ALLOCATE(emis_fac(i01:i02,j01:j02)) ; emis_fac=0.
  180. !<<< TvN
  181. ! no zoom region if region #1 is at 1x1
  182. IF ( (iglbsfc /= 1) .OR. okdebug) THEN
  183. IF (isRoot) THEN
  184. ALLOCATE(glb_arr(nlon360,nlat180))
  185. DO region = 1, nregions
  186. ALLOCATE(emis_glb(region)%surf(im(region), jm(region)))
  187. END DO
  188. ELSE
  189. ALLOCATE(glb_arr(1,1))
  190. DO region = 1, nregions
  191. ALLOCATE(emis_glb(region)%surf(1,1))
  192. END DO
  193. END IF
  194. END IF
  195. !>>> TvN
  196. ! The parameterization of Gong et al. (2003)
  197. ! gives the particle number flux as a function
  198. ! of the radius and the 10-m wind speed u_10:
  199. ! dF/dr80 = f(u_10) x g(r80).
  200. ! The dependence on wind speed is given by
  201. ! the power law f(u_10) = u_10^3.41,
  202. ! and does not affect the size distribution.
  203. ! r80 is the radius at 80% humidity,
  204. ! which is almost exactly 2.0 times the dry radius
  205. ! (Lewis and Schwartz, Sea Salt Aerosol Production, 2004).
  206. ! Note also that in Eq. (2) of Gong et al.
  207. ! B is defined in terms of log(10,r) not ln(r).
  208. !
  209. ! The number mean radii defined in chem_param.F90,
  210. ! i.e. 0.09 um for the accumulation mode,
  211. ! and 0.794 for the coarse mode,
  212. ! were obtained by fitting an accumulation
  213. ! and coarse mode to the dF/dln(r),
  214. ! with r the dry particle radius
  215. ! (see presentation E. Vignati, 21 December 2005).
  216. ! It was verified that the size distribution of Gong et al.
  217. ! can be reasonable well described using these mode radii.
  218. !
  219. ! It is not totally clear how the corresponding prefactors
  220. ! for the two modes were obtained.
  221. ! One way to calculate these prefactors,
  222. ! is to define two size ranges
  223. ! and require that the numbers of particles emitted
  224. ! in both ranges correspond to Gong et al.
  225. ! This procedure is similar to that used in emission_dust.F90,
  226. ! but for particle number instead of mass.
  227. ! Using ranges r1 and r2 for the dry particle radii
  228. ! with r1 from 0.05 to 0.5 um and r2 from 0.5 to 5 um,
  229. ! the resulting prefactors are 9.62013e-3 and 9.05658e-4
  230. ! for the accumulation and coarse mode, respectively.
  231. ! These numbers are very close to the values
  232. ! of 9e-3 and 9e-4, obtained by E. Vignati.
  233. ! They are also insensitive to the
  234. ! value used for the upper bound of r2.
  235. ! (Using a value of 10 instead of 5 um,
  236. ! the prefactors would be 9.62536e-3 and 8.91184e-3.)
  237. !
  238. DO j=j01,j02
  239. norm = 1.e4 * dxy11(j) * sec_month
  240. DO i=i01,i02
  241. ! sea fraction
  242. xsea=1.-lsmask_dat(iglbsfc)%data(i,j,1)/100.
  243. ! sea salt is emitted only over sea without ice cover
  244. area_frac = xsea * (1.-ci_dat(iglbsfc)%data(i,j,1))
  245. if (area_frac .lt. 1.e-10) cycle
  246. emis_fac(i,j) = norm * area_frac
  247. ! Wind speed dependence following W10 parameterization
  248. ! of Salisbury et al. (JGR, 2013; GRL, 2014),
  249. ! which replaces the dependence of Gong et al.
  250. ! Salisbury et al. suggest that "at this stage ...
  251. ! use of W10 is preferable to W37".
  252. ! In TM5, W10 gives better results than W37,
  253. ! which leads to high AOD values compared to MODIS C6,
  254. ! at mid- and high latitudes.
  255. ! The same is true for the wind dependence
  256. ! proposed by Albert et al. (ACPD, 2015).
  257. xwind = SQRT(u10m_dat(iglbsfc)%data(i,j,1)**2+v10m_dat(iglbsfc)%data (i,j,1)**2)
  258. emis_fac(i,j) = emis_fac(i,j) * W10_fac * xwind**W10_exp
  259. ENDDO
  260. ENDDO
  261. !<<< TvN
  262. !===================
  263. ! Accumulation mode
  264. !===================
  265. ! emitted numbers
  266. ! ---------------
  267. DO j=j01,j02
  268. DO i=i01,i02
  269. !>>> TvN
  270. ! sea fraction
  271. !xsea=1.-lsmask_dat(iglbsfc)%data(i,j,1)/100.
  272. !xwind=SQRT(u10m_dat(iglbsfc)%data(i,j,1)**2+v10m_dat(iglbsfc)%data (i,j,1)**2)
  273. !number(i,j) = 0.009*xwind**(3.4269)*1e4*dxy11(j)*xsea*(1.-ci_dat(iglbsfc)%data(i,j,1))*sec_month ! #/gridbox/month
  274. number(i,j) = emis_fac(i,j)*9.62013e-3 ! #/gridbox/month
  275. !<<< TvN
  276. END DO
  277. END DO
  278. CALL fill_target_array( number, 'number acc', 'ss_number acc' ) ! fill emis_temp(region)
  279. IF_NOTOK_RETURN(status=1)
  280. ! vertically distribute according to sector
  281. DO region = 1, nregions
  282. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  283. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  284. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  285. emis_number(region,mode_acs)%d4(:,:,:,4) = emis3d%d3 !#/grid/month
  286. DEALLOCATE(emis3d%d3)
  287. END DO
  288. ! emitted mass
  289. ! ------------
  290. dens = ss_density*0.001 ! in g/cm3
  291. rg1 = radius_ssa*1e6 ! in micron
  292. mass(:,:) = number(:,:)*pi*4./3. *(rg1*1e-4)**3 * EXP(9./2.*alog(sigma_lognormal(3))**2)*dens*1e-3 ! kg/gridbox/month
  293. CALL fill_target_array( mass, 'mass acc', 'ss_mass acc' ) ! fill emis_temp(region)
  294. IF_NOTOK_RETURN(status=1)
  295. ! vertically distribute according to sector
  296. DO region = 1, nregions
  297. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  298. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  299. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  300. emis_mass (region,mode_acs)%d4(:,:,:,4) = emis3d%d3 !kg/grid/month
  301. DEALLOCATE(emis3d%d3)
  302. END DO
  303. !===================
  304. ! Coarse mode
  305. !===================
  306. ! emitted numbers
  307. ! ---------------
  308. DO j=j01,j02
  309. DO i=i01,i02
  310. ! >>> TvN
  311. ! sea fraction
  312. !xsea=1.-lsmask_dat(iglbsfc)%data(i,j,1)/100.
  313. !xwind=SQRT(u10m_dat(iglbsfc)%data(i,j,1)**2+v10m_dat(iglbsfc)%data (i,j,1)**2)
  314. !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
  315. number(i,j) = emis_fac(i,j)*9.05658e-4 !#/cm2/s-->#/gridbox/month
  316. ! <<< TvN
  317. ENDDO
  318. ENDDO
  319. CALL fill_target_array( number, 'number coa', 'ss_number coa' ) ! fill emis_temp(region)
  320. IF_NOTOK_RETURN(status=1)
  321. ! vertically distribute according to sector
  322. DO region = 1, nregions
  323. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  324. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  325. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  326. emis_number(region,mode_cos)%d4(:,:,:,4) = emis3d%d3 !#/grid/month
  327. DEALLOCATE(emis3d%d3)
  328. END DO
  329. ! emitted mass
  330. ! ------------
  331. dens = ss_density*0.001 ! im g/cm3
  332. rg2 = radius_ssc*1e6 ! in micron
  333. mass = number(:,:)*pi*4./3. *(rg2*1e-4)**3 * EXP(9./2.*alog(sigma_lognormal(4))**2)*dens*1e-3 ! kg/gridbox/month
  334. CALL fill_target_array( mass, 'mass coa', 'ss_mass coa' ) ! fill emis_temp(region)
  335. IF_NOTOK_RETURN(status=1)
  336. ! vertically distribute according to sector
  337. DO region = 1, nregions
  338. CALL get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  339. ALLOCATE( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  340. CALL emission_vdist_by_sector( splittype, 'SS', region, emis_temp(region), emis3d, status )
  341. emis_mass(region,mode_cos)%d4(:,:,:,4) = emis3d%d3 !kg/gridbox/month
  342. DEALLOCATE(emis3d%d3)
  343. END DO
  344. !===================
  345. ! Done
  346. !===================
  347. DEALLOCATE(number, mass)
  348. DEALLOCATE(emis_fac)
  349. DO region = 1, nregions
  350. IF (ASSOCIATED(emis_glb(region)%surf)) DEALLOCATE(emis_glb(region)%surf)
  351. END DO
  352. IF (ALLOCATED(glb_arr)) DEALLOCATE(glb_arr)
  353. CONTAINS
  354. !--------------------------------------------------------------------------
  355. ! TM5 !
  356. !--------------------------------------------------------------------------
  357. !BOP
  358. !
  359. ! !IROUTINE: FILL_TARGET_ARRAY
  360. !
  361. ! !DESCRIPTION: Convenient internal program to fill EMIS_TEMP - same as routine in emission_dust.F90
  362. !\\
  363. !\\
  364. ! !INTERFACE:
  365. !
  366. SUBROUTINE fill_target_array( arr_in, str1, str2 )
  367. !
  368. ! !INPUT PARAMETERS:
  369. !
  370. REAL, INTENT(in) :: arr_in(i01:,j01:)
  371. CHARACTER(len=*), INTENT(in) :: str1, str2
  372. !
  373. ! !REVISION HISTORY:
  374. ! 25 Jun 2012 - P. Le Sager - v0
  375. !
  376. !EOP
  377. !------------------------------------------------------------------------
  378. !BOC
  379. ! Test for optimization: if region #1 is at 1x1, assume no zoom region
  380. ! and skip call to coarsen and mpi comm
  381. IF (iglbsfc == 1) THEN
  382. emis_temp(1)%surf = arr_in
  383. IF (okdebug) THEN
  384. CALL gather(dgrid(iglbsfc), arr_in, glb_arr, 0, status)
  385. IF_NOTOK_RETURN(status=1)
  386. !FIXME - really needed? too much messaging
  387. ! IF (isRoot) THEN
  388. ! CALL msg_emis( amonth, ' ', TRIM(str1), 'SS', 1., SUM(glb_arr) )
  389. ! END IF
  390. END IF
  391. ELSE
  392. DO region = 1, nregions
  393. emis_temp(region)%surf = 0.0
  394. END DO
  395. ! gather on 1x1, coarsen to other grid on root, then scatter
  396. !-----------------------------------------------------------
  397. ! FIXME-PERF : Need a coarsen that works independtly on each proc, regardless of zooming... :(
  398. CALL gather(dgrid(iglbsfc), arr_in, glb_arr, 0, status)
  399. IF_NOTOK_RETURN(status=1)
  400. IF (isRoot) THEN
  401. !FIXME - really needed? too much messaging
  402. ! CALL msg_emis( amonth, ' ', TRIM(str1), 'SS', 1., SUM(glb_arr) )
  403. CALL coarsen_emission(TRIM(str2), nlon360, nlat180, glb_arr, emis_glb, add_field, status)
  404. IF_NOTOK_RETURN(status=1)
  405. END IF
  406. DO region = 1, nregions
  407. CALL scatter(dgrid(region), emis_temp(region)%surf, emis_glb(region)%surf, 0, status)
  408. IF_NOTOK_RETURN(status=1)
  409. ENDDO
  410. ENDIF
  411. status=0
  412. END SUBROUTINE FILL_TARGET_ARRAY
  413. !EOC
  414. END SUBROUTINE CALC_EMISSION_SS
  415. !EOC
  416. END MODULE EMISSION_SS