grid_3d.F90 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378
  1. !
  2. ! 3d grid transforms
  3. !
  4. !
  5. ! call Fill3D( lli, levi, ps, field, &
  6. ! lliX, leviX, fieldX, &
  7. ! combkey, status )
  8. !
  9. ! o lli,levi : output field defintions (in)
  10. ! ps : surface pressure in hPa (in)
  11. ! field : output field (out)
  12. !
  13. ! o lliX, leviX : input field definitons (in)
  14. ! fieldX : input field (in)
  15. !
  16. ! o combkey : 'mass-aver', 'sum', 'aver' (in)
  17. !
  18. ! o note that psX is not required ...
  19. !
  20. ! call FillMass ( m , lli, levi, sp , status )
  21. !
  22. ! ! Fill 3D mass array given surface pressure and grid cell aera's.
  23. !
  24. ! call FillMassChange( dm_dt, lli, levi, sp_t1, sp_t2, status )
  25. !
  26. ! ! Fill 3D mass change between two times given different
  27. ! ! surface pressures;
  28. ! ! note that ( a + b sp2 ) - ( a + b sp1 )
  29. ! ! is not the same as a + b (sp2 - sp1) ...
  30. !
  31. !
  32. module grid_3d
  33. implicit none
  34. ! --- in/out -----------------------------
  35. private
  36. public :: Fill3D
  37. public :: Regrid3D
  38. public :: FillMass
  39. public :: FillMassChange
  40. ! --- const ---------------------------------
  41. character(len=*), parameter :: mname = 'grid_3d'
  42. contains
  43. !--------------------------------------------------------------------------
  44. ! TM5 !
  45. !--------------------------------------------------------------------------
  46. !BOP
  47. !
  48. ! !IROUTINE: Regrid3D
  49. !
  50. ! !DESCRIPTION:
  51. ! Performs vertical and horizontal regridding of 3D data sets. This is a
  52. ! convenient wrapper around FillGrid and FillLevels.
  53. !
  54. ! Available methods:
  55. ! for horizontal regridding : 'sum', 'aver', 'area-aver', 'weight'
  56. ! for vertical regridding : 'sum', 'aver', 'mass-aver', 'top', 'bottom'
  57. !
  58. ! See FillGrid and FillLevels for detailed documentation.
  59. !
  60. !\\
  61. !\\
  62. ! !INTERFACE:
  63. !
  64. subroutine Regrid3D( outHInfo, outVInfo, OutData, &
  65. inHInfo, inVInfo, InData, &
  66. hcomb, vcomb, sp, &
  67. nuv, nw, status )
  68. !
  69. ! !USES:
  70. !
  71. use grid_type_ll , only : TllGridInfo, FillGrid
  72. use grid_type_hyb, only : TLevelInfo, FillLevels
  73. !
  74. ! !INPUT PARAMETERS:
  75. !
  76. type(TllGridInfo), intent(in) :: outHInfo ! Output grid : horizontal info (LonLat)
  77. type(TLevelInfo), intent(in) :: outVInfo ! Output grid : vertical info
  78. character(len=*), intent(in) :: nuv ! Horizontal location ('n'=center, 'u'=west edge, 'v'=north edge)
  79. character(len=*), intent(in) :: nw ! Vertical location ('n'=center, 'w'=edge)
  80. real, intent(in) :: sp(:,:) ! Surface pressure in Pa
  81. type(TllGridInfo), intent(in) :: inHInfo ! Input grid : horizontal info
  82. type(TLevelInfo), intent(in) :: inVInfo ! Input grid : vertical info
  83. real, intent(in) :: InData(:,:,:) ! Input data
  84. character(len=*), intent(in) :: hcomb ! combination key for horizontal regridding
  85. character(len=*), intent(in) :: vcomb ! combination key for vertical regridding
  86. !
  87. ! !OUTPUT PARAMETERS:
  88. !
  89. real, intent(out) :: OutData(:,:,:) ! output data
  90. integer, intent(out) :: status ! return code
  91. !
  92. ! !REVISION HISTORY:
  93. ! 15 Dec 2010 - P. Le Sager - written, based on Fill3D.
  94. !
  95. ! !REMARKS:
  96. ! (1) Note that surface pressure on the input grid is not needed.
  97. ! (2) For now, only data on horizontal grid centers can be regriddded.
  98. ! So, you must use 'n' for NUV (pls, 15-12-2010)
  99. !
  100. !EOP
  101. !------------------------------------------------------------------------
  102. !BOC
  103. ! --- local -------------------------
  104. character(len=*), parameter :: name = mname//'/Regrid3D'
  105. real, allocatable :: field_ll(:,:,:)
  106. integer :: l
  107. ! init
  108. ! -----------------
  109. ! output horizontal grid, input levels grid
  110. allocate( field_ll( outHInfo%nlon, outHInfo%nlat, inVInfo%nlev ) )
  111. ! horizontal
  112. ! -----------------
  113. do l = 1, inVInfo%nlev
  114. call FillGrid( outHInfo , nuv, field_ll(:,:,l), &
  115. inHInfo, nuv, InData(:,:,l), hcomb, status )
  116. if (status<0) then
  117. write (*,'("ERROR - only part of target grid filled")')
  118. write (*,'("ERROR in ",a)') name; status=1; return
  119. endif
  120. if (status/=0) then
  121. write (*,'("ERROR in ",a)') name; status=1
  122. return
  123. endif
  124. end do
  125. ! vertical
  126. ! -----------------
  127. call FillLevels( outVInfo, nw, sp, OutData, inVInfo, field_ll, vcomb, status )
  128. if (status/=0) then; write (*,'("ERROR in ",a)') name; status=1; return
  129. end if
  130. ! done
  131. ! -----------------
  132. deallocate( field_ll )
  133. status = 0
  134. end subroutine Regrid3D
  135. !EOC
  136. ! =========================================
  137. subroutine Fill3D( lli, levi, nw, ps, field, &
  138. lliX, leviX, fieldX, &
  139. combkey, status )
  140. use grid_type_ll , only : TllGridInfo, FillGrid
  141. use grid_type_hyb, only : TLevelInfo, FillLevels
  142. ! --- in/out --------------------------------
  143. type(TllGridInfo), intent(in) :: lli
  144. type(TLevelInfo), intent(in) :: levi
  145. character(len=*), intent(in) :: nw
  146. real, intent(in) :: ps(:,:) ! Pa
  147. real, intent(out) :: field(:,:,:)
  148. type(TllGridInfo), intent(in) :: lliX
  149. type(TLevelInfo), intent(in) :: leviX
  150. real, intent(in) :: fieldX(:,:,:)
  151. character(len=*), intent(in) :: combkey
  152. integer, intent(out) :: status
  153. ! --- const ---------------------------------
  154. character(len=*), parameter :: name = mname//'/Fill3D'
  155. ! --- local -------------------------
  156. real, allocatable :: field_ll(:,:,:)
  157. integer :: l
  158. ! --- begin -------------------------
  159. ! output horizontal grid, input levels
  160. allocate( field_ll(lli%nlon,lli%nlat,leviX%nlev) )
  161. select case ( combkey )
  162. !
  163. ! mass average
  164. !
  165. case ( 'mass-aver' )
  166. ! horizontal
  167. do l = 1, leviX%nlev
  168. call FillGrid( lli , 'n', field_ll(:,:,l), &
  169. lliX, 'n', fieldX(:,:,l), 'area-aver', status )
  170. if (status<0) then
  171. write (*,'("ERROR - only part of target grid filled")')
  172. write (*,'("ERROR in ",a)') name; status=1; return
  173. end if
  174. if (status/=0) then; write (*,'("ERROR in ",a)') name; status=1; return; end if
  175. end do
  176. ! vertical
  177. call FillLevels( levi, nw, ps, field, leviX, field_ll, 'mass-aver', status )
  178. if (status/=0) then; write (*,'("ERROR in ",a)') name; status=1; return; end if
  179. !
  180. ! other (should be supported by FillGrid and FillLevels)
  181. !
  182. case default
  183. ! horizontal
  184. do l = 1, leviX%nlev
  185. call FillGrid( lli , 'n', field_ll(:,:,l), &
  186. lliX, 'n', fieldX(:,:,l), combkey, status )
  187. if (status<0) then
  188. write (*,'("ERROR - only part of target grid filled")')
  189. write (*,'("ERROR in ",a)') name; status=1; return
  190. end if
  191. if (status/=0) then; write (*,'("ERROR in ",a)') name; status=1; return; end if
  192. end do
  193. ! vertical
  194. call FillLevels( levi, nw, ps, field, leviX, field_ll, combkey, status )
  195. if (status/=0) then; write (*,'("ERROR in ",a)') name; status=1; return; end if
  196. end select
  197. ! done
  198. deallocate( field_ll )
  199. ! ok
  200. status = 0
  201. end subroutine Fill3D
  202. ! **************************************************************
  203. !
  204. ! p = a + b sp
  205. !
  206. subroutine FillMass( m, lli, levi, sp, status )
  207. use Binas , only : grav
  208. use grid_type_ll , only : TllGridInfo, AreaOper
  209. use grid_type_hyb, only : TLevelInfo
  210. ! --- begin ---------------------------------
  211. real, intent(out) :: m(:,:,:) ! kg
  212. type(TllGridInfo), intent(in) :: lli
  213. type(TLevelInfo), intent(in) :: levi
  214. real, intent(in) :: sp(:,:) ! Pa
  215. integer, intent(out) :: status
  216. ! --- const ---------------------------------
  217. character(len=*), parameter :: rname = mname//'/FillMass'
  218. ! --- local -------------------------
  219. integer :: l
  220. ! --- begin -------------------------
  221. ! check shape of target grid:
  222. if ( (size(m,1) /= lli%nlon ) .or. (size(m,2) /= lli%nlat) .or. &
  223. (size(m,3) /= levi%nlev) ) then
  224. write (*,'("ERROR - target array does not match with grid definition:")')
  225. write (*,'("ERROR - lli : ",i3," x ",i3 )') lli%nlon, lli%nlat
  226. write (*,'("ERROR - levi : ",i3 )') levi%nlev
  227. write (*,'("ERROR - ll : ",i3," x ",i3," x ",i3)') shape(m)
  228. write (*,'("ERROR in ",a)') rname; status=1; return
  229. end if
  230. ! Pa = kg g / A -> kg = A * Pa/g
  231. ! loop over levels
  232. do l = 1, levi%nlev
  233. m(:,:,l) = levi%da(l) + levi%db(l) * sp / grav ! Pa/g = kg/m2
  234. call AreaOper( lli, m(:,:,l), '*', 'm2', status ) ! kg
  235. if (status/=0) then; write (*,'("ERROR in ",a)') rname; status=1; return; end if
  236. end do
  237. ! ok
  238. status = 0
  239. end subroutine FillMass
  240. ! ***
  241. !
  242. ! p1 = a + b sp1
  243. ! p2 = a + b sp2
  244. !
  245. ! p2 - p1 = b (sp2 - sp1)
  246. !
  247. ! m = (p2 - p1)/g * A
  248. !
  249. subroutine FillMassChange( dm, lli, levi, sp1, sp2, status )
  250. use Binas , only : grav
  251. use grid_type_ll , only : TllGridInfo, AreaOper
  252. use grid_type_hyb, only : TLevelInfo
  253. ! --- begin ---------------------------------
  254. real, intent(out) :: dm(:,:,:) ! kg
  255. type(TllGridInfo), intent(in) :: lli
  256. type(TLevelInfo), intent(in) :: levi
  257. real, intent(in) :: sp1(:,:) ! Pa
  258. real, intent(in) :: sp2(:,:) ! Pa
  259. integer, intent(out) :: status
  260. ! --- const ---------------------------------
  261. character(len=*), parameter :: rname = mname//'/FillMassChange'
  262. ! --- local -------------------------
  263. integer :: l
  264. ! --- begin -------------------------
  265. ! check shape of target grid:
  266. if ( (size(dm,1) /= lli%nlon ) .or. (size(dm,2) /= lli%nlat) .or. &
  267. (size(dm,3) /= levi%nlev) ) then
  268. write (*,'("ERROR - target array does not match with grid definition:")')
  269. write (*,'("ERROR - lli : ",i3," x ",i3 )') lli%nlon, lli%nlat
  270. write (*,'("ERROR - levi : ",i3 )') levi%nlev
  271. write (*,'("ERROR - ll : ",i3," x ",i3," x ",i3)') shape(dm)
  272. write (*,'("ERROR in ",a)') rname; status=1; return
  273. end if
  274. ! Pa = kg g / A -> kg = A * Pa/g
  275. ! loop over levels
  276. do l = 1, levi%nlev
  277. dm(:,:,l) = abs(levi%db(l)) * ( sp2 - sp1 ) / grav ! Pa/g = kg/m2
  278. call AreaOper( lli, dm(:,:,l), '*', 'm2', status ) ! kg
  279. if (status/=0) then; write (*,'("ERROR in ",a)') rname; status=1; return; end if
  280. end do
  281. ! ok
  282. status = 0
  283. end subroutine FillMassChange
  284. end module grid_3d