m_prep_4_EnKF.F90 29 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765
  1. ! File: m_prep_4_EnKF.F90
  2. !
  3. ! Created: ???
  4. !
  5. ! Last modified: 29/06/2010
  6. !
  7. ! Purpose: Calculation of HA ("S")
  8. !
  9. ! Description: Calculates HA by going sequentially through each data type.
  10. !
  11. ! Modifications:
  12. ! 09/11/2012 Geir Arne Waagbo:
  13. ! - Added support for OSISAF ice drift obs
  14. ! 29/07/2010 PS:
  15. ! - merged insitu_QC() with generic obs_QC(). Moved
  16. ! insitu_writeforecast() to the point after QC.
  17. ! 29/06/2010 PS:
  18. ! - added generic observation QC: increase the observation
  19. ! error when observation and ensemble mean are much too far
  20. ! away than expected
  21. ! Prior history:
  22. ! Not documented.
  23. module m_prep_4_EnKF
  24. integer, parameter, private :: STRLEN = 512
  25. private read_mean_ssh
  26. contains
  27. ! This subroutine uses the observation and ensembles from the model
  28. ! to prepare the input to the EnKF analysis scheme.
  29. ! The output from this routine is used directly in the global analysis
  30. ! while the output has to be run through a "filter" to be used in the
  31. ! local analysis scheme.
  32. ! S = HA (ensemble observation anomalies)
  33. ! d = d - Hx (innovations)
  34. !
  35. ! S is calculated in two steps:
  36. ! 1. S = HE
  37. ! 2. S = S - repmat(s, 1, m),
  38. ! where s = mean(S')';
  39. ! Note that in reality (with HYCOM) H is different for each member...
  40. ! So that HX must be read "HX" rather than "H * X".
  41. !
  42. subroutine prep_4_EnKF(nrens, enslist, d, S, depths, meandx, nx, ny, nz)
  43. #if defined (QMPI)
  44. use qmpi, only : master, stop_mpi
  45. #else
  46. use qmpi_fake, only : master, stop_mpi
  47. #endif
  48. use mod_measurement
  49. use m_obs
  50. use m_Generate_element_Si
  51. use m_get_mod_fld
  52. use m_read_icemod
  53. use m_parameters
  54. implicit none
  55. integer, intent(in) :: nx, ny, nz ! Model size
  56. integer, intent(in) :: nrens ! Size of ensemble
  57. integer, dimension(:),intent(in) :: enslist ! [CKB,FM] List of existing ens members
  58. real, intent(in) :: depths(nx, ny)
  59. real, intent(in) :: meandx ! mean grid size
  60. real, intent(inout) :: d(nobs)
  61. real, intent(inout) :: S(nobs, nrens)
  62. real :: x(nobs)
  63. integer :: i, j, m, iens
  64. real*4, dimension(nx,ny) :: fldr4
  65. real :: readfld(nx, ny), ai1(nx,ny), ai2(nx,ny), ai3(nx,ny), ai4(nx,ny), ai5(nx,ny), uice(nx,ny), vice(nx,ny)
  66. real :: vi1(nx,ny), vi2(nx,ny), vi3(nx,ny), vi4(nx,ny), vi5(nx,ny)
  67. real :: vs1(nx,ny), vs2(nx,ny), vs3(nx,ny), vs4(nx,ny), vs5(nx,ny)
  68. ! hard-coded for now
  69. integer, parameter :: drnx = 152, drny = 132
  70. real*4, dimension(drnx, drny) :: modzon, modmer
  71. integer, parameter :: drnx_osisaf = 119, drny_osisaf = 177
  72. real*4, dimension(drnx_osisaf, drny_osisaf) :: dX, dY
  73. integer :: reclSLA, ios, reclDRIFT
  74. character*3 :: cmem
  75. character*2 :: day
  76. character*1 :: offset
  77. logical :: ex
  78. character(STRLEN) :: fname
  79. integer :: iuobs
  80. ! FANF: For track assim we launch m_Generate_Si for each day (t=1:Wd)
  81. ! which fills in S at the approriate indices.
  82. ! Wd is is the assimilation cycle (i.e. 7 days)
  83. !
  84. integer :: Wd, t
  85. integer :: tlevel
  86. real :: field2(nx, ny), field3(nx, ny) ! auxiliary fields (e.g. mean SSH,
  87. ! field bias estimate etc.)
  88. integer :: nthisobs, thisobs(nobs)
  89. if (any(obs(:) % id == 'TSLA ')) then
  90. Wd = 6
  91. else
  92. Wd = 0
  93. endif
  94. ! security check
  95. !
  96. if (any(obs(:) % id == 'SSH ') .or. any(obs(:) % id == 'SLA ')) then
  97. if (any(obs(:) % id == 'SLA ')) then
  98. inquire(exist = ex, file = 'model_SLA.uf')
  99. if (.not.ex) then
  100. if (master) print *,'model_SLA.uf does not exist'
  101. call stop_mpi()
  102. end if
  103. end if
  104. if (any(obs(:) % id == 'SSH ')) then
  105. inquire(exist = ex, file = 'model_SSH.uf')
  106. if (.not.ex) then
  107. if (master) print *,'model_SSH.uf does not exist'
  108. call stop_mpi()
  109. end if
  110. end if
  111. end if
  112. ! construct S=HA
  113. !
  114. do iuobs = 1, nuobs
  115. if (master) then
  116. print *, 'prep_4_EnKF: now preparing "', trim(unique_obs(iuobs)), '" observations'
  117. end if
  118. if (trim(unique_obs(iuobs)) == 'ICEC') then
  119. do iens = 1, nrens
  120. write(cmem,'(i3.3)') iens
  121. tlevel = 1
  122. call get_mod_fld_new(trim('forecast'//cmem), readfld, iens,&
  123. 'icec', 0, tlevel, nx, ny)
  124. if (tlevel == -1) then
  125. if (master) then
  126. print *, 'ERROR: get_mod_fld_new(): failed for "icec"'
  127. end if
  128. stop
  129. end if
  130. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  131. readfld, depths, nx, ny, nz, 0)
  132. end do
  133. ! [FM, May 2013: for LIM3 sea ice model]
  134. elseif (trim(unique_obs(iuobs)) == 'AT_I') then
  135. do iens = 1, nrens
  136. write(cmem,'(i3.3)') iens
  137. tlevel = 1
  138. call io_mod_fld(ai1,iens,enslist, &
  139. 'a_i_htc1', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  140. call io_mod_fld(ai2,iens,enslist, &
  141. 'a_i_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  142. call io_mod_fld(ai3,iens,enslist, &
  143. 'a_i_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  144. call io_mod_fld(ai4,iens,enslist, &
  145. 'a_i_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  146. call io_mod_fld(ai5,iens,enslist, &
  147. 'a_i_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  148. if (tlevel == -1) then
  149. if (master) then
  150. print *, 'ERROR: io_mod_fld_new(): failed for "at_i"'
  151. end if
  152. stop
  153. end if
  154. ! Multipply by 100 to match obs conventions
  155. readfld=(ai1+ai2+ai3+ai4+ai5) * 100.0
  156. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  157. readfld, depths, nx, ny, nz, 0)
  158. end do
  159. ! [AD] March 2024: for NEMO-SI3 first test assim. siconc from nemo itself]
  160. elseif (trim(unique_obs(iuobs)) == 'SICON') then
  161. do iens = 1, nrens
  162. write(cmem,'(i3.3)') iens
  163. tlevel = 1
  164. call io_mod_fld(ai1,iens,enslist, &
  165. 'a_i', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  166. !call io_mod_fld(ai2,iens,enslist, &
  167. ! 'a_i_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  168. !call io_mod_fld(ai3,iens,enslist, &
  169. ! 'a_i_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  170. !call io_mod_fld(ai4,iens,enslist, &
  171. ! 'a_i_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  172. !call io_mod_fld(ai5,iens,enslist, &
  173. ! 'a_i_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  174. if (tlevel == -1) then
  175. if (master) then
  176. print *, 'ERROR: io_mod_fld_new(): failed for "sicon"'
  177. end if
  178. stop
  179. end if
  180. ! Multipply by 100 to match obs conventions
  181. !readfld=(ai1+ai2+ai3+ai4+ai5) * 100.0
  182. readfld=(ai1) * 100.0
  183. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  184. readfld, depths, nx, ny, nz, 0)
  185. end do
  186. ! [AD] April 2024: for NEMO-SI3 first test assim. sic from OSI-SAF interpoled on eORCA1]
  187. elseif (trim(unique_obs(iuobs)) == 'SIC') then
  188. do iens = 1, nrens
  189. write(cmem,'(i3.3)') iens
  190. tlevel = 1
  191. call io_mod_fld(ai1,iens,enslist, &
  192. 'a_i', 1, 1, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  193. call io_mod_fld(ai2,iens,enslist, &
  194. 'a_i', 1, 2, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  195. call io_mod_fld(ai3,iens,enslist, &
  196. 'a_i', 1, 3, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  197. call io_mod_fld(ai4,iens,enslist, &
  198. 'a_i', 1, 4, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  199. call io_mod_fld(ai5,iens,enslist, &
  200. 'a_i', 1, 5, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  201. if (tlevel == -1) then
  202. if (master) then
  203. print *, 'ERROR: io_mod_fld_new(): failed for "SIC"'
  204. end if
  205. stop
  206. end if
  207. ! Multipply by 100 to match obs conventions
  208. readfld=(ai1+ai2+ai3+ai4+ai5) * 100.0
  209. ! readfld=(ai1) * 100.0
  210. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  211. readfld, depths, nx, ny, nz, 0)
  212. end do
  213. ! freeboard
  214. elseif(trim(unique_obs(iuobs)) == 'VT_I') then
  215. do iens = 1, nrens
  216. write(cmem, '(i3.3)') iens
  217. tlevel = 1
  218. call io_mod_fld(ai1,iens,enslist, &
  219. 'a_i_htc1', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  220. call io_mod_fld(ai2,iens,enslist, &
  221. 'a_i_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  222. call io_mod_fld(ai3,iens,enslist, &
  223. 'a_i_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  224. call io_mod_fld(ai4,iens,enslist, &
  225. 'a_i_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  226. call io_mod_fld(ai5,iens,enslist, &
  227. 'a_i_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  228. call io_mod_fld(vi1,iens,enslist, &
  229. 'v_i_htc1', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  230. call io_mod_fld(vi2,iens,enslist, &
  231. 'v_i_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  232. call io_mod_fld(vi3,iens,enslist, &
  233. 'v_i_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  234. call io_mod_fld(vi4,iens,enslist, &
  235. 'v_i_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  236. call io_mod_fld(vi5,iens,enslist, &
  237. 'v_i_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  238. if (tlevel == -1) then
  239. if (master) then
  240. print *, 'ERROR: io_mod_fld_nex(): failed for "SIFB"'
  241. end if
  242. stop
  243. end if
  244. readfld=(vi1+vi2+vi3+vi4+vi5)
  245. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  246. readfld, depths, nx, ny, nz, 0)
  247. end do
  248. ! freeboard
  249. elseif(trim(unique_obs(iuobs)) == 'RFB') then
  250. do iens = 1, nrens
  251. write(cmem, '(i3.3)') iens
  252. tlevel = 1
  253. call io_mod_fld(ai1,iens,enslist, &
  254. 'a_i_htc1', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  255. call io_mod_fld(ai2,iens,enslist, &
  256. 'a_i_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  257. call io_mod_fld(ai3,iens,enslist, &
  258. 'a_i_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  259. call io_mod_fld(ai4,iens,enslist, &
  260. 'a_i_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  261. call io_mod_fld(ai5,iens,enslist, &
  262. 'a_i_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  263. call io_mod_fld(vi1,iens,enslist, &
  264. 'v_i_htc1', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  265. call io_mod_fld(vi2,iens,enslist, &
  266. 'v_i_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  267. call io_mod_fld(vi3,iens,enslist, &
  268. 'v_i_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  269. call io_mod_fld(vi4,iens,enslist, &
  270. 'v_i_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  271. call io_mod_fld(vi5,iens,enslist, &
  272. 'v_i_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  273. call io_mod_fld(vs1,iens,enslist, &
  274. 'v_s_htc1', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  275. call io_mod_fld(vs2,iens,enslist, &
  276. 'v_s_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  277. call io_mod_fld(vs3,iens,enslist, &
  278. 'v_s_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  279. call io_mod_fld(vs4,iens,enslist, &
  280. 'v_s_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  281. call io_mod_fld(vs5,iens,enslist, &
  282. 'v_s_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  283. if (tlevel == -1) then
  284. if (master) then
  285. print *, 'ERROR: io_mod_fld_nex(): failed for "SIFB"'
  286. end if
  287. stop
  288. end if
  289. readfld=(((vi1+vi2+vi3+vi4+vi5) * (1024.0 - 899.5) - 330 * (vs1+vs2+vs3+vs4+vs5)) / &
  290. 1024.0-0.25*(vs1 +vs2+vs3+vs4+vs5))
  291. !readfld=(((vi1+vi2+vi3+vi4+vi5) * (1024.0 - 899.5) - 330 * (vs1+vs2+vs3+vs4+vs5)) / 1024.0 - 0.25 * (vs1+vs2+vs3+vs4+vs5)) / (ai1+ai2+ai3+ai4+ai5)
  292. ! Conversion of models' sea ice thickness and snow thickness to
  293. ! model's freeboard using fixed densities for snow (330 kg/m3), ice
  294. ! (899.5 kg/m3 = average of MYI and FYI from Guerreiro et al. 2017
  295. ! and seawater (1024 kg/m3). The model freeboard is then lowered by
  296. ! 25% of the snow depth to account for the fact that the radar
  297. ! measurement underestimates the actual freeboard due to the lower
  298. ! propagation speed of the wave into the snow than in the air.
  299. ! Everything is converted from grid cell mean to in situ by
  300. ! dividing by concentration (if it is not zero). See exchanges
  301. ! e-mail with Sara Fleury 7 December 2020.
  302. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  303. readfld, depths, nx, ny, nz, 0)
  304. end do
  305. elseif (trim(unique_obs(iuobs)) == 'SST') then
  306. do iens = 1, nrens
  307. write(cmem,'(i3.3)') iens
  308. tlevel = 1
  309. call get_mod_fld_new(trim('forecast'//cmem), readfld, iens,&
  310. 'tn', 1, tlevel, nx, ny)
  311. PRINT *, "FRANCOIS"
  312. if (tlevel == -1) then
  313. if (master) then
  314. print *, 'ERROR: get_mod_fld_new(): failed for "SST"'
  315. end if
  316. stop
  317. end if
  318. if (prm_prmestexists('sstb')) then
  319. tlevel = 1
  320. call get_mod_fld_new(trim('forecast'//cmem), field2, iens,&
  321. 'sstb', 0, tlevel, nx, ny)
  322. if (tlevel == -1) then
  323. if (master) then
  324. print *, 'ERROR: get_mod_fld_new(): failed for "sstb"'
  325. end if
  326. stop
  327. end if
  328. readfld = readfld - field2
  329. end if
  330. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  331. readfld, depths, nx, ny, nz, 0)
  332. end do
  333. elseif (trim(unique_obs(iuobs)) == 'SLA' .or.&
  334. trim(unique_obs(iuobs)) == 'TSLA') then
  335. if (trim(unique_obs(iuobs)) == 'TSLA') then
  336. call read_mean_ssh(field2, nx, ny)
  337. end if
  338. inquire(iolength=reclSLA) fldr4
  339. ! FANF loop over each day of the week
  340. do t = 0, Wd
  341. if (trim(unique_obs(iuobs)) == 'TSLA') then
  342. write(day,'(i2.2)') t
  343. fname = trim('model_TSSH_'//day//'.uf')
  344. else
  345. fname = 'model_SLA.uf'
  346. endif
  347. if (master) then
  348. print *, 'TSLA, day', t, ': nobs = ',&
  349. count(obs(uobs_begin(iuobs) : uobs_end(iuobs)) % date == t)
  350. end if
  351. do iens = 1, nrens
  352. open(10, file = trim(fname), access = 'direct',&
  353. status = 'old', recl = reclSLA, action = 'read')
  354. read(10, rec = iens, iostat = ios) fldr4
  355. if (ios /= 0) then
  356. if (master) print *, 'Error reading ', trim(fname), ', member #', iens
  357. call stop_mpi()
  358. end if
  359. close(10)
  360. readfld = fldr4
  361. if (prm_prmestexists('msshb')) then
  362. write(cmem,'(i3.3)') iens
  363. tlevel = 1
  364. call get_mod_fld_new(trim('forecast'//cmem), field3, iens,&
  365. 'msshb', 0, tlevel, nx, ny)
  366. if (tlevel == -1) then
  367. if (master) then
  368. print *, 'ERROR: get_mod_fld_new(): failed for "msshb"'
  369. end if
  370. stop
  371. end if
  372. readfld = readfld - field3 ! mean SSH bias for this member
  373. end if
  374. if (trim(unique_obs(iuobs)) == 'TSLA') then
  375. readfld = readfld - field2 ! mean SSH
  376. end if
  377. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  378. readfld, depths, nx, ny, nz, t)
  379. end do
  380. if (master) then
  381. print *, 'forming S, day', t
  382. print *, ' # of non-zero ens observations = ', count(S /= 0.0)
  383. print *, ' # of zero ens observations = ', count(S == 0.0)
  384. print *, ' # of non-zero observations for member 1 = ', count(S(:, 1) /= 0.0)
  385. end if
  386. end do
  387. elseif (trim(unique_obs(iuobs)) == 'SAL' .or.&
  388. trim(unique_obs(iuobs)) == 'TEM' .or.&
  389. trim(unique_obs(iuobs)) == 'GSAL' .or.&
  390. trim(unique_obs(iuobs)) == 'GTEM') then
  391. if (master) then
  392. print *, ' Interpolating ensemble vectors to the locations of "',&
  393. trim(unique_obs(iuobs)), '" observations'
  394. end if
  395. !
  396. ! for each ensemble member process all profiles "in parallel",
  397. ! reading the fields layer by layer
  398. !
  399. do iens = 1, nrens
  400. call get_S(S(:, iens), trim(unique_obs(iuobs)), nobs, obs, iens)
  401. end do
  402. if (master) then
  403. print *, ' Interpolation completed'
  404. end if
  405. elseif ((trim(unique_obs(iuobs)) == 'U_ICE') .or. trim(unique_obs(iuobs)) == 'V_ICE') then
  406. do iens = 1, nrens
  407. ! [FM] Read the file
  408. !inquire(iolength=reclDRIFT) modzon, modmer
  409. !open(10, file = 'model_ICEDRIFT.uf', access = 'direct',&
  410. ! status = 'old', recl = reclDRIFT, action = 'read')
  411. !read(10, rec = iens, iostat = ios) modzon, modmer
  412. !close(10)
  413. !if (ios /= 0) then
  414. ! if (master) then
  415. ! print *,'ERROR: could not read ensemble ice drift for member ', iens
  416. ! end if
  417. ! call stop_mpi()
  418. !end if
  419. call io_mod_fld(uice,iens,enslist, &
  420. 'u_ice', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  421. call io_mod_fld(vice,iens,enslist, &
  422. 'v_ice', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  423. do m = 1, nobs
  424. i = obs(m) % i_orig_grid
  425. j = obs(m) % j_orig_grid
  426. if (trim(obs(m) % id) == 'U_ICE') then
  427. S(m, iens) = uice(i, j)
  428. elseif (trim(obs(m) % id) == 'V_ICE') then
  429. S(m, iens) = vice(i, j)
  430. end if
  431. end do
  432. end do
  433. elseif ((trim(unique_obs(iuobs)) == 'U2D_I') .OR. trim(unique_obs(iuobs)) == 'V2D_I' ) THEN
  434. ! ADDED BY FM FRANCOIS MASSONNET. u_ice_2d or v_ice_2d is the sea ice u or v-velocity
  435. ! obtained as follows:
  436. ! 1) Rotate OSISAF Low resolution 2-day sea ice drift in a {east,north}
  437. ! reference frame
  438. ! 2) Interpolate to the ORCA grid
  439. ! 3) Rotate to align with the ORCA grid
  440. ! 4) Multiply by 1000 and divide by 2*86400 to convert the 2-day
  441. ! displacement from km to m/s
  442. DO iens=1,nrens
  443. CALL read_icemod(uice,iens,enslist,'iicevelu',nx,ny)
  444. CALL read_icemod(vice,iens,enslist,'iicevelv',nx,ny)
  445. DO m = 1, nobs
  446. i = obs(m) % i_orig_grid
  447. j = obs(m) % j_orig_grid
  448. IF (trim(obs(m) % id) == 'U2D_I') THEN
  449. S(m,iens) = uice(i,j)
  450. ELSEIF (trim(obs(m) % id) == 'V2D_I') THEN
  451. S(m,iens) = vice(i,j)
  452. END IF
  453. END DO ! nobs
  454. END DO ! iens
  455. elseif ((index(unique_obs(iuobs),'DX') > 0 ) .or. (index(unique_obs(iuobs),'DY') > 0)) then
  456. ! OSISAF Ice drift observations (d-2-offset -> d-offset)
  457. !print *, 'Ice drift observation type: ', unique_obs(iuobs)
  458. offset = unique_obs(iuobs)(3:3)
  459. ! Use offset (1,2,3,4 or 5) to open correct model drift file
  460. inquire(iolength=reclDRIFT) dX, dY
  461. open(10, file = 'model_ICEDRIFT_OSISAF'//offset//'.uf', access = 'direct',&
  462. status = 'old', recl = reclDRIFT, action = 'read')
  463. do iens = 1, nrens
  464. read(10, rec = iens, iostat = ios) dX, dY
  465. if (ios /= 0) then
  466. if (master) then
  467. print *,'ERROR: could not read ensemble ice drift for member ', iens
  468. end if
  469. call stop_mpi()
  470. end if
  471. do m = 1, nobs
  472. i = obs(m) % i_orig_grid
  473. j = obs(m) % j_orig_grid
  474. if (index(obs(m)%id,'DX') > 0) then
  475. S(m, iens) = dX(i, j)
  476. elseif (index(obs(m)%id,'DY') > 0) then
  477. S(m, iens) = dY(i, j)
  478. end if
  479. end do
  480. end do
  481. close(10)
  482. else
  483. if (master) then
  484. print *,'ERROR: unknown obs type ' // trim(unique_obs(iuobs))
  485. end if
  486. call stop_mpi()
  487. end if
  488. end do ! iuobs
  489. ! some generic QC - relax fitting if the model and obs are too far apart
  490. !
  491. call obs_QC(nrens, S)
  492. ! add calculated HA to to observations-<type>.nc files for each data type
  493. !
  494. do iuobs = 1, nuobs
  495. if (master) then
  496. nthisobs = 0
  497. do m = 1, nobs
  498. if (trim(unique_obs(iuobs)) == trim(obs(m) % id)) then
  499. nthisobs = nthisobs + 1
  500. thisobs(nthisobs) = m
  501. end if
  502. end do
  503. ! add forecast values to the observation-<TYPE>.nc files
  504. !
  505. call add_forecast(unique_obs(iuobs), S(thisobs(1 : nthisobs), :), obs(thisobs(1 : nthisobs)))
  506. ! append the superobed values (and modified observation error
  507. ! variances) to the file with pre-processed observations (SAL.nc,
  508. ! TEM.nc, GSAL.nc or GTEM.nc)
  509. !
  510. if (trim(unique_obs(iuobs)) == 'SAL' .or.&
  511. trim(unique_obs(iuobs)) == 'TEM' .or.&
  512. trim(unique_obs(iuobs)) == 'GSAL' .or.&
  513. trim(unique_obs(iuobs)) == 'GTEM') then
  514. call insitu_writeforecast(unique_obs(iuobs), nobs, nrens, S, obs)
  515. end if
  516. end if
  517. end do
  518. if (master) then
  519. print *, 'm_prep_4_EnKF: end calculating S = HA'
  520. end if
  521. x = sum(S, DIM = 2) / real(nrens) ! [ FM ] The mean forecast interpolated in the obs.space
  522. if (master) print*,'m_prep_4_EnKF: end calculating Hx'
  523. if (master) then
  524. print *, 'Hx range = ', minval(x), '-', maxval(x)
  525. print *, 'mean(Hx) = ', sum(x) / real(nobs)
  526. end if
  527. if (master) then
  528. print *, 'observation range = ', minval(obs % d), '-', maxval(obs % d)
  529. print *, 'mean(observation) = ', sum(obs % d) / nobs
  530. end if
  531. ! Compute HA = HE - mean(HE)
  532. !
  533. if (master) print*,'prep_4_EnKF(): calculating S = S - x'
  534. do j = 1, nrens
  535. S(:, j) = S(:, j) - x ! [ FM ] This is really where we switch from actual model data to anomalies
  536. enddo
  537. ! Compute innovation
  538. !
  539. d = obs % d - x ! [ FM ] This is exactly was is also done in add_forecast. This is the mean innovation.
  540. if (master) then
  541. print *, ' innovation range = ', minval(d), '-', maxval(d)
  542. if (minval(d) < -1000.0d0) then
  543. print *, 'm_prep_4_EnKF: error: innovation too small detected'
  544. call stop_mpi()
  545. end if
  546. if (maxval(d) > 1000.0d0) then
  547. print *, 'm_prep_4_EnKF: error: innovation too big detected'
  548. call stop_mpi()
  549. end if
  550. end if
  551. end subroutine prep_4_EnKF
  552. subroutine read_mean_ssh(mean_ssh, nx, ny)
  553. #if defined (QMPI)
  554. use qmpi
  555. #else
  556. use qmpi_fake
  557. #endif
  558. use m_parameters
  559. integer, intent(in) :: nx, ny
  560. real, intent(out):: mean_ssh(nx, ny)
  561. logical :: exists
  562. inquire(file = trim(MEANSSHFNAME), exist = exists)
  563. if (.not. exists) then
  564. if (master) then
  565. print *,'ERROR: read_mean_ssh(): file "', trim(MEANSSHFNAME), '" not found'
  566. end if
  567. stop
  568. end if
  569. open (10, file = trim(MEANSSHFNAME), status = 'unknown',form = 'unformatted', action = 'read')
  570. read (10) mean_ssh
  571. close (10)
  572. end subroutine read_mean_ssh
  573. ! This subroutine adds forecast observations (i.e Hx) to the NetCDF
  574. ! observation files for each data type.
  575. !
  576. subroutine add_forecast(obstag, S, obs)
  577. use mod_measurement
  578. use nfw_mod
  579. implicit none
  580. character(OBSTYPESTRLEN), intent(in) :: obstag
  581. real, dimension(:, :), intent(in) :: S
  582. type(measurement), dimension(:) :: obs
  583. character(STRLEN) :: fname
  584. logical :: exists
  585. integer :: ncid
  586. integer :: dids(2), dimlen
  587. logical :: addsobs
  588. integer :: for_id, inn_id, forvar_id, slon_id, slat_id,&
  589. sdepth_id, sipiv_id, sjpiv_id, sd_id, svar_id,&
  590. newvar_id
  591. real, allocatable, dimension(:) :: x, Svar, innovation
  592. integer :: m, p, o
  593. write(fname, '(a, a, a)') 'observations-', trim(obstag), '.nc'
  594. inquire(file = trim(fname), exist = exists)
  595. if (.not. exists) then
  596. print *, 'file "', trim(fname), 'not found, skip adding forecast'
  597. return
  598. else
  599. print *, 'dumping forecast to "', trim(fname), '"'
  600. end if
  601. p = size(S, DIM = 1)
  602. m = size(S, DIM = 2)
  603. allocate(x(p), Svar(p), innovation(p))
  604. x = sum(S, DIM = 2) / real(m); ! [ FM the mean of S=HA ]
  605. Svar = 0.0
  606. do o = 1, p
  607. Svar(o) = sum((S(o, :) - x(o))** 2) ! [ FM thus each row of Svar is the variance (see below) of the forecast]
  608. end do
  609. Svar = Svar / real(m - 1)
  610. innovation = obs % d - x ! [ FM ] the innovation for the mean forecast (or mean of the innovation forecasts)
  611. addsobs = .false.
  612. call nfw_open(fname, nf_write, ncid)
  613. call nfw_inq_dimid(fname, ncid, 'nobs', dids(1))
  614. call nfw_inq_dimlen(fname, ncid, dids(1), dimlen)
  615. call nfw_redef(fname, ncid)
  616. if (dimlen == p) then
  617. dids(2) = dids(1)
  618. elseif (.not. nfw_dim_exists(ncid, 'nsobs')) then
  619. addsobs = .true.
  620. call nfw_def_dim(fname, ncid, 'nsobs', p, dids(2))
  621. call nfw_def_var(fname, ncid, 'slon', nf_float, 1, dids(2), slon_id)
  622. call nfw_def_var(fname, ncid, 'slat', nf_float, 1, dids(2), slat_id)
  623. call nfw_def_var(fname, ncid, 'sdepth', nf_float, 1, dids(2), sdepth_id)
  624. call nfw_def_var(fname, ncid, 'sipiv', nf_int, 1, dids(2), sipiv_id)
  625. call nfw_def_var(fname, ncid, 'sjpiv', nf_int, 1, dids(2), sjpiv_id)
  626. call nfw_def_var(fname, ncid, 'sd', nf_float, 1, dids(2), sd_id)
  627. call nfw_def_var(fname, ncid, 'svar', nf_float, 1, dids(2), svar_id)
  628. end if
  629. if (.not. nfw_var_exists(ncid, 'innovation')) then
  630. call nfw_def_var(fname, ncid, 'innovation', nf_double, 1, dids(2), inn_id)
  631. else
  632. call nfw_inq_varid(fname, ncid, 'innovation', inn_id)
  633. end if
  634. if (.not. nfw_var_exists(ncid, 'forecast')) then
  635. call nfw_def_var(fname, ncid, 'forecast', nf_double, 1, dids(2), for_id)
  636. else
  637. call nfw_inq_varid(fname, ncid, 'forecast', for_id)
  638. end if
  639. if (.not. nfw_var_exists(ncid, 'forecast_variance')) then
  640. call nfw_def_var(fname, ncid, 'forecast_variance', nf_double, 1, dids(2), forvar_id)
  641. else
  642. call nfw_inq_varid(fname, ncid, 'forecast_variance', forvar_id)
  643. end if
  644. if (.not. addsobs) then
  645. if (dimlen == p) then
  646. if (.not. nfw_var_exists(ncid, 'new_var')) then
  647. call nfw_def_var(fname, ncid, 'new_var', nf_double, 1, dids(2), newvar_id)
  648. else
  649. call nfw_inq_varid(fname, ncid, 'new_var', newvar_id)
  650. end if
  651. else
  652. if (.not. nfw_var_exists(ncid, 'new_svar')) then
  653. call nfw_inq_dimid(fname, ncid, 'nsobs', dids(2))
  654. call nfw_def_var(fname, ncid, 'new_svar', nf_double, 1, dids(2), newvar_id)
  655. else
  656. call nfw_inq_varid(fname, ncid, 'new_svar', newvar_id)
  657. end if
  658. end if
  659. end if
  660. call nfw_enddef(fname, ncid)
  661. call nfw_put_var_double(fname, ncid, forvar_id, Svar)
  662. call nfw_put_var_double(fname, ncid, for_id, x)
  663. call nfw_put_var_double(fname, ncid, inn_id, innovation)
  664. if (addsobs) then
  665. call nfw_put_var_double(fname, ncid, slon_id, obs % lon)
  666. call nfw_put_var_double(fname, ncid, slat_id, obs % lat)
  667. call nfw_put_var_double(fname, ncid, sdepth_id, obs % depth)
  668. call nfw_put_var_int(fname, ncid, sipiv_id, obs % ipiv)
  669. call nfw_put_var_int(fname, ncid, sjpiv_id, obs % jpiv)
  670. call nfw_put_var_double(fname, ncid, sd_id, obs % d)
  671. call nfw_put_var_double(fname, ncid, svar_id, obs % var)
  672. else
  673. call nfw_put_var_double(fname, ncid, newvar_id, obs % var)
  674. end if
  675. call nfw_close(fname, ncid)
  676. deallocate(x)
  677. deallocate(Svar)
  678. deallocate(innovation)
  679. end subroutine add_forecast
  680. end module m_prep_4_EnKF