mod_sphere_tools.F90 9.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342
  1. module mod_sphere_tools
  2. contains
  3. ! Routine to get cartesian coordinates from geographical coordinates
  4. function geo2cart(lon,lat)
  5. real, parameter :: rad=1.7453292519943295E-02,deg=57.29577951308232
  6. real, intent(in) :: lon,lat
  7. real, dimension(3) :: geo2cart
  8. real :: lambda
  9. lambda=lat*rad
  10. theta=lon*rad
  11. geo2cart(1)=cos(lambda)*cos(theta)
  12. geo2cart(2)=cos(lambda)*sin(theta)
  13. geo2cart(3)=sin(lambda)
  14. end function geo2cart
  15. ! Routine to calculate cross product of two 3D vectors.
  16. function cross_product(v1,v2)
  17. implicit none
  18. real, intent(in), dimension(3) :: v1,v2
  19. real, dimension(3) :: cross_product
  20. cross_product(1) = v1(2)*v2(3) - v1(3)*v2(2)
  21. cross_product(2) = v1(3)*v2(1) - v1(1)*v2(3)
  22. cross_product(3) = v1(1)*v2(2) - v1(2)*v2(1)
  23. end function cross_product
  24. ! Routine to calculate vector norm (more precise the 2-norm)
  25. function norm2(vector,p)
  26. implicit none
  27. integer, intent(in) :: p
  28. real, intent(in) :: vector(p)
  29. integer :: i
  30. real :: norm2
  31. norm2=0.
  32. do i=1,p
  33. norm2 = norm2 + vector(i)**2
  34. end do
  35. norm2=sqrt(norm2)
  36. end function norm2
  37. ! Routine to calculate wether the point (plon,plat) is in the box
  38. ! defined by crnlon,crnlat. Cnrlon/crnlat must be traversed so that
  39. ! they form a convex polygon in 3D coordinates when following indices.
  40. ! It should work for for all regions defined by crnlon/crnlat...
  41. function inbox(crnlon,crnlat,npt,plon,plat)
  42. implicit none
  43. real, parameter :: rad=1.7453292519943295E-02,deg=57.29577951308232
  44. integer, intent(in) :: npt
  45. real, dimension(npt), intent(in) :: crnlon,crnlat
  46. real, intent(in) :: plon,plat
  47. logical :: inbox
  48. real, dimension(npt,3) :: cvec
  49. real, dimension(3) :: pvec
  50. real, dimension(3,3) :: rvec
  51. real, dimension(3) :: nvec, nvec_prev, cprod
  52. integer :: i,im1,ip1
  53. logical :: lsign
  54. real :: rotsign,old_rotsign
  55. ! point vector from origo
  56. pvec = geo2cart(plon,plat)
  57. !print *,crnlon
  58. !print *,crnlat
  59. ! vector to rectangle corner
  60. do i=1,npt
  61. cvec(i,:) = geo2cart(crnlon(i),crnlat(i))
  62. end do
  63. ! Traverse box boundaries -- Check that traversion is
  64. ! consistent and that point is in box
  65. lsign=.true.
  66. i=1
  67. old_rotsign=0.
  68. do while (i<npt+1 .and. lsign)
  69. ip1= mod(i ,npt)+1
  70. im1= mod(i-2+npt,npt)+1
  71. ! Vectors used to span planes
  72. rvec(3,:) = cvec(ip1,:)
  73. rvec(2,:) = cvec(i ,:)
  74. rvec(1,:) = cvec(im1,:)
  75. ! Normal vector to two spanning planes
  76. nvec = cross_product(rvec(2,:),rvec(3,:))
  77. nvec_prev = cross_product(rvec(1,:),rvec(2,:))
  78. ! As we move to new planes, the cross product rotates in
  79. ! a certain direction
  80. cprod = cross_product(nvec_prev,nvec)
  81. ! For anticlockwise rotation, this should be positive
  82. rotsign=sign(1.,dot_product(cprod,rvec(2,:)))
  83. ! Check that box is consistently traversed
  84. if (i>1 .and. rotsign * old_rotsign < 0) then
  85. print *,'Grid cell not consistently traversed'
  86. print *,'or polygon is not convex'
  87. stop '(inbox2)'
  88. end if
  89. old_rotsign=rotsign
  90. ! If this is true for all four planes, we are in grid box
  91. lsign = lsign .and. (dot_product(nvec,pvec)*rotsign)>0
  92. i=i+1
  93. end do
  94. inbox = lsign
  95. end function inbox
  96. ! Routine to get angle between two vectors defined in
  97. ! geographical coordinates.
  98. function secangle(lon1,lat1,lon2,lat2)
  99. implicit none
  100. real, intent(in), dimension(2) :: lon1,lat1,lon2,lat2
  101. real, dimension(3) :: nx, n2, ny
  102. real :: cos1, cos2
  103. real :: secangle
  104. ! Normal of the planes defined by positions and origo
  105. nx = cross_product(geo2cart(lon1(1),lat1(1)),geo2cart(lon1(2),lat1(2)))
  106. n2 = cross_product(geo2cart(lon2(1),lat2(1)),geo2cart(lon2(2),lat2(2)))
  107. ! Normal to position 1 and vector 1 (x) -- forms rh system
  108. ny = cross_product(geo2cart(lon1(1),lat1(1)), nx)
  109. ny = ny / norm2(ny,3)
  110. ! Angle info 1 -- cosine of angle between planes nx and n2
  111. cos1 = dot_product(nx,n2)
  112. ! Angle info 2 -- Cosine of angle between planes ny and n2
  113. cos2 = dot_product(ny,n2)
  114. ! Angle between vectors 1 and 2
  115. secangle = atan2(cos1,cos2)
  116. end function secangle
  117. ! Intersection routine by Mats Bentsen.
  118. ! --- this routine computes the lat/lon coordinates for the intersection
  119. ! --- of the two geodesic lines which connects the lat/lon pairs a1,a2
  120. ! --- and b1,b2
  121. logical function intersect(lat_a1,lon_a1,lat_a2,lon_a2, &
  122. lat_b1,lon_b1,lat_b2,lon_b2, &
  123. lat_i,lon_i)
  124. implicit none
  125. real lat_a1,lon_a1,lat_a2,lon_a2, &
  126. lat_b1,lon_b1,lat_b2,lon_b2, &
  127. lat_i,lon_i
  128. real lambda,theta, &
  129. x_a1,y_a1,z_a1,x_a2,y_a2,z_a2, &
  130. x_b1,y_b1,z_b1,x_b2,y_b2,z_b2, &
  131. x_na,y_na,z_na,x_nb,y_nb,z_nb, &
  132. x_i,y_i,z_i,l_i, &
  133. x_a,y_a,z_a,x_b,y_b,z_b,l_a,l_b,l_l,a_a,a_b
  134. real rad,deg
  135. parameter(rad=1.7453292519943295E-02,deg=57.29577951308232)
  136. ! --- transforming from spherical to cartesian coordinates
  137. lambda=lat_a1*rad
  138. theta=lon_a1*rad
  139. x_a1=cos(lambda)*cos(theta)
  140. y_a1=cos(lambda)*sin(theta)
  141. z_a1=sin(lambda)
  142. lambda=lat_a2*rad
  143. theta=lon_a2*rad
  144. x_a2=cos(lambda)*cos(theta)
  145. y_a2=cos(lambda)*sin(theta)
  146. z_a2=sin(lambda)
  147. lambda=lat_b1*rad
  148. theta=lon_b1*rad
  149. x_b1=cos(lambda)*cos(theta)
  150. y_b1=cos(lambda)*sin(theta)
  151. z_b1=sin(lambda)
  152. lambda=lat_b2*rad
  153. theta=lon_b2*rad
  154. x_b2=cos(lambda)*cos(theta)
  155. y_b2=cos(lambda)*sin(theta)
  156. z_b2=sin(lambda)
  157. x_na=y_a1*z_a2-y_a2*z_a1
  158. y_na=z_a1*x_a2-z_a2*x_a1
  159. z_na=x_a1*y_a2-x_a2*y_a1
  160. x_nb=y_b1*z_b2-y_b2*z_b1
  161. y_nb=z_b1*x_b2-z_b2*x_b1
  162. z_nb=x_b1*y_b2-x_b2*y_b1
  163. ! --- Let a1 be the vector from the center of the sphere to the point
  164. ! --- (lat_a1,lon_a1) on the sphere. Similar with vectors a2, b1 and b2.
  165. ! --- Then we compute the components and length of a vector i pointing
  166. ! --- along the intersection of the two planes spanned out by the
  167. ! --- vectors a1, a2 and b1, b2 respectively.
  168. x_i=y_na*z_nb-y_nb*z_na
  169. y_i=z_na*x_nb-z_nb*x_na
  170. z_i=x_na*y_nb-x_nb*y_na
  171. l_i=sqrt(x_i*x_i+y_i*y_i+z_i*z_i)
  172. !
  173. ! --- check if i lies between a1 and a2
  174. !
  175. intersect=.true.
  176. !
  177. ! --- first find the vector a between a1 and a2 and its angle a_a to a1
  178. ! --- and a2
  179. !
  180. x_a=x_a1+x_a2
  181. y_a=y_a1+y_a2
  182. z_a=z_a1+z_a2
  183. l_a=sqrt(x_a*x_a+y_a*y_a+z_a*z_a)
  184. l_l=sign(1.,l_i*l_a)*max(1.e-9,abs(l_i*l_a))
  185. a_a=acos(max(-1.,min(1.,x_a1*x_a2+y_a1*y_a2+z_a1*z_a2)))*0.5
  186. ! --- if the angle between i and a is greater than
  187. ! --- a_a, then intersect=.false.
  188. if (acos(max(-1.,min(1.,(x_i*x_a+y_i*y_a+z_i*z_a)/l_l))) &
  189. .gt.a_a) then
  190. ! --- - test the opposite directed intersection vector
  191. x_i=-x_i
  192. y_i=-y_i
  193. z_i=-z_i
  194. if (acos(max(-1.,min(1.,(x_i*x_a+y_i*y_a+z_i*z_a)/l_l))) &
  195. .gt.a_a) intersect=.false.
  196. endif
  197. ! do similar test for b1 and b2
  198. if (intersect) then
  199. x_b=x_b1+x_b2
  200. y_b=y_b1+y_b2
  201. z_b=z_b1+z_b2
  202. l_b=sqrt(x_b*x_b+y_b*y_b+z_b*z_b)
  203. l_l=sign(1.,l_i*l_b)*max(1.e-9,abs(l_i*l_b))
  204. a_b=acos(max(-1.,min(1.,x_b1*x_b2+y_b1*y_b2+z_b1*z_b2)))*0.5
  205. if (acos(max(-1.,min(1.,(x_i*x_b+y_i*y_b+z_i*z_b)/l_l))) &
  206. .gt.a_b) intersect=.false.
  207. endif
  208. ! represent the intersection in lat,lon coordinates
  209. lat_i=atan2(z_i,sqrt(x_i*x_i+y_i*y_i))*deg
  210. lon_i=atan2(y_i,x_i)*deg
  211. end function
  212. elemental real function spherdist(lon1,lat1,lon2,lat2)
  213. ! --- -----------------------------------------
  214. ! --- Computes the distance between geo. pos.
  215. ! --- lon1,lat1 and lon2,lat2.
  216. ! --- INPUT is in degrees.
  217. ! --- -----------------------------------------
  218. implicit none
  219. REAL, intent(in) :: lon1,lat1,lon2,lat2 ! Pos. in degrees
  220. real, parameter :: invradian=0.017453292
  221. real, parameter :: rearth=6371001.0 ! Radius of earth
  222. real rlon1,rlat1,rlon2,rlat2 ! Pos. in radians
  223. real x1,y1,z1,x2,y2,z2 ! Cartesian position
  224. real dx,dy,dz,dr ! Cartesian distances
  225. rlon1=lon1*invradian !lon1 in rad
  226. rlat1=(90.-lat1)*invradian !90-lat1 in rad
  227. rlon2=lon2*invradian !lon2 in rad
  228. rlat2=(90.-lat2)*invradian !90-lat2 in rad
  229. x1= SIN(rlat1)*COS(rlon1) !x,y,z of pos 1.
  230. y1= SIN(rlat1)*SIN(rlon1)
  231. z1= COS(rlat1)
  232. x2= SIN(rlat2)*COS(rlon2) !x,y,z of pos 2.
  233. y2= SIN(rlat2)*SIN(rlon2)
  234. z2= COS(rlat2)
  235. dx=x2-x1 !distances in x, y, z
  236. dy=y2-y1
  237. dz=z2-z1
  238. dr=SQRT(dx*dx+dy*dy+dz*dz) !distance pytagaros
  239. dr=acos(x1*x2+y1*y2+z1*z2) ! Acr length
  240. spherdist=dr*rearth
  241. end function spherdist
  242. end module mod_sphere_tools