phys_gph.F90 9.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382
  1. module Phys_GPH
  2. implicit none
  3. ! --- in/out ---------------------------------
  4. private
  5. public :: PotentialHeight
  6. public :: GeoPotentialHeightB
  7. public :: GeoPotentialHeight
  8. contains
  9. ! ===============================================
  10. !
  11. ! Compute height of air parcel between pressures p0 and p1:
  12. !
  13. ! h+dh +----+ p
  14. ! | |
  15. ! | | temperature constant within parcel
  16. ! | | specific humidity Q (kg h2o / kg air)
  17. ! | |
  18. ! h +----+ p+dp
  19. ! A = 1 m2
  20. !
  21. ! Constants:
  22. !
  23. ! g : gravity
  24. !
  25. ! eps1 = xm_air/xm_h2o - 1
  26. ! = 28.964/18.0 - 1.0 = 0.609
  27. !
  28. ! Mass of parcel:
  29. !
  30. ! dm = A dp/g [kg air]
  31. !
  32. ! Number of mol:
  33. !
  34. ! #mol = #mol_dry_air + #mol_h2o
  35. !
  36. ! = (kg dry air) / (kg dry air/mol dry air) + (kg h2o) / (kg h2o/mol h2o)
  37. !
  38. ! = dm (1-Q) / xm_air + dm Q / xm_h2o
  39. !
  40. ! = dm ( (1-Q)/xm_air + Q/xm_h2o )
  41. !
  42. ! = A dp/g ( 1 - Q + Q xm_air/xm_h2o ) / xm_air
  43. !
  44. ! = A dp/g ( 1 + (xm_air/xm_h2o - 1) Q ) / xm_air
  45. !
  46. ! = A dp/g ( 1 + eps1 Q ) / xm_air
  47. !
  48. ! General gas law:
  49. !
  50. ! p V = R T #mol
  51. !
  52. ! p A dh = R T A dp/g (1 + eps1 Q)/xm_air
  53. !
  54. ! dh = R/xm_air T(1+eQ)/p dp/g
  55. !
  56. ! Integrated height:
  57. !
  58. ! h = int dh = int R/xm_air T(1+eQ)/p dp/g
  59. !
  60. ! p_down
  61. ! = R/xm_air T(1+eQ) ( int 1/p dp ) / g
  62. ! p=p_up
  63. !
  64. ! p_down
  65. ! = R/xm_air T(1+eQ) [ ln(p) ]
  66. ! p=p_up
  67. !
  68. ! = R/M_air T(1+eQ) ln(p_down/p_up) / g
  69. !
  70. !
  71. ! Usefull initial value:
  72. !
  73. ! h0 = orography/g with orography in m*[g]=m2/s2
  74. !
  75. !
  76. ! Result: GeoPotentialHeightB GeoPotentialHeight
  77. !
  78. ! pb(n+1) +-----+ gph(n+1)
  79. ! | |
  80. ! T(n), Q(n) | | gph(n)
  81. ! | |
  82. ! pb(n) +-----+ gph(n)
  83. ! | |
  84. ! T(n-1), Q(n-1) | | gph(n-1)
  85. ! | |
  86. ! pb(n-1) +-----+ gph(n-1)
  87. ! | |
  88. ! : :
  89. ! | |
  90. ! pb(2) +-----+ gph(2)
  91. ! | |
  92. ! | | gph(1)
  93. ! | |
  94. ! pb(1) +-----+ gph(1) = h0
  95. !
  96. !
  97. subroutine PotentialHeight( ptop, pdown, T, Q, dh )
  98. use Binas, only : grav ! 9.80665 m/s2
  99. use Binas, only : xm_air ! 28.964 e-3 kg air/mol
  100. use Binas, only : xm_h2o ! 18.0 e-3 kg h2o/mol
  101. use Binas, only : Rgas ! 8.3144 J/mol/K
  102. ! --- in/out ---------------------------------------
  103. real, intent(in) :: ptop, pdown ! pressure bounds (Pa)
  104. real, intent(in) :: T ! temperature (K)
  105. real, intent(in) :: Q ! specific humidity (kg h2o/kg air)
  106. real, intent(out) :: dh ! geo.pot.height bounds (m)
  107. ! --- const ----------------------------------------
  108. ! eps1 = (kg air/mol) / (kg h2o/mol) - 1 = 0.609
  109. real, parameter :: eps1 = xm_air/xm_h2o - 1.0
  110. ! gas constant for dry air
  111. ! Rgas_air = Rgas / xm_air = 287.0598
  112. real, parameter :: Rgas_air = Rgas / xm_air
  113. ! --- begin -----------------------------------------
  114. ! geo potential height:
  115. dh = T * ( 1.0 + eps1 * Q ) * Rgas_air * log(pdown/ptop) / grav
  116. end subroutine PotentialHeight
  117. ! ***
  118. ! GPH at pressure half levels
  119. subroutine GeoPotentialHeightB( lm, pb, T, Q, h0, gphb )
  120. ! --- in/out ---------------------------------------
  121. integer, intent(in) :: lm
  122. real, intent(in) :: pb(lm+1) ! pressure bounds (Pa)
  123. real, intent(in) :: T(lm) ! temperature (K)
  124. real, intent(in) :: Q(lm) ! specific humidity (kg h2o/kg air)
  125. real, intent(in) :: h0 ! initial height (m)
  126. real, intent(out) :: gphb(lm+1) ! geo.pot.height bounds (m)
  127. ! --- local ----------------------------------------
  128. integer :: topb_gaslaw
  129. integer :: l
  130. logical :: inftop
  131. real :: dh
  132. ! --- begin -----------------------------------------
  133. if ( pb(1) >= pb(lm+1) ) then
  134. !
  135. ! surface -> top
  136. !
  137. ! zero top pressure ? then gph is infinite there ...
  138. inftop = pb(lm+1) < tiny(1.0)
  139. ! maximum boundary for which gph can be computed from general gas law:
  140. if ( inftop ) then
  141. topb_gaslaw = lm
  142. else
  143. topb_gaslaw = lm+1
  144. end if
  145. ! initial height at bottom:
  146. gphb(1) = h0 ! m
  147. ! loop over all boundaries for which gph can be computed
  148. ! from general gas law:
  149. do l = 2, topb_gaslaw
  150. ! potential height of this layer:
  151. call PotentialHeight( pb(l), pb(l-1), T(l-1), Q(l-1), dh )
  152. ! geo potential height at this boundary:
  153. gphb(l) = gphb(l-1) + dh
  154. end do
  155. ! infinte top ? set height to double of mid dh:
  156. if ( inftop ) then
  157. call PotentialHeight( 0.5*pb(lm), pb(lm), T(lm), Q(lm), dh )
  158. gphb(lm+1) = gphb(lm) + 2*dh
  159. end if
  160. else
  161. !
  162. ! top -> surface
  163. !
  164. ! zero top pressure ? then gph is infinite there ...
  165. inftop = pb(1) < tiny(1.0)
  166. ! maximum boundary for which gph can be computed from general gas law:
  167. if ( inftop ) then
  168. topb_gaslaw = 2
  169. else
  170. topb_gaslaw = 1
  171. end if
  172. ! initial height at bottom:
  173. gphb(lm+1) = h0 ! m
  174. ! loop over all boundaries for which gph can be computed
  175. ! from general gas law:
  176. do l = lm, topb_gaslaw, -1
  177. ! potential height of this layer:
  178. call PotentialHeight( pb(l), pb(l+1), T(l), Q(l), dh )
  179. ! geo potential height at this boundary:
  180. gphb(l) = gphb(l+1) + dh
  181. end do
  182. ! infinte top ? set height to double of mid dh:
  183. if ( inftop ) then
  184. call PotentialHeight( 0.5*pb(2), pb(2), T(1), Q(1), dh )
  185. gphb(1) = gphb(2) + 2*dh
  186. end if
  187. end if
  188. end subroutine GeoPotentialHeightB
  189. ! ***
  190. ! GPH at pressure full levels
  191. subroutine GeoPotentialHeight( lm, pb, T, Q, h0, gph )
  192. ! --- in/out ---------------------------------------
  193. integer, intent(in) :: lm
  194. real, intent(in) :: pb(lm+1) ! pressure bounds (Pa)
  195. real, intent(in) :: T(lm) ! temperature (K)
  196. real, intent(in) :: Q(lm) ! specific humidity (kg h2o/kg air)
  197. real, intent(in) :: h0 ! initial height (m)
  198. real, intent(out) :: gph(lm) ! geo.pot.height (m)
  199. ! --- local ----------------------------------------
  200. integer :: l
  201. real :: pf
  202. real :: dh
  203. ! --- begin -----------------------------------------
  204. if ( pb(1) >= pb(lm+1) ) then
  205. !
  206. ! surface -> top
  207. !
  208. ! initial height at bottom:
  209. pf = ( pb(1) + pb(2) )/2.0
  210. call PotentialHeight( pf, pb(1), T(1), Q(1), dh ) ! m
  211. gph(1) = h0 + dh
  212. ! loop over remaining full levels: 2,lm,1 or lm-1,1,-1
  213. do l = 2, lm
  214. ! potential height of second half previous layer:
  215. call PotentialHeight( pb(l), pf, T(l-1), Q(l-1), dh )
  216. gph(l) = gph(l-1) + dh
  217. ! new full pressure:
  218. pf = ( pb(l) + pb(l+1) )/2.0
  219. ! potential height of first half of current layer:
  220. call PotentialHeight( pf, pb(l), T(l), Q(l), dh )
  221. gph(l) = gph(l) + dh
  222. end do
  223. else
  224. !
  225. ! top -> surface
  226. !
  227. ! initial height at bottom:
  228. pf = ( pb(lm+1) + pb(lm) )/2.0
  229. call PotentialHeight( pf, pb(lm+1), T(lm), Q(lm), dh ) ! m
  230. gph(lm) = h0 + dh
  231. ! loop over remaining full levels: 2,lm,1 or lm-1,1,-1
  232. do l = lm-1, 1, -1
  233. ! potential height of second half previous layer:
  234. call PotentialHeight( pb(l+1), pf, T(l+1), Q(l+1), dh )
  235. gph(l) = gph(l+1) + dh
  236. ! new full pressure:
  237. pf = ( pb(l+1) + pb(l) )/2.0
  238. ! potential height of first half of current layer:
  239. call PotentialHeight( pf, pb(l+1), T(l), Q(l), dh )
  240. gph(l) = gph(l) + dh
  241. end do
  242. end if
  243. end subroutine GeoPotentialHeight
  244. end module Phys_GPH
  245. ! ################################################################
  246. !
  247. !
  248. !program test
  249. !
  250. ! use phys_gph
  251. !
  252. ! integer, parameter :: lm = 4
  253. !
  254. ! real :: pb(lm+1)
  255. ! real :: T(lm), Q(lm)
  256. ! real :: gph(lm), gphb(lm+1)
  257. !
  258. ! pb = (/ 1000.0, 500.0, 100.0, 50.0, 0.0 /)
  259. ! T = (/ 293.0, 273.0, 260.0, 280.0 /)
  260. ! Q = 0.0
  261. !
  262. ! print *, ' '
  263. ! print *, 'gph u'
  264. ! call GeoPotentialHeight( lm, pb, T, Q, h0, gph )
  265. ! print *, ' '
  266. ! print *, 'gph u'
  267. ! call GeoPotentialHeightB( lm, pb, T, Q, h0, gphb )
  268. !
  269. ! print *, ' '
  270. ! do l = 1, lm
  271. ! print *, gphb(l)
  272. ! print *, ' ', gph(l)
  273. ! end do
  274. ! print *, gphb(lm+1)
  275. !
  276. !
  277. ! pb = (/ 0.0, 50.0, 100.0, 500.0, 1000.0 /)
  278. ! T = (/ 280.0, 260.0, 273.0, 293.0 /)
  279. ! Q = 0.0
  280. !
  281. ! print *, ' '
  282. ! print *, 'gph d'
  283. ! call GeoPotentialHeight( lm, pb, T, Q, h0, gph )
  284. ! print *, ' '
  285. ! print *, 'gphb d'
  286. ! call GeoPotentialHeightB( lm, pb, T, Q, h0, gphb )
  287. !
  288. ! print *, ' '
  289. ! do l = 1, lm
  290. ! print *, gphb(l)
  291. ! print *, ' ', gph(l)
  292. ! end do
  293. ! print *, gphb(lm+1)
  294. !
  295. !
  296. !end program test
  297. !