calc_pm.F90 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309
  1. #define TRACEBACK write (gol,'("in ",a," (",a,i6,")")') rname, __FILE__, __LINE__ ; call goErr
  2. #define IF_NOTOK_RETURN(action) if (status/=0) then; TRACEBACK; action; return; end if
  3. #define IF_ERROR_RETURN(action) if (status> 0) then; TRACEBACK; action; return; end if
  4. !
  5. #include "tm5.inc"
  6. !-----------------------------------------------------------------------------
  7. ! TM5 !
  8. !-----------------------------------------------------------------------------
  9. !BOP
  10. !
  11. ! !MODULE: CALC_PM
  12. !
  13. ! !DESCRIPTION: PM module to calculate particulate matter concentrations for
  14. ! user defined aerodynamic diameters
  15. !\\
  16. !\\
  17. ! !INTERFACE:
  18. !
  19. MODULE CALC_PM
  20. !
  21. ! !USES:
  22. !
  23. USE GO, only : gol, goErr, goPr
  24. IMPLICIT NONE
  25. PRIVATE
  26. !
  27. ! !PUBLIC MEMBER FUNCTIONS:
  28. !
  29. Public :: PMX_Integrate_3d,PMx_integrate_0d
  30. !
  31. ! !PRIVATE DATA MEMBERS:
  32. !
  33. character(len=*), parameter :: mname = 'calc_pm'
  34. Real, Parameter :: hr2 = 0.5*sqrt(2.0)
  35. !
  36. ! !REVISION HISTORY:
  37. ! 1 Oct 2010 - Achim Strunk - Revised and made suitable for 3D and stations
  38. ! 19 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  39. !
  40. ! !REMARKS:
  41. ! (1) module used only if with_m7 is used
  42. !
  43. !EOP
  44. !------------------------------------------------------------------------
  45. CONTAINS
  46. !--------------------------------------------------------------------------
  47. ! TM5 !
  48. !--------------------------------------------------------------------------
  49. !BOP
  50. !
  51. ! !IROUTINE: PMX_Integrate_0d
  52. !
  53. ! !DESCRIPTION: This routine generates a pmx value for a given grid cell
  54. ! Arbitrary radii may be specified.
  55. ! NO3 and NH4 come from EQSAM and are considered to be in
  56. ! the accumulation mode.
  57. !\\
  58. !\\
  59. ! !INTERFACE:
  60. !
  61. Subroutine PMX_Integrate_0d(region, is, js, ls, isizelimit, pmx, status)
  62. !
  63. ! !USES:
  64. !
  65. Use tracer_data, only : mass_dat
  66. ! I need to convert (/kg air) to (/m3 air)
  67. use chem_param, only : iso4nus,isoanus
  68. use chem_param, only : iso4ais, ibcais, ipomais, isoaais
  69. use chem_param, only : iso4acs, ibcacs, ipomacs, isoaacs, issacs, iduacs
  70. use chem_param, only : iso4cos, ibccos, ipomcos, isoacos, isscos, iducos
  71. use chem_param, only : ibcaii, ipomaii, isoaaii, iduaci
  72. use chem_param, only : iducoi, ino3_a, inh4
  73. Use mo_aero_m7, only : nmod, cmedr2mmedr, sigma
  74. Use m7_data, only : h2o_mode, rw_mode
  75. ! To rewrite ug / cell to ug/m3
  76. Use global_data, only : region_dat
  77. Use Meteodata, only : gph_dat
  78. !
  79. ! !INPUT PARAMETERS:
  80. !
  81. Integer, Intent(In) :: region
  82. real, intent(in) :: isizelimit
  83. integer, Intent(in) :: is, js, ls
  84. !
  85. ! !OUTPUT PARAMETERS:
  86. !
  87. Integer, Intent(out) :: status
  88. real, Intent(out) :: pmx ! units = 'kg/m3'
  89. !
  90. ! !REVISION HISTORY:
  91. ! 7 Oct 2010 - Achim Strunk - v0
  92. ! 19 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  93. !
  94. !EOP
  95. !------------------------------------------------------------------------
  96. !BOC
  97. character(len=*), parameter :: rname = mname//'/PMX_Integrate_0d'
  98. Integer :: imod
  99. Real :: z, rad, sizelimit
  100. Real :: volume ! For converting mass to mass/volume
  101. Real, Dimension(nmod) :: modfrac ! Length: Number of M7 modes (7)
  102. Real, Dimension(nmod) :: lnsigma
  103. ! --- begin -------------------------------------
  104. lnsigma = log(sigma)
  105. ! We can convert micrometers diameter to meters radius. One micrometer diameter is 5.e-7 meter radius.
  106. sizelimit = isizelimit * 5.e-7
  107. ! initialise target value
  108. pmx = 0.0
  109. Do imod = 1, nmod
  110. ! Calculate the median radius of the volume weighted distribution.
  111. rad = rw_mode(region,imod)%d3(is,js,ls) * cmedr2mmedr(imod)
  112. !! if( rad == 0.0 ) then ! changed to a precision depending value
  113. if( rad < 100*TINY(1.0) ) then
  114. modfrac(imod) = 1.0 ! Should not matter, because if the radius is zero,
  115. ! there are no aerosols. But in principle, it would
  116. ! mean that all aerosols are infinitely small, so
  117. ! they would fit in any PM-class.
  118. else
  119. z = ( log(sizelimit) - log(rad) ) / lnsigma(imod)
  120. modfrac(imod) = 0.5 + 0.5 * ERF(z*hr2)
  121. end if
  122. end do
  123. ! And the volume. This is the nicest and cleanest line of code of this whole module.
  124. volume = ( gph_dat(region)%data(is,js,ls+1) - gph_dat(region)%data(is,js,ls) ) * &
  125. region_dat(region)%dxyp(js)
  126. pmx = pmx + mass_dat(region)%rm(is,js,ls, iso4nus) * modfrac(1) / volume
  127. pmx = pmx + mass_dat(region)%rm(is,js,ls, iso4nus) * modfrac(1) / volume
  128. pmx = pmx + mass_dat(region)%rm(is,js,ls, iso4ais) * modfrac(2) / volume
  129. pmx = pmx + mass_dat(region)%rm(is,js,ls, iso4acs) * modfrac(3) / volume
  130. pmx = pmx + mass_dat(region)%rm(is,js,ls, iso4cos) * modfrac(4) / volume
  131. ! preparing for soa
  132. ! pmx = pmx + mass_dat(region)%rm(is,js,ls, isoanus) * modfrac(1) / volume
  133. ! pmx = pmx + mass_dat(region)%rm(is,js,ls, isoanus) * modfrac(1) / volume
  134. ! pmx = pmx + mass_dat(region)%rm(is,js,ls, isoaais) * modfrac(2) / volume
  135. ! pmx = pmx + mass_dat(region)%rm(is,js,ls, isoaacs) * modfrac(3) / volume
  136. ! pmx = pmx + mass_dat(region)%rm(is,js,ls, isoacos) * modfrac(4) / volume
  137. pmx = pmx + mass_dat(region)%rm(is,js,ls, ibcais ) * modfrac(2) / volume
  138. pmx = pmx + mass_dat(region)%rm(is,js,ls, ibcacs ) * modfrac(3) / volume
  139. pmx = pmx + mass_dat(region)%rm(is,js,ls, ibccos ) * modfrac(4) / volume
  140. pmx = pmx + mass_dat(region)%rm(is,js,ls, ibcaii ) * modfrac(5) / volume
  141. pmx = pmx + mass_dat(region)%rm(is,js,ls, ipomais) * modfrac(2) / volume
  142. pmx = pmx + mass_dat(region)%rm(is,js,ls, ipomacs) * modfrac(3) / volume
  143. pmx = pmx + mass_dat(region)%rm(is,js,ls, ipomcos) * modfrac(4) / volume
  144. pmx = pmx + mass_dat(region)%rm(is,js,ls, ipomaii) * modfrac(5) / volume
  145. pmx = pmx + mass_dat(region)%rm(is,js,ls, isoanus) * modfrac(1) / volume
  146. pmx = pmx + mass_dat(region)%rm(is,js,ls, isoaais) * modfrac(2) / volume
  147. pmx = pmx + mass_dat(region)%rm(is,js,ls, isoaacs) * modfrac(3) / volume
  148. pmx = pmx + mass_dat(region)%rm(is,js,ls, isoacos) * modfrac(4) / volume
  149. pmx = pmx + mass_dat(region)%rm(is,js,ls, isoaaii) * modfrac(5) / volume
  150. pmx = pmx + mass_dat(region)%rm(is,js,ls, issacs ) * modfrac(3) / volume
  151. pmx = pmx + mass_dat(region)%rm(is,js,ls, isscos ) * modfrac(4) / volume
  152. pmx = pmx + mass_dat(region)%rm(is,js,ls, iduacs ) * modfrac(3) / volume
  153. pmx = pmx + mass_dat(region)%rm(is,js,ls, iducos ) * modfrac(4) / volume
  154. pmx = pmx + mass_dat(region)%rm(is,js,ls, iduaci ) * modfrac(6) / volume
  155. pmx = pmx + mass_dat(region)%rm(is,js,ls, iducoi ) * modfrac(7) / volume
  156. pmx = pmx + mass_dat(region)%rm(is,js,ls, ino3_a ) * modfrac(3) / volume
  157. pmx = pmx + mass_dat(region)%rm(is,js,ls, inh4 ) * modfrac(3) / volume
  158. ! Count the water
  159. pmx = pmx + h2o_mode(region,1)%d3(is,js,ls)*modfrac(1) / volume
  160. pmx = pmx + h2o_mode(region,2)%d3(is,js,ls)*modfrac(2) / volume
  161. pmx = pmx + h2o_mode(region,3)%d3(is,js,ls)*modfrac(3) / volume
  162. pmx = pmx + h2o_mode(region,4)%d3(is,js,ls)*modfrac(4) / volume
  163. status = 0
  164. End Subroutine PMX_INTEGRATE_0D
  165. !EOC
  166. !--------------------------------------------------------------------------
  167. ! TM5 !
  168. !--------------------------------------------------------------------------
  169. !BOP
  170. !
  171. ! !IROUTINE: PMX_Integrate_3d
  172. !
  173. ! !DESCRIPTION: This routine uses PMX_Integrate_0d for generating a 3d
  174. ! array (over im, jm, lm) of pmx.
  175. ! Arbitrary radii may be specified.
  176. !\\
  177. !\\
  178. ! !INTERFACE:
  179. !
  180. SUBROUTINE PMX_Integrate_3d( region, isizelimit, pmx, status )
  181. !
  182. ! !USES:
  183. !
  184. USE DIMS, ONLY : im, jm, lm
  185. USE TM5_DISTGRID, ONLY : dgrid, Get_DistGrid, gather, scatter_i_band
  186. !
  187. ! !INPUT PARAMETERS:
  188. !
  189. INTEGER, INTENT(in) :: region
  190. REAL, INTENT(in) :: isizelimit
  191. !
  192. ! !OUTPUT PARAMETERS:
  193. !
  194. INTEGER, INTENT(out) :: status
  195. REAL, DIMENSION(:,:,:), INTENT(out) :: pmx
  196. !
  197. ! !REVISION HISTORY:
  198. ! 7 Oct 2010 - Achim Strunk - v0
  199. ! 19 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  200. !
  201. ! !REMARKS:
  202. !
  203. !EOP
  204. !------------------------------------------------------------------------
  205. !BOC
  206. CHARACTER(len=*), PARAMETER :: rname = mname//'/PMX_Integrate_3d'
  207. INTEGER :: is, js, ls, i1, i2, j1, j2
  208. REAL :: pmxl
  209. ! reset value
  210. pmx = 0.0
  211. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  212. DO ls = 1, lm(region)
  213. DO js = j1, j2
  214. DO is = i1, i2
  215. pmxl=0.0
  216. CALL pmx_integrate_0d( region, is, js, ls, isizelimit, pmxl, status )
  217. IF_NOTOK_RETURN( status=1 )
  218. !write(4444,*)pmxl,is,js,ls,ubound(pmx),lbound(pmx)
  219. pmx(is,js,ls) = pmxl
  220. !write(4444,*)pmx(is,js,ls),is,js,ls,ubound(pmx),lbound(pmx)
  221. END DO
  222. END DO
  223. ! write(5555,*)size(pmx),shape(pmx),i2,j2,lm(region),is,js,ls,isizelimit,pmxl
  224. END DO
  225. status = 0
  226. END SUBROUTINE PMX_Integrate_3d
  227. SUBROUTINE PMX_Integrate_2d( region, isizelimit, pmx, status )
  228. !
  229. ! !USES:
  230. !
  231. USE DIMS, ONLY : im, jm, lm
  232. USE TM5_DISTGRID, ONLY : dgrid, Get_DistGrid, gather, scatter_i_band
  233. !
  234. ! !INPUT PARAMETERS:
  235. !
  236. INTEGER, INTENT(in) :: region
  237. REAL, INTENT(in) :: isizelimit
  238. !
  239. ! !OUTPUT PARAMETERS:
  240. !
  241. INTEGER, INTENT(out) :: status
  242. REAL, DIMENSION(:,:), INTENT(out) :: pmx
  243. !
  244. ! !REVISION HISTORY:
  245. ! 7 Oct 2010 - Achim Strunk - v0
  246. ! 19 Jun 2012 - P. Le Sager - adapted for lon-lat MPI domain decomposition
  247. !
  248. ! !REMARKS:
  249. !
  250. !EOP
  251. !------------------------------------------------------------------------
  252. !BOC
  253. CHARACTER(len=*), PARAMETER :: rname = mname//'/PMX_Integrate_3d'
  254. INTEGER :: is, js, ls, i1, i2, j1, j2
  255. REAL :: pmxl
  256. ! reset value
  257. pmx = 0.0
  258. CALL GET_DISTGRID( dgrid(region), I_STRT=i1, I_STOP=i2, J_STRT=j1, J_STOP=j2 )
  259. DO js = j1, j2
  260. DO is = i1, i2
  261. CALL pmx_integrate_0d( region, is, js, ls, isizelimit, pmxl, status )
  262. IF_NOTOK_RETURN( status=1 )
  263. pmx(is,js) = pmxl
  264. END DO
  265. END DO
  266. status = 0
  267. END SUBROUTINE PMX_Integrate_2d
  268. !EOC
  269. END MODULE CALC_PM