m_read_CLS_SSH.F90 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143
  1. module m_read_CLS_SSH
  2. ! Reads CLS SSH data after having read the grid in read_CLS_SST_grid
  3. contains
  4. subroutine read_CLS_SSH(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
  19. ! Array dimensions
  20. integer :: LatLon, NbLatitudes, NbLongitudes
  21. ! Data arrays
  22. real, allocatable :: ssh(:,:), lon(:),lat(:)
  23. ! utilitary
  24. integer ncid, ijmax(2)
  25. real undef,undef_lat, undef_lon
  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. print*,'How far do you go'
  37. ! Get dimension length from id
  38. call nf90_err(nf90_Inquire_Dimension(ncid,LatLon_ID,len=LatLon))
  39. call nf90_err(nf90_Inquire_Dimension(ncid,NbLatitudes_ID,len=NbLatitudes))
  40. call nf90_err(nf90_Inquire_Dimension(ncid,NbLongitudes_ID,len=NbLongitudes))
  41. print*, 'Dimensions:', NbLatitudes, NbLongitudes, LatLon
  42. ! State which variable you want here.. Available vars are shown when you do
  43. ! "ncdump -h " on the netcdf file. This is for SSH
  44. allocate(lon(NbLongitudes))
  45. allocate(lat(NbLatitudes))
  46. ! allocate(latlon0(LatLon))
  47. ! allocate(dlatlon(LatLon))
  48. allocate(ssh(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. ! Variable _FillValue attributes
  54. call nf90_err(nf90_get_att(ncid,vNbLatitudes_ID , '_FillValue',undef_lat))
  55. call nf90_err(nf90_get_att(ncid,vNbLongitudes_ID ,'_FillValue',undef_lon))
  56. call nf90_err(nf90_get_att(ncid,vGrid0001_ID , '_FillValue',undef))
  57. print*, 'Undefined values are ', undef_lat, undef_lon, undef
  58. gr%undef = undef
  59. ! actual variable values (for dimensions of var -- see ncdump, or improve this program)
  60. ! NB: note that index dimensions are different between fortran and C internals.
  61. ! "ncdump" gives C internals.
  62. print *,'test'
  63. call nf90_err(nf90_get_var(ncid,vNbLongitudes_ID ,lon))
  64. !lon = ang180(lon)
  65. print *,'Range Lon', minval(lon), maxval(lon)
  66. call nf90_err(nf90_get_var(ncid,vNbLatitudes_ID ,lat))
  67. print *,'Range Lat', minval(lat), maxval(lat)
  68. call nf90_err(nf90_get_var(ncid,vGrid0001_ID ,ssh))
  69. ssh = ssh - 130. ! LB Correction of offset in MDT
  70. print *,'Range SSH in cm ', minval(ssh), maxval(ssh)
  71. ! print*, 'Latitudes'
  72. ! print '(12f8.2)',lat
  73. ! print*, 'Longitudes'
  74. ! print '(12f8.2)',lon
  75. print '(4a10)','Lat','Lon','SSH[cm]'
  76. ! print*,lat,lon,temp(i),err_temp(i),saln(i),err_saln(i)
  77. ! write(13,*)lat,lon,saln(i),err_saln(i),depth(i)
  78. ! write(14,*)lat,lon,temp(i),err_temp(i),depth(i)
  79. ijmax = minloc(ssh)
  80. do i=ijmax(1)-5, ijmax(1)+5
  81. j = ijmax(2)
  82. print '(4f10.3)', lat(i), lon(j), ssh(i,j)
  83. enddo
  84. call nf90_err (nf90_close(ncid))
  85. !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
  86. ! Fill the data(:) vector
  87. do j=1,NbLongitudes ! gr%ny
  88. do i=1,NbLatitudes ! gr%nx
  89. k=(j-1)*gr%nx+i
  90. data(k)%id = 'SSH'
  91. data(k)%d = ssh(i,j) * 0.01 ! Conversion to meters
  92. data(k)%ipiv = i
  93. data(k)%jpiv = j
  94. data(k)%lat=lat(i)
  95. data(k)%lon=ang180(lon(j))
  96. !LB: Data support is assumed = a square grid cell
  97. !support diameter in meters stored in %a1 (tricky, isn't it ?)
  98. data(k)%a1 = spherdist(lon(j)-.5*gr%dx,lat(i)-.5*gr%dy, &
  99. lon(j)+.5*gr%dx,lat(i)+.5*gr%dy)
  100. data(k)%ns = 1
  101. data(k)%var = 0.05 ! 20cm temporarily, 10 cm by default
  102. print*, 'SSH variance augmented: 22cm'
  103. data(k) % depth = 0.0
  104. valid = (abs( (lon(j)-undef_lon) / undef_lon ) > eps &
  105. .and. abs( (lat(i)-undef_lat) / undef_lat ) > eps &
  106. .and. abs( (ssh(i,j)-undef) / undef ) > eps )
  107. data(k)%status = valid
  108. enddo
  109. enddo
  110. print*, 'Number of data read:', k, gridpoints(gr)
  111. deallocate(lat,lon,ssh)
  112. end subroutine read_CLS_SSH
  113. end module m_read_CLS_SSH