m_read_CLS_SST.F90 5.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146
  1. module m_read_CLS_SST
  2. ! Reads CLS SST data after having read the grid in read_CLS_SST_grid
  3. contains
  4. subroutine read_CLS_SST(fname,gr,data)
  5. use mod_measurement
  6. ! use mod_dimensions
  7. use mod_grid
  8. use m_spherdist
  9. use netcdf
  10. use m_nf90_err
  11. implicit none
  12. type (measurement), intent(inout) :: data(:)
  13. type (grid), intent(inout) :: gr ! CLS measurement grid
  14. character(len=80), intent(in) :: fname
  15. !dimension ids
  16. integer :: NbLatitudes_ID, NbLongitudes_ID, LatLon_ID
  17. ! Variable ids
  18. integer :: vNbLatitudes_ID, vNbLongitudes_ID, vGrid0001_ID, vGrid0004_ID
  19. ! Array dimensions
  20. integer :: LatLon, NbLatitudes, NbLongitudes
  21. ! Data arrays
  22. real,allocatable :: sst(:,:), sst_var(:,:), lon(:),lat(:)!, latlon0(:), dlatlon(:)
  23. ! utilitary
  24. integer ncid, ijmax(2)
  25. real undef,undef_lat, undef_lon, undef_var
  26. integer i, j,k
  27. logical valid
  28. real, parameter :: eps = 0.01 ! test for undefined values
  29. ! Open file
  30. ! filename='sst_topaz_19510.nc'
  31. call nf90_err(NF90_OPEN(trim(fname),NF90_NOCLOBBER,ncid))
  32. ! Get dimension id in netcdf file ...
  33. call nf90_err(nf90_Inq_Dimid(ncid,'LatLon',LatLon_ID))
  34. call nf90_err(nf90_Inq_Dimid(ncid,'NbLatitudes',NbLatitudes_ID))
  35. call nf90_err(nf90_Inq_Dimid(ncid,'NbLongitudes',NbLongitudes_ID))
  36. ! Get dimension length from id
  37. call nf90_err(nf90_Inquire_Dimension(ncid,LatLon_ID,len=LatLon))
  38. call nf90_err(nf90_Inquire_Dimension(ncid,NbLatitudes_ID,len=NbLatitudes))
  39. call nf90_err(nf90_Inquire_Dimension(ncid,NbLongitudes_ID,len=NbLongitudes))
  40. print*, 'Dimensions:', NbLatitudes, NbLongitudes, LatLon
  41. ! State which variable you want here.. Available vars are shown when you do
  42. ! "ncdump -h " on the netcdf file. This is for SST and SDEV of SST
  43. allocate(lon(NbLongitudes))
  44. allocate(lat(NbLatitudes))
  45. ! allocate(latlon0(LatLon))
  46. ! allocate(dlatlon(LatLon))
  47. allocate(sst_var(NbLatitudes,NbLongitudes))
  48. allocate(sst(NbLatitudes,NbLongitudes))
  49. ! Variable ids in netcdf file
  50. call nf90_err(nf90_inq_varid(ncid,'NbLatitudes' ,vNbLatitudes_ID),'NbLatitudes')
  51. call nf90_err(nf90_inq_varid(ncid,'NbLongitudes' ,vNbLongitudes_ID),'NbLongitudes')
  52. call nf90_err(nf90_inq_varid(ncid,'Grid_0001' ,vGrid0001_ID),'Grid_0001')
  53. call nf90_err(nf90_inq_varid(ncid,'Grid_0004' ,vGrid0004_ID),'Grid_0004')
  54. ! Variable _FillValue attributes
  55. call nf90_err(nf90_get_att(ncid,vNbLatitudes_ID , '_FillValue',undef_lat))
  56. call nf90_err(nf90_get_att(ncid,vNbLongitudes_ID ,'_FillValue',undef_lon))
  57. call nf90_err(nf90_get_att(ncid,vGrid0001_ID , '_FillValue',undef))
  58. call nf90_err(nf90_get_att(ncid,vGrid0004_ID , '_FillValue',undef_var))
  59. print*, 'Undefined values are ', undef_lat, undef_lon, undef, undef_var
  60. gr%undef = undef
  61. ! actual variable values (for dimensions of var -- see ncdump, or improve this program)
  62. ! NB: note that index dimensions are different between fortran and C internals.
  63. ! "ncdump" gives C internals.
  64. print *,'test'
  65. call nf90_err(nf90_get_var(ncid,vNbLongitudes_ID ,lon))
  66. print *,'Range Lon', minval(lon), maxval(lon)
  67. call nf90_err(nf90_get_var(ncid,vNbLatitudes_ID ,lat))
  68. print *,'Range Lat', minval(lat), maxval(lat)
  69. call nf90_err(nf90_get_var(ncid,vGrid0001_ID ,sst))
  70. print *,'Range SST', minval(sst), maxval(sst)
  71. call nf90_err(nf90_get_var(ncid,vGrid0004_ID ,sst_var))
  72. print *,'Range Std. Dev.', minval(sst_var), maxval(sst_var)
  73. ! print*, 'Latitudes'
  74. ! print '(12f8.2)',lat
  75. ! print*, 'Longitudes'
  76. ! print '(12f8.2)',lon
  77. print '(4a10)','Lat','Lon','Temp[C]','Error[C]'
  78. ! print*,lat,lon,temp(i),err_temp(i),saln(i),err_saln(i)
  79. ! write(13,*)lat,lon,saln(i),err_saln(i),depth(i)
  80. ! write(14,*)lat,lon,temp(i),err_temp(i),depth(i)
  81. ijmax = minloc(sst)
  82. do i=ijmax(1)-5, ijmax(1)+5
  83. j = ijmax(2)
  84. print '(4f10.3)', lat(i), lon(j), sst(i,j), sst_var(i,j)
  85. enddo
  86. call nf90_err (nf90_close(ncid))
  87. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  88. ! Fill the data(:) vector
  89. do j=1,NbLongitudes ! gr%ny
  90. do i=1,NbLatitudes ! gr%nx
  91. k=(j-1)*gr%nx+i
  92. data(k)%id = 'SST'
  93. data(k)%d = sst(i,j)
  94. data(k)%ipiv = i
  95. data(k)%jpiv = j
  96. data(k)%lat=lat(i)
  97. data(k)%lon=lon(j)
  98. !LB: Data support is assumed = a square grid cell
  99. !support diameter in meters stored in %a1 (tricky, isn't it ?)
  100. data(k)%a1 = spherdist(lon(j)-.5*gr%dx,lat(i)-.5*gr%dy, &
  101. lon(j)+.5*gr%dx,lat(i)+.5*gr%dy)
  102. data(k)%ns = 1
  103. data(k)%var = sst_var(i,j) ! corrected: variance is provided
  104. data(k)%depth = 0.0
  105. valid = (abs( (lon(j)-undef_lon) / undef_lon ) > eps &
  106. .and. abs( (lat(i)-undef_lat) / undef_lat ) > eps &
  107. .and. abs( (sst(i,j)-undef) / undef ) > eps &
  108. .and. abs( (sst_var(i,j)-undef_var)/undef_var) > eps &
  109. .and. sst_var(i,j)> 0 &
  110. .and. sst_var(i,j)< 16 ) ! Sdev too high => perturbations too large
  111. data(k)%status = valid
  112. enddo
  113. enddo
  114. print*, 'Number of data read:', k, gridpoints(gr)
  115. deallocate(lat,lon,sst,sst_var)
  116. end subroutine read_CLS_SST
  117. end module m_read_CLS_SST