m_bilincoeff.F90 3.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107
  1. module m_bilincoeff
  2. use m_oldtonew
  3. implicit none
  4. contains
  5. ! This subroutine uses bilinear interpolation to interpolate the field
  6. ! computed by the model (MICOM) to the position defined by lon, lat
  7. ! The output is the interpolation coeffisients a[1-4]
  8. ! NB NO locations on land.
  9. !
  10. subroutine bilincoeff(glon, glat, nx, ny, lon, lat, ipiv, jpiv, a1, a2, a3, a4)
  11. real, intent(in) :: glon(nx, ny), glat(nx, ny)
  12. integer, intent(in) :: nx ,ny
  13. real, intent(in) :: lon, lat
  14. integer, intent(in) :: ipiv, jpiv
  15. real, intent(out) :: a1, a2, a3, a4
  16. real :: t, u
  17. real :: lat1, lon1, lat2, lon2, latn, lonn
  18. call oldtonew(glat(ipiv, jpiv), glon(ipiv, jpiv), lat1, lon1)
  19. call oldtonew(glat(ipiv + 1, jpiv + 1), glon(ipiv + 1, jpiv + 1), lat2, lon2)
  20. call oldtonew(lat, lon, latn, lonn)
  21. t = (lonn - lon1) / (lon2 - lon1)
  22. u = (latn - lat1) / (lat2 - lat1)
  23. if (t < -0.1 .or. t > 1.1 .or. u < -0.1 .or. u > 1.1) then
  24. print *, 'ERROR: bilincoeff(): t, u = ', t, u, 'for lon, lat =', lon, lat
  25. stop
  26. end if
  27. a1 = (1.0 - t) * (1.0 - u)
  28. a2 = t * (1.0 - u)
  29. a3 = t * u
  30. a4 = (1.0 - t) * u
  31. end subroutine bilincoeff
  32. subroutine bilincoeff1(glon, glat, nx, ny, lon, lat, ipiv, jpiv, a1, a2, a3, a4)
  33. real, intent(in) :: glon(nx, ny), glat(nx, ny)
  34. integer, intent(in) :: nx ,ny
  35. real, intent(in) :: lon, lat
  36. integer, intent(in) :: ipiv, jpiv
  37. real, intent(out) :: a1, a2, a3, a4
  38. real :: xx(4), yy(4)
  39. real :: t, u
  40. xx(1) = glon(ipiv, jpiv)
  41. xx(2) = glon(ipiv + 1, jpiv)
  42. xx(3) = glon(ipiv + 1, jpiv + 1)
  43. xx(4) = glon(ipiv, jpiv + 1)
  44. yy(1) = glat(ipiv, jpiv)
  45. yy(2) = glat(ipiv + 1, jpiv)
  46. yy(3) = glat(ipiv + 1, jpiv + 1)
  47. yy(4) = glat(ipiv, jpiv + 1)
  48. call xy2fij(lon, lat, xx, yy, t, u)
  49. if (t < 0 .or. t > 1 .or. u < 0 .or. u > 1) then
  50. print *, 'ERROR: bilincoeff(): t, u = ', t, u, 'for lon, lat =', lon, lat
  51. ! stop
  52. end if
  53. a1 = (1.0 - t) * (1.0 - u)
  54. a2 = t * (1.0 - u)
  55. a3 = t * u
  56. a4 = (1.0 - t) * u
  57. end subroutine bilincoeff1
  58. subroutine xy2fij(x, y, xx, yy, fi, fj)
  59. real, intent(in) :: x, y
  60. real, intent(in) :: xx(4), yy(4)
  61. real, intent(out) :: fi, fj
  62. real :: a, b, c, d, e, f, g, h
  63. real :: aa, bb, cc
  64. real :: d1, d2
  65. a = xx(1) - xx(2) - xx(4) + xx(3)
  66. b = xx(2) - xx(1)
  67. c = xx(4) - xx(1)
  68. d = xx(1)
  69. e = yy(1) - yy(2) - yy(4) + yy(3)
  70. f = yy(2) - yy(1)
  71. g = yy(4) - yy(1)
  72. h = yy(1)
  73. aa = a * f - b * e;
  74. bb = e * x - a * y + a * h - d * e + c * f - b * g;
  75. cc = g * x - c * y + c * h - d * g;
  76. if (abs(aa) < 1d-5) then
  77. fi = -cc / bb * (1.0d0 + aa * cc / bb / bb);
  78. else
  79. fi = (-bb - sqrt(bb * bb - 4.0d0 * aa * cc)) / (2.0d0 * aa);
  80. end if
  81. d1 = a * fi + c
  82. d2 = e * fi + g
  83. if (abs(d2) > abs(d1)) then
  84. fj = (y - f * fi - h) / d2
  85. else
  86. fj = (x - b * fi - d) / d1
  87. end if
  88. end subroutine xy2fij
  89. end module m_bilincoeff