calc_pm.F90 11 KB

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