m_read_MET_SST.F90 4.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100
  1. module m_read_MET_SST
  2. ! Reads CLS SLA data after having read the grid in read_CLS_SST_grid
  3. integer, parameter, private :: STRLEN = 512
  4. contains
  5. subroutine read_MET_SST(filename,gr,data)
  6. use mod_measurement
  7. use mod_grid
  8. use m_spherdist
  9. use netcdf
  10. use nfw_mod
  11. implicit none
  12. type (measurement), intent(inout) :: data(:)
  13. type (grid), intent(inout) :: gr ! CLS measurement grid
  14. character(len=80), intent(in) :: filename
  15. ! Variable ids
  16. integer :: lon_ID, lat_ID,vsst_ID, vstd_ID, vmask_ID
  17. ! Data arrays
  18. real, allocatable :: sst(:,:), lon(:), lat(:), std(:,:)
  19. integer, allocatable :: mask(:,:)
  20. integer :: ncid ! observations
  21. real, dimension(1) :: undef_sst
  22. integer :: i, j, count1
  23. real, parameter :: eps = 0.01 ! test for undefined values
  24. ! filen name
  25. logical :: ex
  26. print *, 'read_MET_SST:'
  27. inquire(file=trim(filename),exist=ex)
  28. if (ex) then
  29. ! Reading the observation file
  30. call nfw_open(filename, nf_nowrite, ncid)
  31. ! Get dimension id in netcdf file ...
  32. !nb total of data
  33. allocate(lon(gr%nx), lat(gr%ny), sst(gr%nx,gr%ny), std(gr%nx, gr%ny), mask(gr%nx, gr%ny))
  34. ! Variable ids in netcdf file
  35. call nfw_inq_varid(filename, ncid, 'lat', lat_ID)
  36. call nfw_inq_varid(filename, ncid,'lon', lon_ID)
  37. call nfw_inq_varid(filename, ncid,'analysed_sst' ,vsst_ID)
  38. call nfw_inq_varid(filename, ncid,'analysis_error' ,vstd_ID)
  39. call nfw_inq_varid(filename, ncid,'mask' ,vmask_ID)
  40. ! Variable _FillValue attributes
  41. call nfw_get_att_double(filename, ncid, vsst_ID, '_FillValue', undef_sst(1))
  42. gr % undef = undef_sst(1)
  43. ! actual variable values (for dimensions of var -- see ncdump, or improve this program)
  44. ! NB: note that index dimensions are different between fortran and C internals.
  45. ! "ncdump" gives C internals.
  46. call nfw_get_var_double(filename, ncid, lon_ID, lon)
  47. call nfw_get_var_double(filename, ncid, lat_ID, lat)
  48. call nfw_get_var_double(filename, ncid, vsst_ID, sst)
  49. call nfw_get_var_double(filename, ncid, vstd_ID, std)
  50. call nfw_get_var_int(filename, ncid, vmask_ID, mask)
  51. print '(1x, a, 2f10.2)', ' range Lon = ', minval(lon), maxval(lon)
  52. print '(1x, a, 2f10.2)', ' range Lat = ', minval(lat), maxval(lat)
  53. print '(1x, a, 2f10.2)', ' range sst (K) = ', minval(sst), maxval(sst)
  54. print '(1x, a, 2i10)', ' range mask = ', minval(mask), maxval(mask)
  55. call nfw_close(filename, ncid)
  56. count1=1
  57. do i=1,gr%nx
  58. do j=1,gr%ny
  59. !here we only consider:
  60. !data above -30 of lat; valid, within reasonable range,
  61. ! and only open ocean (mask == 1)
  62. !
  63. if (lat(j) > -30 .and.&
  64. abs(sst(i,j)-undef_sst(1)) > eps .and.&
  65. mask(i,j) == 1 .and. &
  66. sst(i,j) > -190 .and. &
  67. sst(i,j) < 4500 .and. &
  68. std(i,j) > 0.0) then
  69. data(count1)%id = 'SST'
  70. data(count1)%d = sst(i,j)*0.01
  71. data(count1)%ipiv = count1 !whatever it is filled afterwards
  72. data(count1)%jpiv = 1 !whaterver it is filled afterwards
  73. data(count1)%lat=lat(j)
  74. data(count1)%lon=lon(i)
  75. data(count1)%a1 = spherdist(lon(i)-.5*gr%dx,lat(j)-.5*gr%dy,lon(i)+.5*gr%dx,lat(j)+.5*gr%dy)
  76. data(count1)%ns = 1 ! 1 for data with a spatial extent
  77. data(count1)%var = (std(i,j) * 0.01 * 8) ** 2 ! Exaggerate, factor 8
  78. data(count1)%date = 0
  79. data(count1)%depth = 0.0
  80. data(count1)%status = .true.
  81. count1=count1+1
  82. endif
  83. enddo !i
  84. enddo !j
  85. print*, ' # of obs read = ', count1
  86. deallocate(lat, lon, sst, mask)
  87. end if ! ex
  88. print *, 'MAX var(SST) = ', maxval(data(1 : count1) % var)
  89. end subroutine read_MET_SST
  90. end module m_read_MET_SST