m_read_jpl_hice.F90 3.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115
  1. module m_read_jpl_hice
  2. contains
  3. subroutine read_jpl_hice(fname, obstype, variance, nx, ny, data)
  4. use mod_measurement
  5. use m_oldtonew
  6. use m_confmap
  7. use m_bilincoeff
  8. use m_pivotp
  9. implicit none
  10. character(*), intent(in) :: fname
  11. character(*), intent(in) :: obstype
  12. real(8), intent(in) :: variance
  13. integer, intent(in) :: nx, ny
  14. type(measurement), allocatable, intent(out) :: data(:)
  15. type(measurement), allocatable :: tmpdata(:)
  16. real(8), dimension(nx, ny) :: modlat, modlon
  17. real(8), dimension(nx, ny) :: depths
  18. integer :: npoints
  19. integer :: i
  20. integer :: nobs
  21. real(8) :: lat, lon, tmp1, tmp2, h
  22. real(8) :: latnew, lonnew
  23. integer :: ipiv, jpiv
  24. open(101, file = trim(fname))
  25. read(101, *) npoints
  26. print *, ' ',trim(fname), ': ', npoints, ' data points'
  27. allocate(tmpdata(npoints))
  28. call confmap_init(nx, ny)
  29. call grid_readxyz(nx, ny, modlat, modlon, depths)
  30. nobs = 0
  31. do i = 1, npoints
  32. read(101, *) lat, lon, tmp1, tmp2, h
  33. if (h > 0.0d0 .and. h < 9990.0d0) then
  34. call oldtonew(lat, lon, latnew, lonnew)
  35. call pivotp(lonnew, latnew, ipiv, jpiv)
  36. if (ipiv < 1 .or. jpiv < 1 .or. ipiv > nx - 1 .or. jpiv > ny - 1) then
  37. cycle
  38. end if
  39. if (depths(ipiv, jpiv) < 10) then
  40. cycle
  41. end if
  42. nobs = nobs + 1
  43. tmpdata(nobs) % d = h / 100.0d0
  44. tmpdata(nobs) % id = obstype
  45. tmpdata(nobs) % var = variance
  46. tmpdata(nobs) % lon = lon
  47. tmpdata(nobs) % lat = lat
  48. tmpdata(nobs) % ipiv = ipiv
  49. tmpdata(nobs) % jpiv = jpiv
  50. tmpdata(nobs) % ns = 0 ! for a point (not gridded) measurement
  51. tmpdata(nobs) % date = 0 ! assimilate synchronously
  52. call bilincoeff(modlon, modlat, nx, ny, lon, lat, ipiv,&
  53. jpiv, tmpdata(nobs) % a1, tmpdata(nobs) % a2, tmpdata(nobs) % a3,&
  54. tmpdata(nobs) % a4)
  55. tmpdata(nobs) % status = .true. ! (active)
  56. tmpdata(nobs) % i_orig_grid = -1 ! not used
  57. tmpdata(nobs) % j_orig_grid = -1 ! not used
  58. end if
  59. end do
  60. close(101)
  61. print *, ' ', nobs, ' valid observations'
  62. allocate(data(nobs))
  63. do i = 1, nobs
  64. data(i) = tmpdata(i)
  65. end do
  66. deallocate(tmpdata)
  67. end subroutine read_jpl_hice
  68. subroutine grid_readxyz(nx, ny, lat, lon, depth)
  69. integer, intent(in) :: nx, ny
  70. real(8), dimension(nx, ny), intent(inout) :: lat, lon, depth
  71. logical :: exists
  72. character(len = 128) :: fname
  73. fname = 'newpos.uf'
  74. inquire(file = fname, exist = exists)
  75. if (.not. exists) then
  76. print *, 'grid_readxyz(): ERROR: "', trim(fname), '" does not exist'
  77. stop
  78. end if
  79. open(10, file = fname, form = 'unformatted', status = 'old')
  80. print *, ' grid_readxyz(): reading "', trim(fname), '"...'
  81. read(10) lat, lon
  82. close(10)
  83. write(fname, '(a, i3.3, a, i3.3, a)') 'depths', nx, 'x', ny, '.uf'
  84. inquire(file = fname, exist = exists)
  85. if (.not. exists) then
  86. print*, 'grid_readxyz(): ERROR: "', trim(fname), '" does not exist'
  87. stop
  88. end if
  89. open (unit = 10, file = fname, status = 'old', form = 'unformatted')
  90. print *, ' grid_readxyz(): reading "', trim(fname), '"...'
  91. read(10) depth
  92. close(10)
  93. end subroutine grid_readxyz
  94. end module m_read_jpl_hice