m_prep_4_EnKF.F90 22 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597
  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. ! hard-coded for now
  67. integer, parameter :: drnx = 152, drny = 132
  68. real*4, dimension(drnx, drny) :: modzon, modmer
  69. integer, parameter :: drnx_osisaf = 119, drny_osisaf = 177
  70. real*4, dimension(drnx_osisaf, drny_osisaf) :: dX, dY
  71. integer :: reclSLA, ios, reclDRIFT
  72. character*3 :: cmem
  73. character*2 :: day
  74. character*1 :: offset
  75. logical :: ex
  76. character(STRLEN) :: fname
  77. integer :: iuobs
  78. ! FANF: For track assim we launch m_Generate_Si for each day (t=1:Wd)
  79. ! which fills in S at the approriate indices.
  80. ! Wd is is the assimilation cycle (i.e. 7 days)
  81. !
  82. integer :: Wd, t
  83. integer :: tlevel
  84. real :: field2(nx, ny), field3(nx, ny) ! auxiliary fields (e.g. mean SSH,
  85. ! field bias estimate etc.)
  86. integer :: nthisobs, thisobs(nobs)
  87. if (any(obs(:) % id == 'TSLA ')) then
  88. Wd = 6
  89. else
  90. Wd = 0
  91. endif
  92. ! security check
  93. !
  94. if (any(obs(:) % id == 'SSH ') .or. any(obs(:) % id == 'SLA ')) then
  95. if (any(obs(:) % id == 'SLA ')) then
  96. inquire(exist = ex, file = 'model_SLA.uf')
  97. if (.not.ex) then
  98. if (master) print *,'model_SLA.uf does not exist'
  99. call stop_mpi()
  100. end if
  101. end if
  102. if (any(obs(:) % id == 'SSH ')) then
  103. inquire(exist = ex, file = 'model_SSH.uf')
  104. if (.not.ex) then
  105. if (master) print *,'model_SSH.uf does not exist'
  106. call stop_mpi()
  107. end if
  108. end if
  109. end if
  110. ! construct S=HA
  111. !
  112. do iuobs = 1, nuobs
  113. if (master) then
  114. print *, 'prep_4_EnKF: now preparing "', trim(unique_obs(iuobs)), '" observations'
  115. end if
  116. if (trim(unique_obs(iuobs)) == 'ICEC') then
  117. do iens = 1, nrens
  118. write(cmem,'(i3.3)') iens
  119. tlevel = 1
  120. call get_mod_fld_new(trim('forecast'//cmem), readfld, iens,&
  121. 'icec', 0, tlevel, nx, ny)
  122. if (tlevel == -1) then
  123. if (master) then
  124. print *, 'ERROR: get_mod_fld_new(): failed for "icec"'
  125. end if
  126. stop
  127. end if
  128. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  129. readfld, depths, nx, ny, nz, 0)
  130. end do
  131. ! [FM, May 2013: for LIM3 sea ice model]
  132. elseif (trim(unique_obs(iuobs)) == 'AT_I') then
  133. do iens = 1, nrens
  134. write(cmem,'(i3.3)') iens
  135. tlevel = 1
  136. call io_mod_fld(ai1,iens,enslist, &
  137. 'a_i_htc1', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  138. call io_mod_fld(ai2,iens,enslist, &
  139. 'a_i_htc2', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  140. call io_mod_fld(ai3,iens,enslist, &
  141. 'a_i_htc3', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  142. call io_mod_fld(ai4,iens,enslist, &
  143. 'a_i_htc4', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  144. call io_mod_fld(ai5,iens,enslist, &
  145. 'a_i_htc5', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  146. if (tlevel == -1) then
  147. if (master) then
  148. print *, 'ERROR: io_mod_fld_new(): failed for "at_i"'
  149. end if
  150. stop
  151. end if
  152. readfld=ai1+ai2+ai3+ai4+ai5
  153. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  154. readfld, depths, nx, ny, nz, 0)
  155. end do
  156. elseif (trim(unique_obs(iuobs)) == 'SST') then
  157. do iens = 1, nrens
  158. write(cmem,'(i3.3)') iens
  159. tlevel = 1
  160. call get_mod_fld_new(trim('forecast'//cmem), readfld, iens,&
  161. 'temp', 1, tlevel, nx, ny)
  162. if (tlevel == -1) then
  163. if (master) then
  164. print *, 'ERROR: get_mod_fld_new(): failed for "SST"'
  165. end if
  166. stop
  167. end if
  168. if (prm_prmestexists('sstb')) then
  169. tlevel = 1
  170. call get_mod_fld_new(trim('forecast'//cmem), field2, iens,&
  171. 'sstb', 0, tlevel, nx, ny)
  172. if (tlevel == -1) then
  173. if (master) then
  174. print *, 'ERROR: get_mod_fld_new(): failed for "sstb"'
  175. end if
  176. stop
  177. end if
  178. readfld = readfld - field2
  179. end if
  180. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  181. readfld, depths, nx, ny, nz, 0)
  182. end do
  183. elseif (trim(unique_obs(iuobs)) == 'SLA' .or.&
  184. trim(unique_obs(iuobs)) == 'TSLA') then
  185. if (trim(unique_obs(iuobs)) == 'TSLA') then
  186. call read_mean_ssh(field2, nx, ny)
  187. end if
  188. inquire(iolength=reclSLA) fldr4
  189. ! FANF loop over each day of the week
  190. do t = 0, Wd
  191. if (trim(unique_obs(iuobs)) == 'TSLA') then
  192. write(day,'(i2.2)') t
  193. fname = trim('model_TSSH_'//day//'.uf')
  194. else
  195. fname = 'model_SLA.uf'
  196. endif
  197. if (master) then
  198. print *, 'TSLA, day', t, ': nobs = ',&
  199. count(obs(uobs_begin(iuobs) : uobs_end(iuobs)) % date == t)
  200. end if
  201. do iens = 1, nrens
  202. open(10, file = trim(fname), access = 'direct',&
  203. status = 'old', recl = reclSLA, action = 'read')
  204. read(10, rec = iens, iostat = ios) fldr4
  205. if (ios /= 0) then
  206. if (master) print *, 'Error reading ', trim(fname), ', member #', iens
  207. call stop_mpi()
  208. end if
  209. close(10)
  210. readfld = fldr4
  211. if (prm_prmestexists('msshb')) then
  212. write(cmem,'(i3.3)') iens
  213. tlevel = 1
  214. call get_mod_fld_new(trim('forecast'//cmem), field3, iens,&
  215. 'msshb', 0, tlevel, nx, ny)
  216. if (tlevel == -1) then
  217. if (master) then
  218. print *, 'ERROR: get_mod_fld_new(): failed for "msshb"'
  219. end if
  220. stop
  221. end if
  222. readfld = readfld - field3 ! mean SSH bias for this member
  223. end if
  224. if (trim(unique_obs(iuobs)) == 'TSLA') then
  225. readfld = readfld - field2 ! mean SSH
  226. end if
  227. call Generate_element_Si(S(:, iens), unique_obs(iuobs),&
  228. readfld, depths, nx, ny, nz, t)
  229. end do
  230. if (master) then
  231. print *, 'forming S, day', t
  232. print *, ' # of non-zero ens observations = ', count(S /= 0.0)
  233. print *, ' # of zero ens observations = ', count(S == 0.0)
  234. print *, ' # of non-zero observations for member 1 = ', count(S(:, 1) /= 0.0)
  235. end if
  236. end do
  237. elseif (trim(unique_obs(iuobs)) == 'SAL' .or.&
  238. trim(unique_obs(iuobs)) == 'TEM' .or.&
  239. trim(unique_obs(iuobs)) == 'GSAL' .or.&
  240. trim(unique_obs(iuobs)) == 'GTEM') then
  241. if (master) then
  242. print *, ' Interpolating ensemble vectors to the locations of "',&
  243. trim(unique_obs(iuobs)), '" observations'
  244. end if
  245. !
  246. ! for each ensemble member process all profiles "in parallel",
  247. ! reading the fields layer by layer
  248. !
  249. do iens = 1, nrens
  250. call get_S(S(:, iens), trim(unique_obs(iuobs)), nobs, obs, iens)
  251. end do
  252. if (master) then
  253. print *, ' Interpolation completed'
  254. end if
  255. elseif ((trim(unique_obs(iuobs)) == 'U_ICE') .or. trim(unique_obs(iuobs)) == 'V_ICE') then
  256. do iens = 1, nrens
  257. ! [FM] Read the file
  258. !inquire(iolength=reclDRIFT) modzon, modmer
  259. !open(10, file = 'model_ICEDRIFT.uf', access = 'direct',&
  260. ! status = 'old', recl = reclDRIFT, action = 'read')
  261. !read(10, rec = iens, iostat = ios) modzon, modmer
  262. !close(10)
  263. !if (ios /= 0) then
  264. ! if (master) then
  265. ! print *,'ERROR: could not read ensemble ice drift for member ', iens
  266. ! end if
  267. ! call stop_mpi()
  268. !end if
  269. call io_mod_fld(uice,iens,enslist, &
  270. 'u_ice', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  271. call io_mod_fld(vice,iens,enslist, &
  272. 'v_ice', 1, 0, 1, nx,ny, 'get',FLOAT(obs(1)%date))
  273. do m = 1, nobs
  274. i = obs(m) % i_orig_grid
  275. j = obs(m) % j_orig_grid
  276. if (trim(obs(m) % id) == 'U_ICE') then
  277. S(m, iens) = uice(i, j)
  278. elseif (trim(obs(m) % id) == 'V_ICE') then
  279. S(m, iens) = vice(i, j)
  280. end if
  281. end do
  282. end do
  283. elseif ((trim(unique_obs(iuobs)) == 'U2D_I') .OR. trim(unique_obs(iuobs)) == 'V2D_I' ) THEN
  284. ! ADDED BY FM FRANCOIS MASSONNET. u_ice_2d or v_ice_2d is the sea ice u or v-velocity
  285. ! obtained as follows:
  286. ! 1) Rotate OSISAF Low resolution 2-day sea ice drift in a {east,north}
  287. ! reference frame
  288. ! 2) Interpolate to the ORCA grid
  289. ! 3) Rotate to align with the ORCA grid
  290. ! 4) Multiply by 1000 and divide by 2*86400 to convert the 2-day
  291. ! displacement from km to m/s
  292. DO iens=1,nrens
  293. CALL read_icemod(uice,iens,enslist,'iicevelu',nx,ny)
  294. CALL read_icemod(vice,iens,enslist,'iicevelv',nx,ny)
  295. DO m = 1, nobs
  296. i = obs(m) % i_orig_grid
  297. j = obs(m) % j_orig_grid
  298. IF (trim(obs(m) % id) == 'U2D_I') THEN
  299. S(m,iens) = uice(i,j)
  300. ELSEIF (trim(obs(m) % id) == 'V2D_I') THEN
  301. S(m,iens) = vice(i,j)
  302. END IF
  303. END DO ! nobs
  304. END DO ! iens
  305. elseif ((index(unique_obs(iuobs),'DX') > 0 ) .or. (index(unique_obs(iuobs),'DY') > 0)) then
  306. ! OSISAF Ice drift observations (d-2-offset -> d-offset)
  307. !print *, 'Ice drift observation type: ', unique_obs(iuobs)
  308. offset = unique_obs(iuobs)(3:3)
  309. ! Use offset (1,2,3,4 or 5) to open correct model drift file
  310. inquire(iolength=reclDRIFT) dX, dY
  311. open(10, file = 'model_ICEDRIFT_OSISAF'//offset//'.uf', access = 'direct',&
  312. status = 'old', recl = reclDRIFT, action = 'read')
  313. do iens = 1, nrens
  314. read(10, rec = iens, iostat = ios) dX, dY
  315. if (ios /= 0) then
  316. if (master) then
  317. print *,'ERROR: could not read ensemble ice drift for member ', iens
  318. end if
  319. call stop_mpi()
  320. end if
  321. do m = 1, nobs
  322. i = obs(m) % i_orig_grid
  323. j = obs(m) % j_orig_grid
  324. if (index(obs(m)%id,'DX') > 0) then
  325. S(m, iens) = dX(i, j)
  326. elseif (index(obs(m)%id,'DY') > 0) then
  327. S(m, iens) = dY(i, j)
  328. end if
  329. end do
  330. end do
  331. close(10)
  332. else
  333. if (master) then
  334. print *,'ERROR: unknown obs type ' // trim(unique_obs(iuobs))
  335. end if
  336. call stop_mpi()
  337. end if
  338. end do ! iuobs
  339. ! some generic QC - relax fitting if the model and obs are too far apart
  340. !
  341. call obs_QC(nrens, S)
  342. ! add calculated HA to to observations-<type>.nc files for each data type
  343. !
  344. do iuobs = 1, nuobs
  345. if (master) then
  346. nthisobs = 0
  347. do m = 1, nobs
  348. if (trim(unique_obs(iuobs)) == trim(obs(m) % id)) then
  349. nthisobs = nthisobs + 1
  350. thisobs(nthisobs) = m
  351. end if
  352. end do
  353. ! add forecast values to the observation-<TYPE>.nc files
  354. !
  355. call add_forecast(unique_obs(iuobs), S(thisobs(1 : nthisobs), :), obs(thisobs(1 : nthisobs)))
  356. ! append the superobed values (and modified observation error
  357. ! variances) to the file with pre-processed observations (SAL.nc,
  358. ! TEM.nc, GSAL.nc or GTEM.nc)
  359. !
  360. if (trim(unique_obs(iuobs)) == 'SAL' .or.&
  361. trim(unique_obs(iuobs)) == 'TEM' .or.&
  362. trim(unique_obs(iuobs)) == 'GSAL' .or.&
  363. trim(unique_obs(iuobs)) == 'GTEM') then
  364. call insitu_writeforecast(unique_obs(iuobs), nobs, nrens, S, obs)
  365. end if
  366. end if
  367. end do
  368. if (master) then
  369. print *, 'm_prep_4_EnKF: end calculating S = HA'
  370. end if
  371. x = sum(S, DIM = 2) / real(nrens) ! [ FM ] The mean forecast interpolated in the obs.space
  372. if (master) print*,'m_prep_4_EnKF: end calculating Hx'
  373. if (master) then
  374. print *, 'Hx range = ', minval(x), '-', maxval(x)
  375. print *, 'mean(Hx) = ', sum(x) / real(nobs)
  376. end if
  377. if (master) then
  378. print *, 'observation range = ', minval(obs % d), '-', maxval(obs % d)
  379. print *, 'mean(observation) = ', sum(obs % d) / nobs
  380. end if
  381. ! Compute HA = HE - mean(HE)
  382. !
  383. if (master) print*,'prep_4_EnKF(): calculating S = S - x'
  384. do j = 1, nrens
  385. S(:, j) = S(:, j) - x ! [ FM ] This is really where we switch from actual model data to anomalies
  386. enddo
  387. ! Compute innovation
  388. !
  389. d = obs % d - x ! [ FM ] This is exactly was is also done in add_forecast. This is the mean innovation.
  390. if (master) then
  391. print *, ' innovation range = ', minval(d), '-', maxval(d)
  392. if (minval(d) < -1000.0d0) then
  393. print *, 'm_prep_4_EnKF: error: innovation too small detected'
  394. call stop_mpi()
  395. end if
  396. if (maxval(d) > 1000.0d0) then
  397. print *, 'm_prep_4_EnKF: error: innovation too big detected'
  398. call stop_mpi()
  399. end if
  400. end if
  401. end subroutine prep_4_EnKF
  402. subroutine read_mean_ssh(mean_ssh, nx, ny)
  403. #if defined (QMPI)
  404. use qmpi
  405. #else
  406. use qmpi_fake
  407. #endif
  408. use m_parameters
  409. integer, intent(in) :: nx, ny
  410. real, intent(out):: mean_ssh(nx, ny)
  411. logical :: exists
  412. inquire(file = trim(MEANSSHFNAME), exist = exists)
  413. if (.not. exists) then
  414. if (master) then
  415. print *,'ERROR: read_mean_ssh(): file "', trim(MEANSSHFNAME), '" not found'
  416. end if
  417. stop
  418. end if
  419. open (10, file = trim(MEANSSHFNAME), status = 'unknown',form = 'unformatted', action = 'read')
  420. read (10) mean_ssh
  421. close (10)
  422. end subroutine read_mean_ssh
  423. ! This subroutine adds forecast observations (i.e Hx) to the NetCDF
  424. ! observation files for each data type.
  425. !
  426. subroutine add_forecast(obstag, S, obs)
  427. use mod_measurement
  428. use nfw_mod
  429. implicit none
  430. character(OBSTYPESTRLEN), intent(in) :: obstag
  431. real, dimension(:, :), intent(in) :: S
  432. type(measurement), dimension(:) :: obs
  433. character(STRLEN) :: fname
  434. logical :: exists
  435. integer :: ncid
  436. integer :: dids(2), dimlen
  437. logical :: addsobs
  438. integer :: for_id, inn_id, forvar_id, slon_id, slat_id,&
  439. sdepth_id, sipiv_id, sjpiv_id, sd_id, svar_id,&
  440. newvar_id
  441. real, allocatable, dimension(:) :: x, Svar, innovation
  442. integer :: m, p, o
  443. write(fname, '(a, a, a)') 'observations-', trim(obstag), '.nc'
  444. inquire(file = trim(fname), exist = exists)
  445. if (.not. exists) then
  446. print *, 'file "', trim(fname), 'not found, skip adding forecast'
  447. return
  448. else
  449. print *, 'dumping forecast to "', trim(fname), '"'
  450. end if
  451. p = size(S, DIM = 1)
  452. m = size(S, DIM = 2)
  453. allocate(x(p), Svar(p), innovation(p))
  454. x = sum(S, DIM = 2) / real(m); ! [ FM the mean of S=HA ]
  455. Svar = 0.0
  456. do o = 1, p
  457. Svar(o) = sum((S(o, :) - x(o))** 2) ! [ FM thus each row of Svar is the variance (see below) of the forecast]
  458. end do
  459. Svar = Svar / real(m - 1)
  460. innovation = obs % d - x ! [ FM ] the innovation for the mean forecast (or mean of the innovation forecasts)
  461. addsobs = .false.
  462. call nfw_open(fname, nf_write, ncid)
  463. call nfw_inq_dimid(fname, ncid, 'nobs', dids(1))
  464. call nfw_inq_dimlen(fname, ncid, dids(1), dimlen)
  465. call nfw_redef(fname, ncid)
  466. if (dimlen == p) then
  467. dids(2) = dids(1)
  468. elseif (.not. nfw_dim_exists(ncid, 'nsobs')) then
  469. addsobs = .true.
  470. call nfw_def_dim(fname, ncid, 'nsobs', p, dids(2))
  471. call nfw_def_var(fname, ncid, 'slon', nf_float, 1, dids(2), slon_id)
  472. call nfw_def_var(fname, ncid, 'slat', nf_float, 1, dids(2), slat_id)
  473. call nfw_def_var(fname, ncid, 'sdepth', nf_float, 1, dids(2), sdepth_id)
  474. call nfw_def_var(fname, ncid, 'sipiv', nf_int, 1, dids(2), sipiv_id)
  475. call nfw_def_var(fname, ncid, 'sjpiv', nf_int, 1, dids(2), sjpiv_id)
  476. call nfw_def_var(fname, ncid, 'sd', nf_float, 1, dids(2), sd_id)
  477. call nfw_def_var(fname, ncid, 'svar', nf_float, 1, dids(2), svar_id)
  478. end if
  479. if (.not. nfw_var_exists(ncid, 'innovation')) then
  480. call nfw_def_var(fname, ncid, 'innovation', nf_double, 1, dids(2), inn_id)
  481. else
  482. call nfw_inq_varid(fname, ncid, 'innovation', inn_id)
  483. end if
  484. if (.not. nfw_var_exists(ncid, 'forecast')) then
  485. call nfw_def_var(fname, ncid, 'forecast', nf_double, 1, dids(2), for_id)
  486. else
  487. call nfw_inq_varid(fname, ncid, 'forecast', for_id)
  488. end if
  489. if (.not. nfw_var_exists(ncid, 'forecast_variance')) then
  490. call nfw_def_var(fname, ncid, 'forecast_variance', nf_double, 1, dids(2), forvar_id)
  491. else
  492. call nfw_inq_varid(fname, ncid, 'forecast_variance', forvar_id)
  493. end if
  494. if (.not. addsobs) then
  495. if (dimlen == p) then
  496. if (.not. nfw_var_exists(ncid, 'new_var')) then
  497. call nfw_def_var(fname, ncid, 'new_var', nf_double, 1, dids(2), newvar_id)
  498. else
  499. call nfw_inq_varid(fname, ncid, 'new_var', newvar_id)
  500. end if
  501. else
  502. if (.not. nfw_var_exists(ncid, 'new_svar')) then
  503. call nfw_inq_dimid(fname, ncid, 'nsobs', dids(2))
  504. call nfw_def_var(fname, ncid, 'new_svar', nf_double, 1, dids(2), newvar_id)
  505. else
  506. call nfw_inq_varid(fname, ncid, 'new_svar', newvar_id)
  507. end if
  508. end if
  509. end if
  510. call nfw_enddef(fname, ncid)
  511. call nfw_put_var_double(fname, ncid, forvar_id, Svar)
  512. call nfw_put_var_double(fname, ncid, for_id, x)
  513. call nfw_put_var_double(fname, ncid, inn_id, innovation)
  514. if (addsobs) then
  515. call nfw_put_var_double(fname, ncid, slon_id, obs % lon)
  516. call nfw_put_var_double(fname, ncid, slat_id, obs % lat)
  517. call nfw_put_var_double(fname, ncid, sdepth_id, obs % depth)
  518. call nfw_put_var_int(fname, ncid, sipiv_id, obs % ipiv)
  519. call nfw_put_var_int(fname, ncid, sjpiv_id, obs % jpiv)
  520. call nfw_put_var_double(fname, ncid, sd_id, obs % d)
  521. call nfw_put_var_double(fname, ncid, svar_id, obs % var)
  522. else
  523. call nfw_put_var_double(fname, ncid, newvar_id, obs % var)
  524. end if
  525. call nfw_close(fname, ncid)
  526. deallocate(x)
  527. deallocate(Svar)
  528. deallocate(innovation)
  529. end subroutine add_forecast
  530. end module m_prep_4_EnKF