p_prep_obs_ORCA1.F90 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285
  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', 'siconc' if model forced by itself)
  33. !!!!!TAG-DEV-AD : change mask-ORCA1.nc into mask-eORCA1.nc
  34. ! character(len=80), parameter :: maskfile = 'mask-ORCA1.nc' !hc!
  35. character(len=80), parameter :: maskfile = 'mask-eORCA1.nc' !hc!
  36. ! Data vars
  37. !!!!!TAG-DEV-AD : adapted coordinate for eORCA1
  38. integer, parameter :: nx=360, ny=331 ! Unfortunately, still hard coded.
  39. ! integer, parameter :: nx=362, ny=292 ! Unfortunately, still hard coded.
  40. real, dimension(nx,ny) :: lons, lats, obsfld, errorfld, obsfldu, obsfldv, errorfldu, errorfldv
  41. REAL :: obsmin, obsmax, errmin, errmax
  42. REAL :: latmin_NH = 40.0
  43. REAL :: latmax_NH = 90.0
  44. REAL :: latmin_SH = 40.0 ! Same for SH (sign will be added)
  45. REAL :: latmax_SH = 90.0
  46. integer, dimension(nx,ny) :: mask
  47. integer :: obscnt ! Counts nr of obs's.
  48. ! Other vars
  49. character(len=99) dummy ! To read cmd line args
  50. ! for loops (haha)
  51. integer :: i, j, varID, icomp
  52. ! Ice thickness category stuff
  53. integer, parameter :: nhtc=5 !hc! nr of ice thickn. cat's
  54. real(KIND=8) :: rdate
  55. ! Obs stuff
  56. type (measurement), allocatable :: obs(:)
  57. ! Need to fill:
  58. ! d (val), var (error var), id (iiceconc etc.), lon, lat, depths,
  59. ! ipic, jpic (i/j-pivot point in grid), ns (support, 0: point meas.),
  60. ! a1-4 (bilin. coeff), status (not used)
  61. narg= iargc()
  62. PRINT *, narg
  63. if (narg<=1) then
  64. ! Write info
  65. write(*,*)
  66. write(*,*) " (prep_obs) takes a real obs, extracts the desired variable and outputs"
  67. write(*,*) " it in a format that the EnKF can read & treat ('observations.uf')."
  68. write(*,*)
  69. write(*,*) " A file named mask.nc containing the variables tmaskutil, nav_lon and nav_lat"
  70. write(*,*) " is expected to be in the current directory (ORCA-file)"
  71. write(*,*)
  72. write(*,*) " Three command line arguments are expected:"
  73. write(*,*) " 1. Path to the nc file of which the data is to be extracted."
  74. write(*,*) " 2. Variable name that can be found in there, 'h_i_htc1' or"
  75. write(*,*) " 'sic'. or dxdy_ice"
  76. write(*,*) " 3. A tag with the date, e.g. 19790520"
  77. write(*,*)
  78. write(*,*) " Hope to see you again soon."
  79. write(*,*)
  80. stop "(prep_obs): Stopped."
  81. end if
  82. ! Command line arguments
  83. nvar=narg-1
  84. ! Get them
  85. call getarg(1, dummy)
  86. read(dummy,'(A)') filename
  87. ! Some parameter checking
  88. inquire(file=trim(filename), exist=ex)
  89. if (.not.ex) then
  90. print *, '(p_prep_obs): file does not exist: '// trim(filename)
  91. stop
  92. end if
  93. ! Get mask, lons & lats
  94. ! open the maskfile
  95. ierr = nf90_open(trim(maskfile),nf90_NoWrite,ncid); if (ierr.ne.nf90_noerr) call handle_err(ierr, "opening mask file")
  96. ! Find VarID of tmask & get values
  97. ierr = nf90_inq_varid(ncid, 'tmask', varID); if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID tmask")
  98. ierr = nf90_get_var(ncid, varID, mask) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable tmaks")
  99. ! Find VarID of longitude & get vals
  100. ierr = nf90_inq_varid(ncid, 'nav_lon', varID); if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID nav_lon")
  101. ierr = nf90_get_var(ncid, varID, lons) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable nav_lon")
  102. ! Find VarID of latitude & get vals
  103. ierr = nf90_inq_varid(ncid, 'nav_lat', varID); if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID nav_lat")
  104. ierr = nf90_get_var(ncid, varID, lats) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable nav_lat")
  105. ! Close maskfile
  106. ierr = nf90_close(ncid)
  107. if (ierr.ne.nf90_noerr) call handle_err(ierr, "closing")
  108. allocate( obs(nvar*sum(mask)), STAT=ierr )
  109. if (ierr.ne.0) call handle_err(ierr, "allocating obs") !no netcdf-error! ohwell.
  110. obscnt = 0 ! Keeps track of nr of added obs's.
  111. call getarg(2, dummy)
  112. read(dummy,'(A)') varname
  113. call getarg(3, dummy)
  114. read(dummy,*) rdate
  115. IF ( TRIM(varname) .eq. 'rfb' ) THEN
  116. WRITE(*,*) "Handling ", TRIM(varname)
  117. ! Min and max values for error used to screen the data (any data with
  118. ! standard deviation in between those values will be selected
  119. obsmin = 0.01
  120. obsmax = 10.0
  121. errmin = 0.01
  122. errmax = 1.0
  123. ELSEIF ( TRIM(varname) .eq. 'vt_i') THEN
  124. WRITE(*,*) "Handling ", TRIM(varname)
  125. obsmin = 0.01
  126. obsmax = 50.0
  127. errmin = 0.01
  128. errmax = 1.0
  129. ELSEIF ( TRIM(varname) .eq. 'at_i') THEN
  130. WRITE(*,*) "Handling ", TRIM(varname)
  131. obsmin = 0.001
  132. obsmax = 1.0
  133. errmin = 0.001
  134. errmax = 0.5
  135. ELSEIF ( TRIM(varname) .eq. 'siconc') THEN
  136. WRITE(*,*) "Handling ", TRIM(varname)
  137. obsmin = 0.001
  138. obsmax = 1.0
  139. errmin = 0.001
  140. errmax = 0.5
  141. ELSE
  142. WRITE(*,*) " (prep_obs) Currently only the variables rfb (sea ice freeboard),"
  143. WRITE(*,*) " vt_i (total sea ice volume)"
  144. WRITE(*,*) " at_i (total sea ice concentration)"
  145. WRITE(*,*) " can be processed "
  146. STOP "(prep_obs) Aborted"
  147. ENDIF
  148. ! User info
  149. WRITE(*,*) "(prep_obs) Extracting "//TRIM(varname)//" from "//TRIM(filename)
  150. ! Some preparations
  151. obsfld=0.
  152. ! open the netCDF file
  153. ierr = nf90_open(trim(filename),nf90_NoWrite,ncid) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "opening")
  154. ! Read observation data
  155. ! Find VarID of varname
  156. ierr = nf90_inq_varid(ncid, trim(varname), varID) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID")
  157. ! get values
  158. ierr = nf90_get_var(ncid, varID, obsfld) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable")
  159. ! Find VarID of varname
  160. cvarerror=TRIM(varname)//'_sd'
  161. ierr = nf90_inq_varid(ncid, cvarerror, varID) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "inquiring varID")
  162. ! get values
  163. ierr = nf90_get_var(ncid, varID, errorfld) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "getting variable")
  164. ! Close file
  165. ierr = nf90_close(ncid) ; if (ierr.ne.nf90_noerr) call handle_err(ierr, "closing")
  166. ! User info - ADAPT ACCORDINGLY
  167. WRITE (*,*) " (prep_obs) Using data >40N and <45S"
  168. ! Loop over space
  169. DO i = 1, SIZE(mask, 1)
  170. DO j = 1, SIZE(mask, 2)
  171. ! Pick out ocean points where data is existent
  172. IF ( (errorfld(i,j)) ** 2 .GT. errmin ** 2 &
  173. .AND. (errorfld(i,j)) ** 2 .LT. errmax ** 2 &
  174. .AND. obsfld(i,j) .GT. obsmin &
  175. .AND. obsfld(i,j) .LT. obsmax ) THEN
  176. ! Limit 'obs' spatially
  177. IF ( ( lats(i,j) .GE. latmin_NH &
  178. .AND. lats(i,j) .LE. latmax_NH ) &
  179. .OR.( lats(i,j) .LE. (-latmin_SH) &
  180. .AND. lats(i,j) .GE. (-latmax_SH) ) &
  181. ) THEN
  182. obscnt = obscnt + 1
  183. obs(obscnt)%d = obsfld(i,j)
  184. obs(obscnt)%lon = lons(i,j)
  185. obs(obscnt)%lat = lats(i,j)
  186. obs(obscnt)%ipiv = i ! the i-point of the grid of the model
  187. obs(obscnt)%jpiv = j ! the j-point of the grid of the model
  188. ! Put other data into obs type array
  189. obs(obscnt)%var = (errorfld(i,j))**2 ! The variance
  190. obs(obscnt)%id = TRIM(varname)
  191. obs(obscnt)%depths = 0
  192. obs(obscnt)%ns = 0
  193. obs(obscnt)%a1 = 1
  194. obs(obscnt)%a2 = 0
  195. obs(obscnt)%a3 = 0
  196. obs(obscnt)%a4 = 0
  197. obs(obscnt)%status = .TRUE.
  198. obs(obscnt)%i_orig_grid = -1
  199. obs(obscnt)%j_orig_grid = -1
  200. obs(obscnt)%h = -1.0
  201. obs(obscnt)%date = rdate
  202. obs(obscnt)%orig_id = 0
  203. END IF ! Spatial selection
  204. END IF ! if valid point
  205. END DO ! j
  206. END DO ! i
  207. !Write data:
  208. INQUIRE(iolength=i)obs(1)
  209. OPEN (11, file='observations.uf', status='replace',form='unformatted', access='direct', recl=i)
  210. DO j = 1, obscnt
  211. WRITE(11, rec=j)obs(j)
  212. ENDDO
  213. CLOSE(11)
  214. ! Write data to textfile, for visual checking
  215. OPEN(12, file = 'observations.txt')
  216. DO j = 1, obscnt
  217. WRITE(12, FMT = 53) obs(j)
  218. 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)
  219. ENDDO
  220. CLOSE(12)
  221. ! Tell user how many obs there were
  222. WRITE(*,*) " (prep_obs) Wrote out this many obs: ", obscnt
  223. WRITE(*,*) " (prep_obs) Number of ocean points : ", sum(mask)
  224. ! Cleanup
  225. IF (allocated(obs)) deallocate(obs,STAT=ierr)
  226. WRITE(*,*) ' (prep_obs) End successfully reached'; WRITE(*,*)
  227. contains
  228. subroutine handle_err(status, infomsg)
  229. integer, intent ( in) :: status
  230. character(len = *), intent ( in), optional :: infomsg
  231. if(status /= nf90_noerr) then
  232. if (present(infomsg)) then
  233. print *, '(prep_obs) Error while '//infomsg//' - '//trim(nf90_strerror(status))
  234. else
  235. print *, trim(nf90_strerror(status))
  236. endif ! opt arg
  237. print *,'(prep_obs)'
  238. stop
  239. end if ! check error status
  240. end subroutine handle_err
  241. end program prep_obs