MGridTools.F90 6.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182
  1. ! First include the set of model-wide compiler flags
  2. #include "tm5.inc"
  3. module MGridTools
  4. ! Tool subroutines related to the horizontal and vertical grid
  5. implicit none
  6. private
  7. ! subroutine GetTMCellIndex( xlon,xlat,ix,iy,dx,dy )
  8. ! subroutine GetTMCellIndex4(xlon,xlat,ix,ix2,iy,iy2,gdx,gdy)
  9. public :: getTMCellIndex, getTMCellIndex4
  10. contains
  11. subroutine GetTMCellIndex( im, jm, xlon, xlat, ix, iy, dx, dy )
  12. !------------------------------------------------------------------------
  13. ! Determine TM5 cell indices (ix,iy) and the relative position inside
  14. ! the cell (dx,dy) from the cell centre of the location (xlon,xlat).
  15. !
  16. ! im, jm : TM array dimensions [in]
  17. ! xlon,xlat : longitude, latitude in degrees [in]
  18. ! ix,iy : indices of the TM5 grid cell containing (xlon,xlat) [out]
  19. ! dx,dy : relative position with respect to the cell centre,
  20. ! el. (-0.5,0.5) [out]
  21. !------------------------------------------------------------------------
  22. use physconstants, only : pi
  23. implicit none
  24. integer,intent(in) :: im, jm
  25. real,intent(in) :: xlon, xlat
  26. integer,intent(out) :: ix, iy
  27. real,intent(out) :: dx, dy
  28. !
  29. real :: dlon,dlat,dlondeg,dlatdeg,phim1,tetam1
  30. real :: rx,ry,todegree,firstcellmid,xlon2
  31. !
  32. todegree=180.0/pi
  33. !
  34. xlon2=xlon
  35. if( xlon < -180.0 ) xlon2=xlon+360.0
  36. if( xlon >= 180.0 ) xlon2=xlon-360.0 ! xlon2 in [-180,180)
  37. !
  38. if ( im == 0 ) then
  39. print*,'ERROR in ass_gridtools.f90; im = 0 '
  40. stop
  41. end if
  42. !
  43. dlon = 2.0*pi/real(im)
  44. dlondeg = todegree*dlon ! deg
  45. phim1 = 0.5*dlon - pi ! lon first grid cell
  46. firstcellmid = todegree*phim1 ! deg
  47. rx = (xlon2-firstcellmid)/dlondeg + 1.0 ! old grid: [+1.0,im+1.0)
  48. ! new grid: [+0.5,im+0.5)
  49. !old dlon=360.0/real(im)
  50. !old rx = xlon/dlon+0.5*im+1.0
  51. ix = nint(rx) ! old: [1,im+1] new: [1,im]
  52. dx = rx - real(ix)
  53. if ( ix == im+1 ) ix = 1
  54. if ( ix == 0 ) ix = im ! for safety ...
  55. !
  56. if ( jm == 0 ) then
  57. print*,'ERROR in ass_gridtools.f90; jm = 0 '
  58. stop
  59. end if
  60. !
  61. dlat = pi/real(jm)
  62. dlatdeg = todegree*dlat ! deg
  63. tetam1 = 0.5*dlat-0.5*pi ! lat first grid cell
  64. firstcellmid = todegree*tetam1 ! deg
  65. ry = (xlat-firstcellmid)/dlatdeg+1.0 ! [0.5,jm+0.5]
  66. iy = nint(ry) ! [0,jm+1] for rounding error
  67. if ( iy < 1 ) iy = 1 ! [1,jm+1]
  68. if ( iy > jm ) iy = jm ! [1,jm]
  69. !old ry = 0.5+(xlat+90.0)*real(jm)/180.0
  70. !old iy = 1+int((xlat+90.0)*real(jm)/180.0)
  71. dy = ry - real(iy)
  72. !
  73. end subroutine GetTMCellIndex
  74. subroutine GetTMCellIndex4 ( im, jm, xlon, xlat, ix4a, iy4a, w4a )
  75. !------------------------------------------------------------------------
  76. !
  77. ! Find the four grid points that
  78. ! are nearest to the position (lon,lat) on the globe.
  79. !
  80. ! im, jm : TM array dimensions [in]
  81. ! (xlon,xlat) : position on the globe, in degrees (in)
  82. ! ix4a,iy4a : east-north indices of the four TM5 grid points (out)
  83. ! w4a : weight factors of the four points (out)
  84. ! Henk Eskes, KNMI
  85. !------------------------------------------------------------------------
  86. use physconstants, only : pi
  87. implicit none
  88. ! in/out
  89. integer,intent(in) :: im, jm
  90. real,intent(in) :: xlon, xlat
  91. integer,dimension(4),intent(out) :: ix4a, iy4a
  92. real,dimension(4),intent(out) :: w4a
  93. ! local
  94. integer :: ix,ix2,iy,iy2
  95. real :: gdx,gdy
  96. real :: dlon,dlat,xlon2,phim1,tetam1
  97. real :: rx,ry,todegree,firstcellmid,dlondeg,dlatdeg
  98. !
  99. ! determine south-west model cell indices (ix,iy)
  100. ! and inside-cell cordinate (gdx,gdy), both el. [0,1]
  101. ! corresponding to the pixel coordinates (longitude=xx,latitude=yy)
  102. ! ( lon el. [-180,180] )
  103. !
  104. todegree=180.0/pi
  105. !
  106. xlon2=xlon
  107. if( xlon < -180.0 ) xlon2=xlon+360.0
  108. if( xlon >= 180.0 ) xlon2=xlon-360.0 ! xlon2 in [-180,180)
  109. !
  110. dlon = 2.0*pi/real(im)
  111. dlondeg = todegree*dlon ! deg
  112. phim1 = 0.5*dlon-pi ! lon first grid cell, in rad
  113. firstcellmid = todegree*phim1 ! deg
  114. rx = (xlon2-firstcellmid)/dlondeg + 1.0 ! old grid: [+1.0,im+1.0)
  115. ! new grid: [+0.5,im+0.5)
  116. ix = int(rx)
  117. gdx = rx - real(ix)
  118. if ( ix == im+1 ) ix = 1
  119. if ( ix == 0 ) ix = im ! for safety ...
  120. !
  121. dlat = pi/real(jm)
  122. dlatdeg = todegree*dlat ! deg
  123. tetam1 = 0.5*dlat - 0.5*pi ! lat first grid cell
  124. firstcellmid = todegree*tetam1 ! deg
  125. ry = (xlat-firstcellmid)/dlatdeg + 1.0 ! [0.5,jm+0.5]
  126. iy = int(ry)
  127. gdy = ry - real(iy)
  128. !
  129. ! determine coordinates of other corner points (ix,iy2),(ix2,iy),(ix2,iy2)
  130. ix2 = ix + 1
  131. if( ix2 == im + 1 ) then ! periodic east-west boundary
  132. ix2 = 1
  133. end if
  134. iy2 = iy + 1
  135. if( iy >= jm ) then ! north pole
  136. iy = jm - 1
  137. iy2 = jm
  138. gdy = 1.0
  139. end if
  140. if( iy <= 0 ) then ! south pole
  141. iy = 1
  142. iy2 = 2
  143. gdy = 0.0
  144. end if
  145. !
  146. ix4a(:) = (/ ix, ix2, ix, ix2 /)
  147. iy4a(:) = (/ iy, iy, iy2, iy2 /)
  148. !if ( (iy == 1) .or. (iy == jm) ) then ! polar cap exception !JDMP
  149. ! ix4a(1) = 1 !JDMP
  150. ! ix4a(2) = 1 !JDMP
  151. !end if !JDMP
  152. !if ( (iy2 == 1) .or. (iy2 == jm) ) then ! polar cap exception !JDMP
  153. ! ix4a(3) = 1 !JDMP
  154. ! ix4a(4) = 1 !JDMP
  155. !end if !JDMP
  156. w4a(1) = (1.0-gdx)*(1.0-gdy)
  157. w4a(2) = gdx *(1.0-gdy)
  158. w4a(3) = (1.0-gdx)* gdy
  159. w4a(4) = gdx * gdy
  160. !
  161. end subroutine GetTMCellIndex4
  162. end module MGridTools