m_read_OSISAF_data.F90 4.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143
  1. module m_read_OSISAF_data
  2. contains
  3. ! Reads OSISAF ice drift data
  4. ! 2012-11-13: Geir Arne Waagbø (met.no)
  5. subroutine read_OSISAF_data(driftfile, gr, data, numdata, var, offset)
  6. use nfw_mod
  7. use mod_measurement
  8. use mod_grid
  9. implicit none
  10. character(*), intent(in) :: driftfile
  11. integer, intent(in) :: numdata
  12. type(measurement), dimension(numdata),intent(inout) :: data
  13. type(grid), intent(in) :: gr
  14. real, intent(in) :: var
  15. character(len=1), intent(in) :: offset
  16. integer :: dimids(2)
  17. integer , dimension(2) :: dimsizes
  18. integer :: lon_id, lat_id, dx_id, dy_id, qflag_id
  19. real, dimension(:,:), allocatable :: drlon, drlat, drdX, drdY
  20. integer, dimension(:,:), allocatable :: qflag
  21. integer :: ncid
  22. real, dimension(1) :: fillval
  23. integer :: i,j,k,icomp
  24. integer :: drnx, drny
  25. logical :: valid
  26. ! Get dimensions of drift file
  27. call nfw_open(driftfile, nf_nowrite, ncid)
  28. call nfw_inq_varid(driftfile, ncid, 'dX', dx_id)
  29. call nfw_inq_vardimid(driftfile, ncid, dx_id, dimids)
  30. do i = 1, 2
  31. call nfw_inq_dimlen(driftfile, ncid, dimids(i), dimsizes(i))
  32. end do
  33. drnx=dimsizes(1)
  34. drny=dimsizes(2)
  35. if (gr % reg) then
  36. print *,'NB: OSISAF data should be specified as irregular !'
  37. print *,' Currently it is set as regular..'
  38. print *,'(read_OSISAF_data)'
  39. call exit(1)
  40. end if
  41. ! Which should match numdata dimension
  42. ! NB !!! Mult by 2 for two vector components
  43. if (2*drnx*drny /= numdata .or. 2*drnx*drny /= gr%nx) then
  44. print *,'Different dimensions - data file and specified'
  45. print *,'dimsizes(1)=',drnx
  46. print *,'dimsizes(2)=',drny
  47. print *,'nx =',gr%nx
  48. print *,'(read_OSISAF_data)'
  49. call exit(1)
  50. end if
  51. ! Read data from drift file
  52. allocate(drlon(drnx,drny))
  53. allocate(drlat(drnx,drny))
  54. allocate(drdX(drnx,drny))
  55. allocate(drdY(drnx,drny))
  56. allocate(qflag(drnx,drny))
  57. call nfw_inq_varid(driftfile, ncid, 'lon', lon_id)
  58. call nfw_get_var_double(driftfile, ncid, lon_id, drlon)
  59. call nfw_inq_varid(driftfile, ncid, 'lat', lat_id)
  60. call nfw_get_var_double(driftfile, ncid, lat_id, drlat)
  61. call nfw_inq_varid(driftfile, ncid, 'dX', dx_id)
  62. call nfw_get_var_double(driftfile, ncid, dx_id, drdX)
  63. ! Change dY_v1p4 to dY from version 1.4 of OSISAF product file
  64. call nfw_inq_varid(driftfile, ncid, 'dY_v1p4', dy_id)
  65. call nfw_get_var_double(driftfile, ncid, dy_id, drdY)
  66. call nfw_get_att_double(driftfile, ncid, dx_id, '_FillValue', fillval)
  67. where (abs(drdX - fillval(1)) < 1e-4 * fillval(1))
  68. drdX = gr % undef
  69. end where
  70. call nfw_get_att_double(driftfile, ncid, dy_id, '_FillValue', fillval)
  71. where (abs(drdY - fillval(1)) < 1e-4 * fillval(1))
  72. drdY = gr % undef
  73. end where
  74. ! Change status_flag_v1p4 to status_flag from version 1.4 of OSISAF product file
  75. call nfw_inq_varid(driftfile, ncid, 'status_flag_v1p4', qflag_id)
  76. call nfw_get_var_int(driftfile, ncid, qflag_id, qflag)
  77. call nfw_close(driftfile, ncid)
  78. k = 0
  79. do icomp = 1, 2
  80. do j = 1, drny
  81. do i = 1, drnx
  82. k = k + 1
  83. valid = qflag(i,j) == 30 ! Only use observations with quality_flag==30
  84. if (icomp==1) then
  85. data(k)%id = 'DX'//offset
  86. data(k)%d = drdX(i,j)
  87. valid = valid .and. abs((drdX(i,j)-gr%undef) / gr%undef) > 1e-4
  88. else
  89. data(k)%id = 'DY'//offset
  90. data(k)%d = drdY(i,j)
  91. valid = valid .and. abs((drdY(i,j)-gr%undef) / gr%undef) > 1e-4
  92. end if
  93. if (.not. valid .or. mod(i,2)==0 .or. mod(j,2)==0) then
  94. ! Skip invalid observations, or observations on grid points with
  95. ! even i- or j-indices (to avoid over assimilation)
  96. data(k)%d = gr%undef
  97. end if
  98. data(k)%ipiv = i ! Not really used for ice drift
  99. data(k)%jpiv = j ! Not really used for ice drift
  100. data(k)%i_orig_grid = i ! Used for ice drift
  101. data(k)%j_orig_grid = j ! Used for ice drift
  102. data(k)%lat=drlat(i,j)
  103. data(k)%lon=ang180(drlon(i,j))
  104. ! Each vector represents the average drift of a 120kmx120km area of sea ice
  105. ! The a1 value should be in meters, although other values are in km
  106. data(k)%a1 = 1.4*60000 ! 1.4 represents square root of 2
  107. data(k)%ns = 1
  108. data(k)%var = var ! fom idrft.hdr specification
  109. data(k)%depth = 0.0
  110. data(k)%status = valid
  111. enddo
  112. enddo
  113. enddo
  114. print *, 'Number of data read:', k, gridpoints(gr)
  115. end subroutine
  116. end module