p_prep_obs_ORCA25.F90 8.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240
  1. program prep_obs
  2. ! This takes a model restart file, extracts the desired variable
  3. ! and brings into a format that the EnKF V2 can read & treat ('observations.uf').
  4. !
  5. ! !!! AGAIN: THIS USES EnKF VERSION 2 !!!
  6. !
  7. ! Two command line arguments are expected:
  8. ! 1. Path to the nc file of which the data is to be extracted.
  9. ! 2. Variable name that can be found in there.
  10. !
  11. ! Output is written to 'observations.uf'. If lots of files are to be treated,
  12. ! use a shell script to call me and rename output.
  13. !
  14. ! Warning: Currently only the data from the first time step are being read.
  15. ! (No problem with restart files, normally)
  16. !
  17. !
  18. ! TO DO: Add possibility of treating several obs types.
  19. !
  20. !
  21. ! (c) September 2008, Christof König Beatty, Christof.KonigBeatty@uclouvain.be
  22. ! (c) May 2009, generalized by same.
  23. ! (c) Jun 2011, simplified by F. Massonnet
  24. ! (c) April 2016, cleaned by F. Massonnet
  25. use mod_measurement
  26. use netcdf
  27. implicit none
  28. ! NetCDF vars
  29. character(len=99) :: filename, cvarerror, cvarerroru, cvarerrorv
  30. integer :: ncid, ierr, nvar, narg
  31. logical :: ex
  32. character(len=16) :: varname, varnameu, varnamev ! Type of measurement ('a_i_htcX', 'u_ice', 'v_ice', maybe 'v_i_htcX')
  33. character(len=80), parameter :: maskfile = 'mask-ORCA25.nc' !hc!
  34. ! Data vars
  35. integer, parameter :: nx=1442, ny=1050 ! Unfortunately, still hard coded.
  36. real, dimension(nx,ny) :: lons, lats, obsfld, errorfld, obsfldu, obsfldv, errorfldu, errorfldv
  37. integer, dimension(nx,ny) :: mask
  38. integer :: obscnt ! Counts nr of obs's.
  39. ! Other vars
  40. character(len=99) dummy ! To read cmd line args
  41. ! for loops (haha)
  42. integer :: i, j, varID, icomp
  43. ! Ice thickness category stuff
  44. integer, parameter :: nhtc=5 !hc! nr of ice thickn. cat's
  45. real(KIND=8) :: rdate
  46. ! Obs stuff
  47. type (measurement), allocatable :: obs(:)
  48. ! Need to fill:
  49. ! d (val), var (error var), id (iiceconc etc.), lon, lat, depths,
  50. ! ipic, jpic (i/j-pivot point in grid), ns (support, 0: point meas.),
  51. ! a1-4 (bilin. coeff), status (not used)
  52. narg= iargc()
  53. PRINT *, narg
  54. if (narg<=1) then
  55. ! Write info
  56. write(*,*)
  57. write(*,*) " (prep_obs) takes a real obs, extracts the desired variable and outputs"
  58. write(*,*) " it in a format that the EnKF can read & treat ('observations.uf')."
  59. write(*,*)
  60. write(*,*) " A file named mask.nc containing the variables tmaskutil, nav_lon and nav_lat"
  61. write(*,*) " is expected to be in the current directory (ORCA-file)"
  62. write(*,*)
  63. write(*,*) " Three command line arguments are expected:"
  64. write(*,*) " 1. Path to the nc file of which the data is to be extracted."
  65. write(*,*) " 2. Variable name that can be found in there, 'h_i_htc1' or"
  66. write(*,*) " 'sic'. or dxdy_ice"
  67. write(*,*) " 3. A tag with the date, e.g. 19790520"
  68. write(*,*)
  69. write(*,*) " Hope to see you again soon."
  70. write(*,*)
  71. stop "(prep_obs): Stopped."
  72. end if
  73. ! Command line arguments
  74. nvar=narg-1
  75. ! Get them
  76. call getarg(1, dummy)
  77. read(dummy,'(A)') filename
  78. ! Some parameter checking
  79. inquire(file=trim(filename), exist=ex)
  80. if (.not.ex) then
  81. print *, '(p_prep_obs): file does not exist: '// trim(filename)
  82. stop
  83. end if
  84. ! Get mask, lons & lats
  85. ! open the maskfile
  86. ierr = nf90_open(trim(maskfile),nf90_NoWrite,ncid); if (ierr.ne.nf90_noerr) call handle_err(ierr, "opening mask file")
  87. ! Find VarID of tmask & get values
  88. ierr = nf90_inq_varid(ncid, 'tmask', varID); if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID tmask")
  89. ierr = nf90_get_var(ncid, varID, mask) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable tmaks")
  90. ! Find VarID of longitude & get vals
  91. ierr = nf90_inq_varid(ncid, 'nav_lon', varID); if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID nav_lon")
  92. ierr = nf90_get_var(ncid, varID, lons) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable nav_lon")
  93. ! Find VarID of latitude & get vals
  94. ierr = nf90_inq_varid(ncid, 'nav_lat', varID); if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID nav_lat")
  95. ierr = nf90_get_var(ncid, varID, lats) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable nav_lat")
  96. ! Close maskfile
  97. ierr = nf90_close(ncid)
  98. if (ierr.ne.nf90_noerr) call handle_err(ierr, "closing")
  99. allocate( obs(nvar*sum(mask)), STAT=ierr )
  100. if (ierr.ne.0) call handle_err(ierr, "allocating obs") !no netcdf-error! ohwell.
  101. obscnt = 0 ! Keeps track of nr of added obs's.
  102. call getarg(2, dummy)
  103. read(dummy,'(A)') varname
  104. call getarg(3, dummy)
  105. read(dummy,*) rdate
  106. if ( trim(varname).ne.'sic' ) then
  107. write(*,*) ' (prep_obs) Currently only the total sea ice concentration'
  108. write(*,*) ' sic can be processed '
  109. stop " (prep_obs) Aborted"
  110. endif
  111. ! User info
  112. write(*,*) ' (prep_obs) Extracting '//trim(varname)//' from '//trim(filename)
  113. ! Some preparations
  114. obsfld=0.
  115. ! open the netCDF file
  116. ierr = nf90_open(trim(filename),nf90_NoWrite,ncid) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "opening")
  117. ! Read observation data
  118. ! Find VarID of varname
  119. ierr = nf90_inq_varid(ncid, trim(varname), varID) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID")
  120. ! get values
  121. ierr = nf90_get_var(ncid, varID, obsfld) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable")
  122. ! Find VarID of varname
  123. cvarerror='error_'//varname
  124. ierr = nf90_inq_varid(ncid, cvarerror, varID) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID")
  125. ! get values
  126. ierr = nf90_get_var(ncid, varID, errorfld) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable")
  127. ! Close file
  128. ierr = nf90_close(ncid) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "closing")
  129. ! User info - ADAPT ACCORDINGLY
  130. write (*,*) " (prep_obs) Using data >40N and <45S"
  131. ! Loop over space
  132. do i=1,size(mask,1)
  133. do j=1,size(mask,2)
  134. ! Pick out ocean points where data is existent
  135. if ( errorfld(i,j) .GE. 1.0 ) then
  136. ! Limit 'obs' to >40°N or <45°S
  137. if ( (lats(i,j).ge.40) .or. (lats(i,j).le.-45) ) then
  138. obscnt = obscnt + 1
  139. obs(obscnt)%d = obsfld(i,j)/100.0
  140. obs(obscnt)%lon = lons(i,j)
  141. obs(obscnt)%lat = lats(i,j)
  142. obs(obscnt)%ipiv = i ! the i-point of the grid of the model
  143. obs(obscnt)%jpiv = j ! the j-point of the grid of the model
  144. ! Put other data into obs type array
  145. obs(obscnt)%var = (errorfld(i,j)/100.0)**2 ! The variance
  146. obs(obscnt)%id = TRIM(varname)
  147. obs(obscnt)%depths = 0
  148. obs(obscnt)%ns = 0
  149. obs(obscnt)%a1 = 1
  150. obs(obscnt)%a2 = 0
  151. obs(obscnt)%a3 = 0
  152. obs(obscnt)%a4 = 0
  153. obs(obscnt)%status = .TRUE.
  154. obs(obscnt)%i_orig_grid = -1
  155. obs(obscnt)%j_orig_grid = -1
  156. obs(obscnt)%h = -1.0
  157. obs(obscnt)%date = rdate
  158. obs(obscnt)%orig_id = 0
  159. end if ! >40°N or <50°S
  160. end if ! ocean point
  161. end do ! j
  162. end do ! i
  163. !Write data:
  164. inquire(iolength=i)obs(1)
  165. open (11, file='observations.uf', status='replace',form='unformatted', access='direct', recl=i)
  166. do j=1,obscnt
  167. write(11,rec=j)obs(j)
  168. enddo
  169. close(11)
  170. PRINT *, rdate
  171. ! !XXX:
  172. ! ! Write data to textfile, for visual checking
  173. open(12, file='observations.txt')
  174. do j=1,obscnt
  175. write(12,FMT=53) obs(j)
  176. 53 FORMAT(f8.4,X,f8.4,X,a8,X,2(f10.5,X),f4.2,X,2(I3,X),I1,X,4(f5.2,X),L,X,2(I3,X),f5.2,X,I8,X,I1)
  177. enddo
  178. close(12)
  179. ! !:XXX
  180. ! Tell user how many obs there were
  181. write(*,*) " (prep_obs) Wrote out this many obs: ", obscnt
  182. write(*,*) " (prep_obs) Number of ocean points : ", sum(mask)
  183. ! Cleanup
  184. if (allocated(obs)) deallocate(obs,STAT=ierr)
  185. write(*,*) ' (prep_obs) End successfully reached'; write(*,*)
  186. contains
  187. subroutine handle_err(status, infomsg)
  188. integer, intent ( in) :: status
  189. character(len = *), intent ( in), optional :: infomsg
  190. if(status /= nf90_noerr) then
  191. if (present(infomsg)) then
  192. print *, '(prep_obs) Error while '//infomsg//' - '//trim(nf90_strerror(status))
  193. else
  194. print *, trim(nf90_strerror(status))
  195. endif ! opt arg
  196. print *,'(prep_obs)'
  197. stop
  198. end if ! check error status
  199. end subroutine handle_err
  200. end program prep_obs