m_obs.F90 11 KB

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