emission_dms.F90 28 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776
  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_DMS
  14. !
  15. ! !DESCRIPTION:
  16. !\\
  17. !\\
  18. ! !INTERFACE:
  19. !
  20. MODULE EMISSION_DMS
  21. !
  22. ! !USES:
  23. !
  24. use GO, only : gol, goErr, goPr
  25. use tm5_distgrid, only : dgrid, get_distgrid, scatter, gather
  26. use partools, only : isRoot
  27. use dims, only : nregions
  28. use global_types, only : emis_data, d3_data
  29. IMPLICIT NONE
  30. PRIVATE
  31. !
  32. ! !PUBLIC MEMBER FUNCTIONS:
  33. !
  34. public :: emission_dms_init ! allocate
  35. public :: emission_dms_declare ! read monthly input
  36. public :: emission_dms_apply ! distribute & add emissions to tracer array
  37. public :: emission_dms_done ! deallocate
  38. public :: getDMS ! compute DMS from ocean
  39. !
  40. ! !PRIVATE DATA MEMBERS:
  41. !
  42. character(len=*), parameter :: mname = 'emission_dms'
  43. type(emis_data), dimension(nregions), target :: dms_land, dms_sea ! land and sea flux
  44. type(emis_data), dimension(nregions), target :: dms_tool ! work array
  45. real, dimension(:,:), allocatable :: dms_conc ! 1x1 ocean concentrations (nmol/L)
  46. ! The original calculation of the ocean DMS source was based on:
  47. ! - Climatology of sea water concentrations from Kettle et al., GBC, 1999
  48. ! - Exchange rate based on piston velocity from Liss and Merlivat, 1986
  49. ! (Air-sea gas exchange rates: Introduction and synthesis,
  50. ! The Role of Air-Sea Exchange in Geochemical Cycling,
  51. ! P. B. Menard, 113–127, D. Reidel, Norwell, Mass.)
  52. ! - The assumption that below -20 degrees C, there is 100% ice cover
  53. ! The updated scheme instead uses:
  54. ! - Climatology of sea water concentrations from Lana et al., GBC, 2011
  55. ! - Exchange rate based on piston velocity from Wanninkhof et al.,
  56. ! Limnol. Oceanogr. Methods, 2014.
  57. ! - Masking of emissions based on fractional ice cover directly
  58. !
  59. ! Optionally, it is possible to calculate the Schmidt number
  60. ! based on the actual sea surface temperature
  61. ! instead of the 2-m air temperature.
  62. ! However, approximating SST by t2m
  63. ! has only a small impact on the results:
  64. ! In a simulation driven by ERA-Interim for the year 2010,
  65. ! it reduces the global annual source of DMS by -2.3%,
  66. ! while the pattern remains nearly identical.
  67. ! For practical reasons (in particular because the current set of
  68. ! fields received in EC-Earth does not include SST),
  69. ! the 2-m temperature is therefore used by default,
  70. ! also in the new implementation.
  71. !
  72. ! We have also kept another simplification of the original implementations:
  73. ! - The sea water concentrations change discontinuously
  74. ! from month to month, following the input climatology.
  75. ! This could be improved by applying a linear interpolation
  76. ! between consecutive month.
  77. ! Given the short lifetimes of DMS and SO2 (on the order of a day),
  78. ! this may have a (small) impact on sulfate production.
  79. !
  80. ! The new scheme is turned on by default.
  81. ! For testing purposes, one can revert to the old scheme by setting use_old to true.
  82. logical, parameter :: use_old = .false. ! default value
  83. !logical, parameter :: use_old = .true. ! testing only
  84. !
  85. ! !REVISION HISTORY:
  86. ! 1 Oct 2010 - Achim Strunk - standardized routines name, new apply method
  87. ! 22 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  88. !
  89. ! !REMARKS:
  90. ! (1) see getDMS
  91. !
  92. !EOP
  93. !------------------------------------------------------------------------
  94. CONTAINS
  95. !--------------------------------------------------------------------------
  96. ! TM5 !
  97. !--------------------------------------------------------------------------
  98. !BOP
  99. !
  100. ! !IROUTINE: EMISSION_DMS_INIT
  101. !
  102. ! !DESCRIPTION: Allocate space needed to handle the emissions.
  103. !\\
  104. !\\
  105. ! !INTERFACE:
  106. !
  107. SUBROUTINE EMISSION_DMS_INIT( status )
  108. !
  109. ! !USES:
  110. !
  111. use dims, only : im, jm, iglbsfc
  112. use meteo, only : Set
  113. use meteodata, only : wspd_dat, ci_dat, lsmask_dat, sst_dat
  114. use meteodata, only : t2m_dat
  115. !
  116. ! !OUTPUT PARAMETERS:
  117. !
  118. integer, intent(out) :: status
  119. !
  120. ! !REVISION HISTORY:
  121. ! 1 Oct 2010 - Achim Strunk - extracted from old 'declare' routine
  122. ! 22 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  123. !
  124. !EOP
  125. !------------------------------------------------------------------------
  126. !BOC
  127. character(len=*), parameter :: rname = mname//'/Emission_DMS_Init'
  128. integer :: region
  129. integer :: imr, jmr, i1, i2, j1, j2
  130. ! --- begin --------------------------------------
  131. CALL GET_DISTGRID( dgrid(iglbsfc), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  132. allocate( dms_conc(i1:i2, j1:j2) )
  133. do region = 1, nregions
  134. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  135. imr = im(region) ; jmr = jm(region)
  136. allocate(dms_land(region)%surf(i1:i2,j1:j2))
  137. allocate(dms_sea (region)%surf(i1:i2,j1:j2))
  138. ! work array used in coarsen/scatter
  139. if (isRoot) then
  140. allocate(dms_tool(region)%surf(imr, jmr ))
  141. else
  142. allocate(dms_tool(region)%surf(1, 1))
  143. end if
  144. enddo
  145. ! Activate required meteo
  146. call Set( wspd_dat(iglbsfc), status, used=.true. )
  147. call Set( ci_dat(iglbsfc), status, used=.true. )
  148. call Set(lsmask_dat(iglbsfc), status, used=.true. )
  149. ! t2m is only used in old scheme inplace of sst:
  150. if (use_old) then
  151. call Set(t2m_dat(iglbsfc), status, used=.true. )
  152. else
  153. call Set(sst_dat(iglbsfc), status, used=.true. )
  154. endif
  155. status = 0
  156. END SUBROUTINE EMISSION_DMS_INIT
  157. !EOC
  158. !--------------------------------------------------------------------------
  159. ! TM5 !
  160. !--------------------------------------------------------------------------
  161. !BOP
  162. !
  163. ! !IROUTINE: EMISSION_DMS_DONE
  164. !
  165. ! !DESCRIPTION: Deallocate space needed to handle the emissions.
  166. !\\
  167. !\\
  168. ! !INTERFACE:
  169. !
  170. SUBROUTINE EMISSION_DMS_DONE( status )
  171. !
  172. ! !OUTPUT PARAMETERS:
  173. !
  174. integer, intent(out) :: status
  175. !
  176. ! !REVISION HISTORY:
  177. ! 1 Oct 2010 - Achim Strunk - rename old 'free_emission_dms'
  178. !
  179. !EOP
  180. !------------------------------------------------------------------------
  181. !BOC
  182. character(len=*), parameter :: rname = mname//'/Emission_DMS_Done'
  183. integer :: region
  184. ! --- begin ---------------------------------
  185. do region = 1, nregions
  186. deallocate( dms_land (region)%surf)
  187. deallocate( dms_tool (region)%surf)
  188. deallocate( dms_sea (region)%surf)
  189. end do
  190. deallocate(dms_conc)
  191. status = 0
  192. END SUBROUTINE EMISSION_DMS_DONE
  193. !EOC
  194. !--------------------------------------------------------------------------
  195. ! TM5 !
  196. !--------------------------------------------------------------------------
  197. !BOP
  198. !
  199. ! !IROUTINE: EMISSION_DMS_DECLARE
  200. !
  201. ! !DESCRIPTION: Opens, reads and evaluates input files (per month).
  202. ! Land surface fluxes from vegetation and soils
  203. ! are based on the dataset from Spiro et al. (JGR, 1992).
  204. ! Ocean concentrations are used to retrieve fluxes by
  205. ! temperature and wind speed.
  206. ! Emissions from biomass burning are provided for CMIP6
  207. ! but have not been included, because their contribution
  208. ! is likely very small. It would require extra coding
  209. ! to include these.
  210. !\\
  211. !\\
  212. ! !INTERFACE:
  213. !
  214. SUBROUTINE EMISSION_DMS_DECLARE( status )
  215. !
  216. ! !USES:
  217. !
  218. USE MDF, ONLY : MDF_Open, MDF_READ, MDF_Get_Var, MDF_Close
  219. USE MDF, ONLY : MDF_NETCDF, MDF_Inq_VarID
  220. use toolbox, only : coarsen_emission
  221. use dims, only : idate, sec_month, sec_year, iglbsfc, nlat180, nlon360
  222. use chem_param, only : xms
  223. use emission_data, only : msg_emis, emis_input_dir_dms
  224. !
  225. ! !OUTPUT PARAMETERS:
  226. !
  227. integer, intent(out) :: status
  228. !
  229. ! !REVISION HISTORY:
  230. ! 1 Oct 2010 - Achim Strunk - updated to emission standard
  231. ! 22 Jun 2012 - Ph. Le Sager - adapted for lon-lat MPI domain decomposition
  232. !
  233. !EOP
  234. !------------------------------------------------------------------------
  235. !BOC
  236. character(len=*), parameter :: rname = mname//'/declare_emission_dms'
  237. integer :: region, i1, i2, j1, j2, nlatsrc, nlonsrc, mo
  238. integer, parameter :: add_field=0
  239. integer, parameter :: amonth=2
  240. real :: year2month
  241. integer :: hid, varid
  242. character(len=8) :: vname8
  243. character(len=9) :: vname9
  244. ! work array on 1x1 resolution
  245. real, dimension(:,:), allocatable :: glb_dms_conc
  246. ! --- begin ---------------------------------
  247. year2month = sec_month/sec_year ! scale factor for yearly total to specific month
  248. nlatsrc = nlat180 ! or dgrid(iglbsfc)%jm_region, accessible thru get_distgrid(...)
  249. nlonsrc = nlon360 ! or dgrid(iglbsfc)%im_region, accessible thru get_distgrid(...)
  250. ! to deal with sea water at 1x1
  251. if(isRoot)then
  252. allocate(glb_dms_conc(nlonsrc,nlatsrc))
  253. else
  254. allocate(glb_dms_conc(1,1))
  255. end if
  256. !
  257. ! DMS emissions from landsurface
  258. ! DMS seawater concentrations. The actual flux is a function of windspeed and temperature
  259. ! and is calculated each time step. The resolution dependency needs to be checked.
  260. ! We will use the 10m winds on 1x1 resolution.....
  261. !
  262. write(gol,'(" EMISS-INFO ------------- read DMS emissions -------------")'); call goPr
  263. if (isRoot) then
  264. ! monthly emissions from Land Spiro dataset
  265. CALL MDF_Open( trim(emis_input_dir_dms)//'/DMSland.nc4', MDF_NETCDF, MDF_READ, hid, status )
  266. IF_NOTOK_RETURN(status=1)
  267. mo = idate(2) - 1
  268. if ( mo < 10 ) then
  269. write(vname8,'("spec1__",i1)')mo
  270. CALL MDF_Inq_VarID( hid, vname8, varid, status )
  271. IF_NOTOK_RETURN(status=1)
  272. else
  273. write(vname9,'("spec1__",i2)')mo
  274. CALL MDF_Inq_VarID( hid, vname9, varid, status )
  275. IF_NOTOK_RETURN(status=1)
  276. endif
  277. CALL MDF_Get_Var( hid, varid, glb_dms_conc, status )
  278. IF_NOTOK_RETURN(status=1)
  279. CALL MDF_Close( hid, status )
  280. IF_NOTOK_RETURN(status=1)
  281. call msg_emis(amonth, 'SPIRO', 'vegetation/soil', 'DMS', xms, sum(glb_dms_conc))
  282. call coarsen_emission('dms_land', nlonsrc, nlatsrc, glb_dms_conc, dms_tool, add_field, status)
  283. IF_NOTOK_RETURN(status=1)
  284. ! seawater concentrations
  285. if (use_old) then
  286. CALL MDF_Open( trim(emis_input_dir_dms)//'/DMSconc.nc4', MDF_NETCDF, MDF_READ, hid, status )
  287. IF_NOTOK_RETURN(status=1)
  288. mo = idate(2) - 1
  289. if ( mo < 10 ) then
  290. write(vname8,'("spec1__",i1)')mo
  291. CALL MDF_Inq_VarID( hid, vname8, varid, status )
  292. IF_NOTOK_RETURN(status=1)
  293. else
  294. write(vname9,'("spec1__",i2)')mo
  295. CALL MDF_Inq_VarID( hid, vname9, varid, status )
  296. IF_NOTOK_RETURN(status=1)
  297. endif
  298. CALL MDF_Get_Var( hid, varid, glb_dms_conc, status )
  299. IF_NOTOK_RETURN(status=1)
  300. else
  301. CALL MDF_Open( trim(emis_input_dir_dms)//'/DMS_ocean_conc.nc', MDF_NETCDF, MDF_READ, hid, status )
  302. IF_NOTOK_RETURN(status=1)
  303. CALL MDF_Inq_VarID( hid, trim('DMS'), varid, status )
  304. IF_ERROR_RETURN(status=1)
  305. CALL MDF_Get_Var( hid, varid, glb_dms_conc, status, start=(/1,1,idate(2)/) )
  306. IF_NOTOK_RETURN(status=1)
  307. endif
  308. CALL MDF_Close( hid, status )
  309. IF_NOTOK_RETURN(status=1)
  310. end if
  311. do region = 1, nregions
  312. call scatter(dgrid(region), dms_land(region)%surf, dms_tool(region)%surf, 0, status)
  313. IF_NOTOK_RETURN(status=1)
  314. end do
  315. call scatter(dgrid(iglbsfc), dms_conc, glb_dms_conc, 0, status)
  316. IF_NOTOK_RETURN(status=1)
  317. deallocate(glb_dms_conc)
  318. ! ok
  319. status = 0
  320. END SUBROUTINE EMISSION_DMS_DECLARE
  321. !EOC
  322. !--------------------------------------------------------------------------
  323. ! TM5 !
  324. !--------------------------------------------------------------------------
  325. !BOP
  326. !
  327. ! !IROUTINE: EMISSION_DMS_APPLY
  328. !
  329. ! !DESCRIPTION: Take monthly emissions, and
  330. ! - split them vertically
  331. ! - apply time splitting factors
  332. ! - add them up (add_3d)
  333. !\\
  334. !\\
  335. ! !INTERFACE:
  336. !
  337. SUBROUTINE EMISSION_DMS_APPLY( region, status )
  338. !
  339. ! !USES:
  340. !
  341. use dims, only : okdebug, itaur, nsrce, tref
  342. use dims, only : im, jm, lm
  343. use datetime, only : tau2date
  344. use emission_data, only : emission_vdist_by_sector
  345. use emission_data, only : do_add_3d
  346. use chem_param, only : idms, xmdms, xms
  347. use emission_data, only : vd_class_name_len
  348. !
  349. ! !INPUT PARAMETERS:
  350. !
  351. integer, intent(in) :: region
  352. !
  353. ! !OUTPUT PARAMETERS:
  354. !
  355. integer, intent(out) :: status
  356. !
  357. ! !REVISION HISTORY:
  358. ! 1 Oct 2010 - Achim Strunk - updated to use new vertical distribution method
  359. ! 28 Mar 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  360. !
  361. !EOP
  362. !------------------------------------------------------------------------
  363. !BOC
  364. character(len=*), parameter :: rname = mname//'/emission_dms_apply'
  365. integer,dimension(6) :: idater
  366. real :: dtime, fraction
  367. integer :: imr, jmr, lmr, i1, i2, j1, j2
  368. type(d3_data) :: emis3d
  369. character(len=vd_class_name_len) :: splittype
  370. ! --- begin -----------------------------------------
  371. if( okdebug ) then
  372. write(gol,*) 'start of emission_dms_apply'; call goPr
  373. endif
  374. call tau2date(itaur(region),idater)
  375. dtime=float(nsrce)/(2*tref(region)) !emissions are added in two steps...XYZECCEZYX.
  376. if( okdebug ) then
  377. write(gol,*) 'emission_dms_apply in region ',region,' at date: ',idater, ' with time step:', dtime ; call goPr
  378. endif
  379. ! get a working structure for 3d emissions
  380. call get_distgrid( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  381. allocate( emis3d%d3(i1:i2,j1:j2,lm(region)) ) ; emis3d%d3 = 0.0
  382. ! default: no additional splitting
  383. fraction = 1.0
  384. ! ----------------------------------------------------------------------------------------
  385. ! distinguish here between sectors and whether they should have additional splitting
  386. ! if( ar5_sectors(lsec)%catname == 'biomassburning' ) fraction = fraction * bb_frac etc...
  387. ! ----------------------------------------------------------------------------------------
  388. splittype = 'surface'
  389. ! -------
  390. ! dms sea
  391. ! -------
  392. ! vertically distribute according to sector
  393. call emission_vdist_by_sector( splittype, 'DMS', region, dms_sea(region), emis3d, status )
  394. IF_NOTOK_RETURN(status=1;deallocate(emis3d%d3))
  395. ! add dataset
  396. call do_add_3d( region, idms, i1, j1, emis3d%d3, xmdms, xms, status, fraction )
  397. IF_NOTOK_RETURN(status=1)
  398. ! -------
  399. ! dms land
  400. ! -------
  401. emis3d%d3 = 0.0
  402. ! vertically distribute according to sector
  403. call emission_vdist_by_sector( splittype, 'DMS', region, dms_land(region), emis3d, status )
  404. IF_NOTOK_RETURN(status=1;deallocate(emis3d%d3))
  405. ! add dataset
  406. call do_add_3d( region, idms, i1, j1, emis3d%d3, xmdms, xms, status, fraction )
  407. IF_NOTOK_RETURN(status=1)
  408. if( okdebug ) then
  409. write(gol,*) 'end of emission_dms_apply'; call goPr
  410. endif
  411. deallocate( emis3d%d3 )
  412. ! OK
  413. status = 0
  414. END SUBROUTINE EMISSION_DMS_APPLY
  415. !EOC
  416. !--------------------------------------------------------------------------
  417. ! TM5 !
  418. !--------------------------------------------------------------------------
  419. !BOP
  420. !
  421. ! !IROUTINE: GETDMS
  422. !
  423. ! !DESCRIPTION: calculate the DMS flux based on high resolution ocean
  424. ! concentration and wind speeds and temperatures,
  425. ! and calculate the 'coarser' resolution fluxes from it
  426. !
  427. ! input needed: U2 and T2 on 1x1
  428. !\\
  429. !\\
  430. ! !INTERFACE:
  431. !
  432. SUBROUTINE GETDMS( status )
  433. !
  434. ! !USES:
  435. !
  436. use toolbox, only : coarsen_emission
  437. use binas, only : T0
  438. use meteodata, only : t2m_dat, wspd_dat
  439. use meteodata, only : ci_dat, lsmask_dat, sst_dat
  440. use datetime, only : tau2date
  441. use Dims, only : itau, staggered, sec_year, sec_month, dxy11, iglbsfc, im, jm
  442. use chem_param, only : xms
  443. !
  444. ! !OUTPUT PARAMETERS:
  445. !
  446. integer, intent(out) :: status
  447. !
  448. ! !REVISION HISTORY:
  449. ! 22 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  450. !
  451. ! !REMARKS:
  452. ! (1) Required inputs: U2 and T2 on 1x1
  453. !
  454. !EOP
  455. !------------------------------------------------------------------------
  456. !BOC
  457. character(len=*), parameter :: rname = mname//'/getdms'
  458. real, dimension(:,:), allocatable :: emi_dms, emis_glb, wspd, temperature
  459. real :: xwind, zschmidt, resfac, tot_flux, kw, tt, xlon, xlat, piston_dms
  460. real :: prefac, norm, xsea, area_frac, emis_fac, k660, tt2
  461. integer :: i, j, l, region, i1, j1, i2, j2
  462. integer, dimension(6) :: idates
  463. integer, parameter :: level1=1, add_field=0
  464. integer :: nlatsrc, nlonsrc
  465. ! determine correction factor to account for less effective mixing
  466. ! at lower averaged windspeeds occuring at lower resolution
  467. CALL GET_DISTGRID( dgrid(iglbsfc), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  468. nlatsrc = jm(iglbsfc) ! or nlat180, or dgrid(iglbsfc)%jm_region, accessible thru get_distgrid(...)
  469. nlonsrc = im(iglbsfc) ! or nlon360, or dgrid(iglbsfc)%jm_region, accessible thru get_distgrid(...)
  470. if (iglbsfc /= 1) then ! assume no zoom region if region #1 is at 1x1 (see below)
  471. if (isRoot)then
  472. allocate(emis_glb (nlonsrc,nlatsrc))
  473. else
  474. allocate(emis_glb (1,1))
  475. end if
  476. end if
  477. allocate(emi_dms (i1:i2, j1:j2))
  478. allocate(wspd (i1:i2, j1:j2))
  479. allocate(temperature(i1:i2, j1:j2))
  480. call tau2date(itau+staggered,idates) !CMK at end 3-hourly interval
  481. if (use_old) then
  482. ! would be better to use SST in the old scheme also
  483. ! but using air temperature allows to reproduce the old implementation
  484. temperature = t2m_dat(iglbsfc)%data(:,:,1)
  485. else
  486. temperature = sst_dat(iglbsfc)%data(:,:,1)
  487. endif
  488. ! convert to celsius
  489. temperature = temperature - T0
  490. wspd = wspd_dat(iglbsfc)%data(:,:,1)
  491. ! Resolution dependent tuning factor
  492. ! that can be used to account for variability in wind speed
  493. ! Currently set to 1.
  494. resfac=1.
  495. if (use_old) then
  496. ! THE RESOLUTION DEPENDENCY WITH PREVIOUS TM3 RUNS WAS AS FOLLOWS
  497. ! NEEDS TO BE CHECKED FD
  498. !FD if (jm.eq.24) resfac=1.32
  499. !FD if (jm.eq.48) resfac=1.065
  500. emi_dms=0.
  501. do j=j1,j2
  502. do i=i1,i2
  503. !
  504. ! ocean surface temperature is approximated by air surface temperature
  505. ! we assume that if temperature<-20. ! sea ice prevents DMS emissions (sea-ice data from ECMWF would be usefull)
  506. ! max ocean temperature is 28.0 C
  507. ! minumum ocean temperature 0. C
  508. !fd tt=min(28.,ts(i,j)-273.15) ! the original parameterisation needs the seawater temperature
  509. ! as a proxy we use the air temperature, but restrict it to 28 C.
  510. ! surface or 2m temperature of air is probably the best proxy.
  511. tt=min(28.0,temperature(i,j))
  512. if (dms_conc(i,j).gt.0..and.tt.gt.-20.) then
  513. ! TvN: don't know where this expression comes from
  514. zschmidt=max(0.0,3652.047271-246.99*TT+8.536397*TT*TT-0.124397*TT*TT*TT)
  515. xwind=wspd(i,j)
  516. !
  517. ! Liss and Merlivat windspeed dependency of transfer velocity kw
  518. ! including correction for Schmidt number
  519. ! windspeed at middle of gridbox is used should be 10 m windspeed [m/s]
  520. if (xwind.le.3.6) then
  521. kw=0.17*xwind
  522. piston_dms=kw*(zschmidt/600.)**(-0.666)
  523. else !xwind>3.6)
  524. if (xwind.le.13.0) then
  525. kw=2.85*xwind-9.65
  526. piston_dms=kw*(zschmidt/600.)**(-0.5)
  527. else !xwind>13.0
  528. kw=5.9*xwind-49.3
  529. piston_dms=kw*(zschmidt/600.)**(-0.5)
  530. endif !x<13.0
  531. endif
  532. !
  533. ! cm/hr=> m/s * nmol/l=>mol/m3 * m2 * g S/mol=>g S/s
  534. !
  535. emi_dms(i,j)=resfac*piston_dms/360000.*dms_conc(i,j)*1e-6*dxy11(j)*xms
  536. endif !(dms_conc gt 0)
  537. enddo
  538. enddo
  539. else
  540. ! New parameterization based on Lana et al. and Wanninkhof et al.
  541. emi_dms=0.
  542. ! 1.e-6 converts nmol/L to mol/m3
  543. ! 1.e-2/3600. converts cm/hr to m/s
  544. ! xms converts mol to g S
  545. prefac = resfac * 1.e-6 * xms * 1.e-2 / 3600.
  546. do j=j1,j2
  547. norm = prefac * dxy11(j)
  548. do i=i1,i2
  549. ! Undefined values in the original input files
  550. ! of the DMS sea water concentrations have been set to -1.e9.
  551. ! These should be excluded.
  552. ! Locations with zero concentrations can be excluded as well.
  553. if (dms_conc(i,j) <= 0. ) cycle
  554. ! sea fraction
  555. xsea=1.-lsmask_dat(iglbsfc)%data(i,j,1)/100.
  556. ! DMS is emitted only over sea without ice cover
  557. area_frac = xsea * (1.-ci_dat(iglbsfc)%data(i,j,1))
  558. if (area_frac .lt. 1.e-10 ) cycle
  559. emis_fac = norm * area_frac
  560. ! Wind speed dependence from Wanninkhof et al.,
  561. ! Limnology and Oceanography: Methods, 351-362, 2014.
  562. ! "It should provide good estimates for most insoluble
  563. ! gases at intermediate wind speed ranges (3-15 m s–1).
  564. ! At low winds, non-wind effects such as chemical enhancement
  565. ! and thermal boundary layer processes influence gas transfer,
  566. ! and this quadratic relationship will underestimate gas transfer.
  567. ! At high winds (≈> 15 m s–1), bubble-enhanced exchange will
  568. ! affect gases differently depending on their solubility,
  569. ! and the relationship is only suitable for CO2
  570. ! under these conditions.
  571. ! The differences in physical and chemical processes
  572. ! in boundary layers and their impact on gases at high and low winds
  573. ! need to be taken into consideration when estimating uncertainty.
  574. ! Since over 94% of the winds over the ocean are in the 3-15 m s–1 range,
  575. ! the regional and global gas transfer velocities
  576. ! for gases listed in Table 1 can be determined
  577. ! using the above relationship."
  578. k660 = 0.251 * wspd(i,j)**2
  579. ! where k660 is in units cm/h.
  580. ! We prefer this relation above that proposed
  581. ! by Nightingale et al. (Glob. Biogeochem. Cycl., 2000):
  582. ! k600 = 0.333 * u_10m + 0.222 * (u_10m)^2
  583. ! The presence of the linear term enhances the transfer coefficient
  584. ! at low wind speeds, but this regime is of minor importance.
  585. ! The relation by Wanninkhof gives a higher transfer coefficient
  586. ! at moderate to high wind speeds.
  587. ! Because wind speed tends to be underestimated
  588. ! because of the coarse horizontal resolution,
  589. ! this could be an argument to use the formulation
  590. ! by Wanninkhof et al., which is also more recent.
  591. ! In both formulations, the transfer coefficient depends on temperature
  592. ! as the square root of the Schmidt number.
  593. ! We calculate the Schmidt number based on the fourth-order polynomial fit
  594. ! given in Wanninkhof et al. (2014) for sea water at 3.5% salinity.
  595. ! It is based on the data by Saltzman et al. (JGR, 1993),
  596. ! but is applicable over over a wider range of temperatures (–2°C to 40°C)
  597. ! than the third-order polynomial proposed by Saltzman et al. (0°C and 30°C).
  598. ! At the high end, the fourth order polynomial only weakly depends
  599. ! on temperature, and reaches a minimum at 43.78°C.
  600. ! Since the increase thereafter falls outside the range of the fit,
  601. ! we set the Schmidt number at higher temperatures
  602. ! equal to the value at 43.78°C.
  603. !
  604. ! Using the temperature dependence from Wanninkhof et al.,
  605. ! the transfer coefficient increases by less than about 3.2%
  606. ! per degree temperature increase.
  607. tt=min(43.78,temperature(i,j))
  608. tt2=tt*tt
  609. zschmidt=2855.7-177.63*tt+6.0438*tt2-0.11645*tt2*tt+0.00094743*tt2*tt2
  610. kw = k660 * (660./zschmidt)**0.5
  611. ! "DMS fluxes are generally parameterized assuming water side resitance only,
  612. ! but air side resistance can also be significant at cold temperatures
  613. ! and high wind speeds." (Lana et al.)
  614. ! Here the air side resistance is neglected:
  615. piston_dms = kw
  616. emi_dms(i,j)=emis_fac*piston_dms*dms_conc(i,j)
  617. enddo
  618. enddo
  619. endif
  620. !-----------------------------------------------------------
  621. ! Fill target array
  622. !-----------------------------------------------------------
  623. ! Test for optimization: if region #1 is at 1x1, assume no zoom region and skip call to coarsen and mpi comm
  624. if (iglbsfc == 1) then
  625. dms_sea(1)%surf = emi_dms*sec_month*1e-3 !expected in kgS/month
  626. else
  627. !-----------------------------------------------------------
  628. ! gather on 1x1, coarsen to other grid on root, then scatter
  629. !-----------------------------------------------------------
  630. ! FIXME-PERF : Need a coarsen that works independtly on each proc
  631. call gather(dgrid(iglbsfc), emi_dms, emis_glb, 0, status)
  632. IF_NOTOK_RETURN(status=1)
  633. if (isRoot) then
  634. ! FIXME : if really needed, write into budget - here, it gives too much output
  635. !write(gol,*)'the yearly sum of GETDMS would be ',sum(emis_glb) *1e-12*sec_year,' Tg S/a'; call goPr
  636. ! assume equal emissions over season
  637. call coarsen_emission('emi_dms', nlonsrc, nlatsrc, emis_glb, dms_tool, add_field, status)
  638. IF_NOTOK_RETURN(status=1)
  639. end if
  640. do region = 1, nregions
  641. call scatter(dgrid(region), dms_sea(region)%surf, dms_tool(region)%surf, 0, status)
  642. IF_NOTOK_RETURN(status=1)
  643. dms_sea(region)%surf = dms_sea(region)%surf*sec_month*1e-3 !expected in kgS/month
  644. enddo
  645. deallocate(emis_glb)
  646. end if
  647. !-----------------------------------------------------------
  648. deallocate(emi_dms)
  649. deallocate(temperature)
  650. deallocate(wspd)
  651. END SUBROUTINE GETDMS
  652. !EOC
  653. END MODULE EMISSION_DMS