m_obs.F90 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378
  1. ! File: m_obs.F90
  2. !
  3. ! Created: 6 Feb 2008
  4. !
  5. ! Last modified: 21 Feb 2008
  6. !
  7. ! Author: Pavel Sakov*
  8. ! NERSC
  9. !
  10. ! Purpose: Generic code to deal with observations.
  11. !
  12. ! Description: This module contains the following functions and subroutines:
  13. ! - obs_setobs
  14. ! reads the observations into allocatable array obs(nobs)
  15. ! of type(measurement)
  16. ! - obs_prepareobs
  17. ! conducts state-dependent pre-processing of observations
  18. ! - obs_prepareuobs
  19. ! conducts state-dependent pre-processing of observations
  20. ! of a given type
  21. ! It also contains the following data:
  22. ! - obs
  23. ! allocatable array of type(measurement)
  24. ! - nobs
  25. ! number of observations (may differ from the size of the
  26. ! array)
  27. !
  28. ! * This file contains some modified code of unknown origin
  29. ! from EnKF package. In particular, the code here supersedes
  30. ! the code from:
  31. ! m_get_nrobs_d.F90
  32. ! m_get_obs_d.F90
  33. !
  34. ! Modifications:
  35. ! 09/11/2012 Geir Arne Waagbo:
  36. ! -- Added support for OSISAF ice drift obs
  37. ! 29/07/2010 PS:
  38. ! -- modified obs_QC(). The maximal increment now does not go to
  39. ! 0 as the innovation increases, but rather is limited by
  40. ! KMAX * sigma_ens
  41. ! 29/06/2010 PS:
  42. ! -- added obs_QC()
  43. ! 26/02/2008 PS:
  44. ! -- put "obs" and "nobs" as public data in this module
  45. module m_obs
  46. #if defined (QMPI)
  47. use qmpi
  48. #else
  49. use qmpi_fake
  50. #endif
  51. use mod_measurement
  52. use m_uobs
  53. use m_insitu
  54. implicit none
  55. !
  56. ! public stuff
  57. !
  58. integer, public :: nobs = -1
  59. type(measurement), allocatable, dimension(:), public :: obs
  60. public obs_readobs
  61. public obs_prepareobs
  62. public obs_QC
  63. !
  64. ! private stuff
  65. !
  66. private obs_testrange
  67. integer, parameter, private :: STRLEN = 512
  68. real, parameter, private :: TEM_MIN = -2.0d0
  69. real, parameter, private :: TEM_MAX = 50.0d0
  70. real, parameter, private :: SAL_MIN = 2.0d0
  71. real, parameter, private :: SAL_MAX = 40.0d0
  72. real, parameter, private :: SSH_MIN = -3.0d0
  73. real, parameter, private :: SSH_MAX = 3.0d0
  74. real, parameter, private :: ICEC_MIN = 0.0d0
  75. real, parameter, private :: ICEC_MAX = 0.999d0 ! [FM] Changed from 0.996 to 0.999
  76. real, parameter, private :: RFB_MIN = 0.0d0 ! FM 2020
  77. real, parameter, private :: RFB_MAX = 10.0d0
  78. real, parameter, private :: VT_I_MIN = 0.0d0 ! FM 2020
  79. real, parameter, private :: VT_I_MAX = 10.0d0
  80. real, parameter, private :: UVICE_MIN = -100.0
  81. real, parameter, private :: UVICE_MAX = 100.0
  82. private obs_prepareuobs, obs_realloc
  83. contains
  84. ! Obtain observations to be used for assimilation from the file
  85. ! "observation.uf". Store the number of observations in "nobs" and the data
  86. ! in the array "obs".
  87. !
  88. subroutine obs_readobs
  89. use m_parameters
  90. logical :: exists = .false.
  91. type(measurement) :: record
  92. integer :: rsize
  93. integer :: ios
  94. integer :: o
  95. CHARACTER(LEN=*), PARAMETER :: &
  96. FMT2 = "(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)"
  97. real :: myX
  98. real :: myY
  99. !========== TEST
  100. ! inquire(iolength = rsize) record
  101. ! !open(10, file = 'test.txt', form = 'unformatted',&
  102. ! ! access = 'direct', recl = rsize, status = 'old')
  103. ! allocate(obs(2))
  104. ! open(10, file = 'observations.txt')!, form = 'unformatted',&
  105. ! !access = 'direct', recl = rsize, status = 'old')
  106. ! !read(10, *) obs(1)
  107. !
  108. ! do o = 1, 2
  109. ! read(10, *) obs(o)
  110. ! PRINT *, obs(o)
  111. ! end do
  112. ! close(10)
  113. ! stop
  114. !==========
  115. if (nobs >= 0) then
  116. return
  117. end if
  118. ! Testing existence of file
  119. inquire(file = 'observations.txt', exist = exists)
  120. !inquire(file = 'observations.uf', exist = exists)
  121. if (.not. exists) then
  122. if (master) then
  123. print *, 'ERROR: obs_getnobs(): file "observations.txt" does not exist'
  124. end if
  125. stop
  126. end if
  127. inquire(iolength = rsize) record
  128. open(10, file = 'observations.txt')!, form = 'unformatted',&
  129. ! EXPERIMENTAL
  130. !open(10, file = 'observations.uf', form = 'unformatted',&
  131. ! access = 'direct', recl = rsize, status = 'old')!, form = 'unformatted',&
  132. !access = 'direct', recl = rsize, status = 'old')
  133. ! END EXPERIMENTAL
  134. ! I guess there is no other way to work out the length other than read the
  135. ! file in fortran - PS
  136. !
  137. o = 1
  138. do while (.true.)
  139. read(10, *, iostat = ios) record
  140. if (ios /= 0) then
  141. nobs = o - 1
  142. exit
  143. end if
  144. o = o + 1
  145. enddo
  146. allocate(obs(nobs))
  147. ! PS - there were problem with using rewind(): g95 reported:
  148. ! "Cannot REWIND a file opened for DIRECT access". Therefore reopen.
  149. !
  150. close(10)
  151. open(10, file = 'observations.txt')!, form = 'unformatted',&
  152. ! BEGIN EXPERIMENTAL
  153. !open(10, file = 'observations.uf', form = 'unformatted',&
  154. ! access = 'direct', recl = rsize, status = 'old')
  155. ! -- END EXPERIMENTAL
  156. do o = 1, nobs
  157. read(10, *) obs(o)
  158. call ucase(obs(o) % id)
  159. !PRINT *, obs(o)
  160. enddo
  161. close(10)
  162. if (RFACTOR1 /= 1.0d0) then
  163. do o = 1, nobs
  164. obs(o) % var = obs(o) % var * RFACTOR1
  165. end do
  166. end if
  167. call uobs_get(obs % id, nobs, master)
  168. call obs_testrange
  169. end subroutine obs_readobs
  170. subroutine obs_testrange
  171. integer :: o, uo, nbad
  172. real :: dmin, dmax
  173. if (master) then
  174. print '(a)', ' EnKF: testing range for each type of obs '
  175. end if
  176. do uo = 1, nuobs
  177. if (trim(unique_obs(uo)) == 'SST' .or. trim(unique_obs(uo)) == 'TEM'&
  178. .or. trim(unique_obs(uo)) == 'GTEM') then
  179. dmin = TEM_MIN
  180. dmax = TEM_MAX
  181. elseif (trim(unique_obs(uo)) == 'SAL'&
  182. .or. trim(unique_obs(uo)) == 'GSAL') then
  183. dmin = SAL_MIN
  184. dmax = SAL_MAX
  185. elseif (trim(unique_obs(uo)) == 'SLA'&
  186. .or. trim(unique_obs(uo)) == 'TSLA'&
  187. .or. trim(unique_obs(uo)) == 'SSH') then
  188. dmin = SSH_MIN
  189. dmax = SSH_MAX
  190. elseif (trim(unique_obs(uo)) == 'ICEC') then
  191. dmin = ICEC_MIN
  192. dmax = ICEC_MAX
  193. elseif (trim(unique_obs(uo)) == 'AT_I') then ! [FM] Added as we assimilate total ice conc. (opposed to indiv. category
  194. dmin = ICEC_MIN
  195. dmax = ICEC_MAX
  196. elseif (trim(unique_obs(uo)) == 'RFB') then ! FM added 2020
  197. dmin = RFB_MIN
  198. dmax = RFB_MAX
  199. elseif (trim(unique_obs(uo)) == 'VT_I') then ! FM added 2021
  200. dmin = VT_I_MIN
  201. dmax = VT_I_MAX
  202. elseif (trim(unique_obs(uo)) == 'V_ICE'&
  203. .or. trim(unique_obs(uo)) == 'U_ICE') then
  204. dmin = UVICE_MIN
  205. dmax = UVICE_MAX
  206. elseif (trim(unique_obs(uo)) == 'U2D_I'& ! [FM] OSISAF 2-day sea ice drift converted to m/s and interpolated onto ORCA
  207. .OR. trim(unique_obs(uo)) == 'V2D_I') THEN
  208. dmin = UVICE_MIN
  209. dmax = UVICE_MAX
  210. elseif ((index(trim(unique_obs(uo)),'DX') .gt. 0) &
  211. .or. (index(trim(unique_obs(uo)),'DY') .gt. 0)) then
  212. ! The type can be DX1,DX2,..,DX5,DY1,..DY5
  213. dmin = UVICE_MIN
  214. dmax = UVICE_MAX
  215. else
  216. dmin = -1.0d6
  217. dmax = 1.0d6
  218. print *, 'ERROR: obs_testrange(): "', trim(unique_obs(uo)), '": unknown type'
  219. stop
  220. end if
  221. nbad = 0
  222. do o = uobs_begin(uo), uobs_end(uo)
  223. if (obs(o) % status .and.&
  224. (obs(o) % d < dmin .or. obs(o) % d > dmax)) then
  225. obs(o) % status = .false.
  226. nbad = nbad + 1
  227. end if
  228. end do
  229. if (master) then
  230. print '(a, a, a, i6, a)', ' ', trim(unique_obs(uo)), ': ', nbad, ' outliers'
  231. end if
  232. end do
  233. if (master) then
  234. print *
  235. end if
  236. end subroutine obs_testrange
  237. ! Prepare observations before allocating matrices S, D, and A in EnKF().
  238. ! This invloves mainly thinning, superobing, or sorting.
  239. !
  240. ! Note that generically this processing can not be completely outsourced
  241. ! to the preprocessing stage, at least for in-situ data, because its thinning
  242. ! involves reading ensemble members for layer depth information.
  243. !
  244. subroutine obs_prepareobs()
  245. implicit none
  246. integer :: iuobs
  247. if (master) then
  248. print '(a)', ' EnKF: preparing observations'
  249. end if
  250. do iuobs = 1, nuobs
  251. call obs_prepareuobs(trim(unique_obs(iuobs)))
  252. end do
  253. ! calculate again the number of observation of each type (that could change
  254. ! in prepare_obs)
  255. call uobs_get(obs % id, nobs, master)
  256. end subroutine obs_prepareobs
  257. ! Prepare (thin, superob) observations of type "obstag".
  258. !
  259. subroutine obs_prepareuobs(obstag)
  260. character(*), intent(in) :: obstag
  261. character(STRLEN) :: fname
  262. if (trim(obstag) == 'SAL' .or. trim(obstag) == 'TEM' .or.&
  263. trim(obstag) == 'GSAL' .or. trim(obstag) == 'GTEM') then
  264. call insitu_prepareobs(trim(obstag), nobs, obs)
  265. if (master) then
  266. write(fname, '(a, ".nc")') trim(obstag)
  267. print *, 'Writing "', trim(obstag), '" obs to be assimilated to "',&
  268. trim(fname), '"'
  269. call insitu_writeprofiles(fname, trim(obstag), nobs, obs);
  270. end if
  271. else
  272. ! do nothing for obs of other types for now
  273. end if
  274. call obs_realloc
  275. end subroutine obs_prepareuobs
  276. subroutine obs_realloc()
  277. type(measurement), allocatable :: newobs(:)
  278. if (nobs < 0 .or. nobs == size(obs)) then
  279. return
  280. end if
  281. allocate(newobs(nobs))
  282. newobs = obs(1 : nobs)
  283. deallocate(obs)
  284. allocate(obs(nobs))
  285. obs = newobs
  286. deallocate(newobs)
  287. end subroutine obs_realloc
  288. subroutine obs_QC(m, S)
  289. use m_parameters
  290. implicit none
  291. integer :: m
  292. real :: S(nobs, m)
  293. integer :: nmodified(nuobs)
  294. real :: so(m), smean, svar
  295. integer :: o, uo
  296. real :: ovar, inn, newovar
  297. if (master) then
  298. print *, 'Starting generic observation QC'
  299. end if
  300. nmodified = 0
  301. do uo = 1, nuobs
  302. do o = uobs_begin(uo), uobs_end(uo)
  303. so = S(o, :);
  304. smean = sum(so) / m ! must be 0...
  305. svar = sum((so - smean) ** 2) / real(m - 1)
  306. ovar = obs(o) % var
  307. inn = obs(o) % d - smean
  308. obs(o) % var = sqrt((svar + ovar) ** 2 +&
  309. svar * (inn / KFACTOR) ** 2) - svar
  310. if (svar > 0 .and. obs(o) % var / ovar > 2.0d0) then
  311. nmodified(uo) = nmodified(uo) + 1
  312. end if
  313. end do
  314. end do
  315. if (master) then
  316. do uo = 1, nuobs
  317. print *, ' ', trim(unique_obs(uo)), ':'
  318. print *, ' # of observations:', uobs_end(uo) - uobs_begin(uo) + 1
  319. print *, ' (of them) substantially modified:', nmodified(uo)
  320. end do
  321. end if
  322. end subroutine obs_QC
  323. end module m_obs