m_get_def_wet_point.F90 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162
  1. ! Bilinear coeffisients are calculated witin this program for all the
  2. ! observation points.
  3. ! Only wet points are stored in the observation.uf file
  4. module m_get_def_wet_point
  5. implicit none
  6. integer, parameter, private :: STRLEN = 512
  7. character(STRLEN), parameter, private :: MEANSSHFNAME = "meanssh.uf"
  8. private read_mean_ssh
  9. private land_nearby
  10. contains
  11. subroutine get_def_wet_point(obs, data, gr, depths, modlat, modlon, nrobs, nx, ny)
  12. ! Program converts to a general format readable by ENKF (observations.uf)
  13. use mod_measurement
  14. use mod_grid
  15. ! Functions to be used
  16. use m_oldtonew
  17. use m_bilincoeff
  18. use m_pivotp
  19. use m_confmap
  20. use m_spherdist
  21. integer, intent(in) :: nx, ny
  22. type (measurement), intent(in) :: data(:)
  23. type (measurement), intent(inout) :: obs(:)
  24. type (grid), intent(in) :: gr ! observations grid
  25. real, dimension(nx, ny), intent(in) :: depths, modlat, modlon
  26. integer, intent(out) :: nrobs
  27. integer, parameter :: maxobs = 1441 * 722 !2*400*600 ! maximum number of observations
  28. real, dimension(nx, ny) :: mean_ssh
  29. integer k, imin, imax, jmin, jmax
  30. integer ipiv, jpiv, nsupport, nsmin, nsmax
  31. real :: x0, y0
  32. real wetsill, griddiag, mingridsize, minobssize
  33. logical wet
  34. ! gr = default_grid
  35. nrobs = 0;
  36. nsmin = maxobs;
  37. nsmax = 0
  38. mingridsize = 1.E+10;
  39. minobssize = 1.E+10 ! in meters
  40. call confmap_init(nx,ny) ! Initialize conformal mapping before calling
  41. !Calculate pivot points
  42. !Find wet points (all neigbours in water)
  43. !Find the points with defined data value
  44. !Put the data into the obs data structture
  45. !Compute bilinear coefficients
  46. call read_mean_ssh(mean_ssh, nx, ny)
  47. do k = 1, gridpoints(gr)
  48. if (data(k) % id .eq. 'SLA' .or. data(k) % id .eq. 'sla' .or. &
  49. data(k) % id.eq. 'SSH' .or. data(k)%id .eq. 'ssh' .or.&
  50. data(k)%id.eq.'TSLA') then
  51. wetsill = 200. ! Discarding data in shallow waters
  52. else
  53. wetsill=10.
  54. endif
  55. call oldtonew(data(k) % lat, data(k) % lon, y0, x0)
  56. call pivotp(x0, y0, ipiv, jpiv)
  57. ! Discard obs on model boundaries (TODO: cyclic domains)
  58. ! Also valid if ns=0
  59. imin = ipiv - data(k) % ns
  60. imax = ipiv + data(k) % ns + 1
  61. jmin = jpiv - data(k) % ns
  62. jmax = jpiv + data(k) % ns + 1
  63. if ((imin .le. 0) .or. (jmin .le. 0) .or. (imax .ge. nx) .or. &
  64. (jmax .ge. ny)) cycle
  65. ! Is observation surrounded by wet grid points?
  66. if (any(depths(imin:imax, jmin:jmax) < wetsill .or. depths(imin:imax, jmin:jmax) == depths(imin:imax, jmin:jmax) + 1.0)) cycle
  67. wet = data(k) % status ! Discards inconsistent/Fill values
  68. if (data(k) % id .eq. 'SLA' .or. data(k) % id .eq. 'sla' .or.&
  69. data(k) % id .eq. 'TSLA') then
  70. wet = wet .and. (mean_ssh(ipiv, jpiv) < 990)
  71. wet = wet .and. .not. land_nearby(nx, ny, mean_ssh, modlon, modlat,&
  72. ipiv, jpiv, data(k) % lon, data(k) % lat)
  73. endif
  74. if(.not. undefined(data(k) % d, gr) .and. wet) then
  75. nrobs = nrobs + 1
  76. obs(nrobs) = data(k)
  77. obs(nrobs) % ipiv = ipiv
  78. obs(nrobs) % jpiv= jpiv
  79. obs(nrobs) % status = .true. ! Wet obs
  80. if (data(k) % ns > 0) then ! large support data: a1 is the obs support(m)
  81. griddiag = spherdist(modlon(ipiv, jpiv), modlat(ipiv, jpiv), &
  82. modlon(ipiv + 1, jpiv + 1), modlat(ipiv + 1, jpiv + 1))
  83. !FC: 0.5 because m_Generate_element_Si runs from -ns to +ns
  84. nsupport = floor(0.5 * data(k) % a1 / griddiag) ! both in meters
  85. obs(nrobs)%ns = nsupport ! number of grid cells in the diagonal
  86. nsmin = min(nsmin, nsupport)
  87. nsmax = max(nsmax, nsupport)
  88. mingridsize = min(mingridsize, griddiag)
  89. minobssize = min(minobssize, data(k) % a1)
  90. else
  91. obs(nrobs) % ns = 0 ! point measurements have zero support
  92. endif
  93. call bilincoeff(modlon, modlat, nx, ny, obs(nrobs)%lon, &
  94. obs(nrobs) % lat, obs(nrobs) % ipiv, obs(nrobs) % jpiv, &
  95. obs(nrobs) % a1, obs(nrobs) % a2, obs(nrobs) % a3, &
  96. obs(nrobs) % a4)
  97. endif
  98. end do
  99. print*, 'Number of defined and wet observations: nrobs ', nrobs
  100. print*, 'Support (in nb of cells) between: ', nsmin, ' and ', nsmax
  101. print '(2(a,f8.3),a)', ' Minimum obs support: ', 0.001 * minobssize, &
  102. 'km, min grid diagonal: ', 0.001 * mingridsize, ' km'
  103. end subroutine get_def_wet_point
  104. subroutine read_mean_ssh(mean_ssh, nx, ny)
  105. integer, intent(in) :: nx, ny
  106. real, intent(out):: mean_ssh(nx, ny)
  107. logical :: exists
  108. inquire(file = trim(MEANSSHFNAME), exist = exists)
  109. if (.not. exists) then
  110. print *,'ERROR: read_mean_ssh(): file "', trim(MEANSSHFNAME), '" not found'
  111. stop
  112. end if
  113. open (10, file = trim(MEANSSHFNAME), status = 'unknown',form = 'unformatted')
  114. read (10) mean_ssh
  115. close (10)
  116. end subroutine read_mean_ssh
  117. logical function land_nearby(nx, ny, mean_ssh, modlon, modlat, ipiv, jpiv, obslon, obslat)
  118. use m_spherdist
  119. implicit none
  120. real, parameter :: Dis0 = 50.0d0
  121. integer, intent (in) :: nx, ny, ipiv, jpiv
  122. real, dimension(nx,ny), intent(in) :: mean_ssh, modlon, modlat
  123. real, intent (in) :: obslon,obslat
  124. integer :: ii, jj, ncells
  125. real :: griddist
  126. land_nearby = .false.
  127. ncells = ceiling(Dis0 / spherdist(modlon(ipiv, jpiv), modlat(ipiv, jpiv),&
  128. modlon(ipiv, jpiv + 1), modlat(ipiv, jpiv + 1)))
  129. do jj = max(jpiv - ncells, 1), min(jpiv + ncells, ny)
  130. do ii = max(ipiv - ncells, 1), min(ipiv + ncells, nx)
  131. griddist = spherdist(modlon(ii, jj), modlat(ii, jj), obslon, obslat)
  132. if (mean_ssh(ipiv,jpiv) < 990 .and. griddist < Dis0) then
  133. land_nearby = .true.
  134. return
  135. end if
  136. enddo
  137. enddo
  138. end function land_nearby
  139. end module m_get_def_wet_point