m_insitu.F90 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788
  1. ! File: m_insitu.F90
  2. !
  3. ! Created: 6 Feb 2008
  4. !
  5. ! Last modified: 13 Feb 2008
  6. !
  7. ! Author: Pavel Sakov
  8. ! NERSC
  9. !
  10. ! Purpose: The code to deal with insitu observations.
  11. !
  12. ! Description: This module contains the following subroutines:
  13. ! - insitu_setprofiles
  14. ! breaks the measurements into profiles and returns
  15. ! arrays of start and end indices for each profile
  16. ! - insitu_writeprofiles
  17. ! writes profiles to a netCDF file
  18. ! - insitu_prepareobs
  19. ! sorts out the measurements within profiles so they
  20. ! go in surface to bottom order and thins the measurements
  21. ! by keeping max 1 measurements per layer of the first
  22. ! ensemble member
  23. ! It also contains the following data:
  24. ! nprof
  25. ! - the number of profiles
  26. ! pstart(nprof)
  27. ! - start indices for each profile in the array "obs" of
  28. ! type(measurement) stored in module m_obs
  29. ! pend(nprof)
  30. ! - end indices for each profile
  31. !
  32. ! Modifications:
  33. ! 30/7/2010 PS: added profile pivot points to profile output
  34. ! files (SAL.nc etc.)
  35. ! 29/7/2010 PS: some rather minor changes, including interface
  36. ! of insitu_writeforecast()
  37. ! 13/02/2008 PS: added insitu_writeprofiles()
  38. ! 26/02/2008 PS: put "nprof", "pstart" and "pend" as public data
  39. ! in this module
  40. ! 20/04/2008 PS: added insitu_QC() and insitu_writeforecast()
  41. ! 29/07/2010 PS: removed insitu_QC(). There is a generic obs QC
  42. ! procedure in m_obs.F90 now.
  43. module m_insitu
  44. use mod_measurement
  45. !use m_parse_blkdat
  46. use m_get_mod_xyz
  47. use m_get_mod_fld
  48. use m_io_mod_fld
  49. #if defined (QMPI)
  50. use qmpi
  51. #else
  52. use qmpi_fake
  53. #endif
  54. implicit none
  55. !
  56. ! public stuff
  57. !
  58. integer, allocatable, dimension(:), public :: pstart
  59. integer, allocatable, dimension(:), public :: pend
  60. integer, public :: nprof
  61. public insitu_setprofiles
  62. public insitu_prepareobs
  63. public insitu_writeprofiles
  64. !
  65. ! private stuff
  66. !
  67. real, parameter, private :: ONEMETER = 9806.0
  68. integer, parameter, private :: STRLEN = 512
  69. ! The portion of the layer thickness at which the variability in
  70. ! vertical data will be used for estimating the vertical representativeness
  71. ! error.
  72. !
  73. real, parameter, private :: VARCOEFF1 = 0.15
  74. ! A factor by which a calculated vertical representativeness error variance
  75. ! will be reduced if the data is in different layers
  76. !
  77. real, parameter, private :: VARCOEFF2 = 2.0
  78. ! Write information about this profile. Set to < 1 to switch off.
  79. !
  80. integer, parameter, private :: PDEBUGINFO = 0
  81. ! Integers used to tag the fields (to avoid parsing the string tags)
  82. !
  83. integer, parameter, private :: NONE = 0
  84. integer, parameter, private :: TEMPERATURE = 1
  85. integer, parameter, private :: SALINITY = 2
  86. real, parameter, private :: TEM_MIN = -2.0
  87. real, parameter, private :: TEM_MAX = 35.0
  88. real, parameter, private :: SAL_MIN = 5.0
  89. real, parameter, private :: SAL_MAX = 41.0
  90. ! Maximum allowed deviation between the observation and ensemble mean in
  91. ! terms of combined standard deviation.
  92. !
  93. real, parameter, private :: SAL_MAXRATIO = 10.0
  94. real, parameter, private :: TEM_MAXRATIO = 5.0
  95. ! If an observation is not considered an outlier, the observation error
  96. ! variance is modified so that the distance between the observation and the
  97. ! endemble mean is within DIST_MAX * sqrt(sigma_obs^2 + sigma_ens^2).
  98. ! Bigger values of DIST_MAX result in a more violent assimilation.
  99. !
  100. real, parameter, private :: DIST_MAX = 2.0
  101. contains
  102. ! Work out the number of profiles, each identified by "obs % i_orig_grid"
  103. ! and return start id of the first and the last obs in the profile in
  104. ! arrays "pstart" and "pend". "pstart" and "pend" are publicly available
  105. ! arrays stored in this module.
  106. !
  107. subroutine insitu_setprofiles(obstag, nobs, obs)
  108. character(*), intent(in) :: obstag
  109. integer, intent(in) :: nobs
  110. type(measurement), dimension(:), intent(inout) :: obs
  111. integer, allocatable, dimension(:) :: tmp, tmp1
  112. integer :: o, o1, o2, p, nobsp
  113. type(measurement), allocatable, dimension(:) :: tmpobs
  114. if (nobs == 0) then
  115. return
  116. end if
  117. if (allocated(pstart)) then
  118. deallocate(pstart)
  119. deallocate(pend)
  120. end if
  121. ! find the very first obs of the right kind
  122. !
  123. o1 = 1
  124. do while (trim(obs(o1) % id) /= trim(obstag) .and. o1 <= nobs)
  125. o1 = o1 + 1
  126. end do
  127. if (o1 > nobs) then
  128. return
  129. end if
  130. ! find the very last obs of the right kind
  131. !
  132. o2 = nobs
  133. do while (trim(obs(o2) % id) /= trim(obstag) .and. o2 >= 0)
  134. o2 = o2 - 1
  135. end do
  136. nprof = 1
  137. do o = 2, o2
  138. if (obs(o) % ipiv /= obs(o - 1) % ipiv .or.&
  139. obs(o) % jpiv /= obs(o - 1) % jpiv .or.&
  140. obs(o) % date /= obs(o - 1) % date) then
  141. nprof = nprof + 1
  142. end if
  143. end do
  144. allocate(pstart(nprof))
  145. allocate(pend(nprof))
  146. ! identify profiles
  147. !
  148. ! PS: This is a tricky cycle but it seems it is doing the job. Do not
  149. ! meddle with it.
  150. !
  151. pend = 0
  152. nprof = 1
  153. pstart(1) = o1
  154. do o = o1, o2
  155. ! find obs from the same profile
  156. !
  157. if (trim(obs(o) % id) == trim(obstag) .and.&
  158. ((obs(o) % i_orig_grid > 0 .and.&
  159. obs(o) % i_orig_grid == obs(pstart(nprof)) % i_orig_grid) .or.&
  160. (obs(o) % i_orig_grid <= 0 .and.&
  161. obs(o) % ipiv == obs(pstart(nprof)) % ipiv .and.&
  162. obs(o) % jpiv == obs(pstart(nprof)) % jpiv .and.&
  163. obs(o) % date == obs(pstart(nprof)) % date))) then
  164. pend(nprof) = o
  165. cycle
  166. end if
  167. if (trim(obs(o) % id) /= trim(obstag)) then
  168. print *, 'ERROR: insitu_setprofiles(): obs id does not match processed obs tag'
  169. stop
  170. end if
  171. ! if there were no obs of the right type in this profile yet,
  172. ! then pend(nprof) has not been set yet and therefore the condition
  173. ! below will yield "false"
  174. !
  175. if (pend(nprof) >= pstart(nprof)) then
  176. nprof = nprof + 1
  177. end if
  178. if (PDEBUGINFO > 0) then
  179. print *, ' DEBUG: new profile #', nprof, ', o =', o, ', id =', obs(o) % i_orig_grid
  180. end if
  181. pstart(nprof) = o
  182. pend(nprof) = o
  183. end do
  184. if (pend(nprof) < pstart(nprof)) then
  185. nprof = nprof - 1
  186. end if
  187. ! truncate "pstat" and "pend" to length "nprof"
  188. !
  189. allocate(tmp(nprof))
  190. tmp = pstart(1 : nprof)
  191. deallocate(pstart)
  192. allocate(pstart(nprof))
  193. pstart = tmp
  194. tmp = pend(1 : nprof)
  195. deallocate(pend)
  196. allocate(pend(nprof))
  197. pend = tmp
  198. deallocate(tmp)
  199. ! for glider data - sort observations in each profile by increasing depth
  200. !
  201. if (trim(obstag) == 'GSAL'.or. trim(obstag) == 'GTEM') then
  202. allocate(tmp(nobs))
  203. allocate(tmp1(nobs))
  204. allocate(tmpobs(nobs))
  205. do p = 1, nprof
  206. nobsp = pend(p) - pstart(p) + 1
  207. do o = 1, nobsp
  208. tmp(o) = o
  209. end do
  210. !
  211. ! (using procedure from pre_local_analysis())
  212. !
  213. call order(dble(nobsp), obs(pstart(p) : pend(p)) % depth,&
  214. dble(nobsp), tmp, tmp1)
  215. tmpobs(1 : nobsp) = obs(pstart(p) : pend(p))
  216. do o = 1, nobsp
  217. obs(pstart(p) + o - 1) = tmpobs(tmp1(o))
  218. end do
  219. end do
  220. deallocate(tmp, tmp1, tmpobs)
  221. end if
  222. end subroutine insitu_setprofiles
  223. ! 1. Sort out the obs within profiles so that they are stored in order of
  224. ! increasing depth.
  225. ! 2. Thin observations by keeping a single obs within a layer using the
  226. ! layers from the first ensemble member
  227. !
  228. subroutine insitu_prepareobs(obstag, nobs, obs)
  229. character(*), intent(in) :: obstag
  230. integer, intent(inout) :: nobs
  231. type(measurement), dimension(:), intent(inout) :: obs
  232. ! profiles
  233. !
  234. integer, allocatable, dimension(:) :: pnow
  235. integer :: nobs_max
  236. integer :: p, o
  237. type(measurement), allocatable, dimension(:) :: profile
  238. integer, allocatable, dimension(:) :: ipiv, jpiv
  239. real, allocatable, dimension(:) :: a1, a2, a3, a4
  240. real, allocatable, dimension(:) :: z1, z2
  241. integer :: nrev
  242. integer :: ndel
  243. integer :: oo
  244. real :: rdummy
  245. integer :: k, nk, ni, nj
  246. character(80) :: fname
  247. integer :: tlevel
  248. real, allocatable, dimension(:, :) :: dz2d
  249. real, dimension(2, 2) :: dz_cell
  250. real :: dz, zcentre
  251. integer :: best
  252. logical :: isrogue
  253. ! As we thin the measurements within each layer, it still may be a good
  254. ! idea to update the obs error variance if the variability within the layer
  255. ! is big enough. `dmin' and `dmax' give the min and max measured values
  256. ! within the layer.
  257. !
  258. real :: dmin, dmax
  259. real :: var1, var2
  260. integer :: nobsnew, nobs_thistype, nobs_othertype
  261. if (master) then
  262. print '(a, a, a)', ' insitu_prepareobs(', trim(obstag), '):'
  263. print '(a, i6)', ' total # of obs = ', nobs
  264. end if
  265. if (nobs == 0) then
  266. return
  267. end if
  268. call insitu_setprofiles(trim(obstag), nobs, obs)
  269. if (master) then
  270. print '(a, a, a, i6)', ' # of obs of type "', trim(obstag), '" = ',&
  271. sum(pend(1 : nprof) - pstart(1 : nprof)) + nprof
  272. print '(a, i4)', ' # of profiles = ', nprof
  273. end if
  274. ! find the maximal # of obs in a single profile
  275. !
  276. nobs_max = 0
  277. do p = 1, nprof
  278. nobs_max = max(nobs_max, pend(p) - pstart(p) + 1)
  279. end do
  280. if (master) then
  281. print '(a, i4)', ' max # of obs in a profile before thinning = ', nobs_max
  282. end if
  283. ! reverse the obs in profiles that go from bottom to surface
  284. !
  285. allocate(profile(nobs_max))
  286. nrev = 0
  287. do p = 1, nprof
  288. if (obs(pstart(p)) % depth > obs(pend(p)) % depth) then
  289. profile(1 : pend(p) - pstart(p) + 1) = obs(pstart(p) : pend(p))
  290. do o = 0, pend(p) - pstart(p)
  291. obs(pstart(p) + o) = profile(pend(p) - o)
  292. end do
  293. nrev = nrev + 1
  294. end if
  295. end do
  296. deallocate(profile)
  297. if (nrev > 0 .and. master) then
  298. print *, ' ', nrev, ' profile(s) reversed'
  299. end if
  300. ! check for rogue obs
  301. !
  302. ndel = 0
  303. do p = 1, nprof
  304. isrogue = .false.
  305. do o = pstart(p) + 1, pend(p)
  306. ! shift the remaining obs in this profile one obs down
  307. !
  308. if (obs(o) % depth <= obs(o - 1) % depth) then
  309. isrogue = .true.
  310. do oo = o + 1, pend(p)
  311. obs(oo - 1) = obs(oo)
  312. end do
  313. ndel = ndel + 1
  314. pend(p) = pend(p) - 1
  315. end if
  316. end do
  317. if (isrogue .and. master) then
  318. print *, ' a rogue obs detected in profile # ', p
  319. end if
  320. end do
  321. if (ndel > 0 .and. master) then
  322. print *, ' ', ndel, 'rogue obs deleted'
  323. end if
  324. !
  325. ! Now to the thinning of the profiles.
  326. !
  327. allocate(ipiv(nprof))
  328. allocate(jpiv(nprof))
  329. allocate(a1(nprof))
  330. allocate(a2(nprof))
  331. allocate(a3(nprof))
  332. allocate(a4(nprof))
  333. ipiv = obs(pstart(1 : nprof)) % ipiv
  334. jpiv = obs(pstart(1 : nprof)) % jpiv
  335. a1 = obs(pstart(1 : nprof)) % a1
  336. a2 = obs(pstart(1 : nprof)) % a2
  337. a3 = obs(pstart(1 : nprof)) % a3
  338. a4 = obs(pstart(1 : nprof)) % a4
  339. ! get the grid dimensions
  340. !
  341. !call parse_blkdat('kdm ','integer', rdummy, nk)
  342. !call parse_blkdat('idm ','integer', rdummy, ni)
  343. !call parse_blkdat('jdm ','integer', rdummy, nj)
  344. call get_mod_xyz(ni, nj, nk) ! CKB,FM Changed from using m_parse_blkdat
  345. ! get the data file name
  346. !
  347. if (trim(obstag) /= 'SAL' .and. trim(obstag) /= 'TEM' .and.&
  348. trim(obstag) /= 'GSAL'.and. trim(obstag) /= 'GTEM') then
  349. print *, 'ERROR: get_S(): unknown observation tag "', trim(obstag), '"'
  350. stop
  351. end if
  352. fname = 'forecast001'
  353. allocate(z1(nprof))
  354. allocate(z2(nprof))
  355. allocate(pnow(nprof))
  356. allocate(dz2d(ni, nj))
  357. ! data thinning cycle
  358. !
  359. if (master) then
  360. print *, ' maximum one observation per layer will be retained after thinning'
  361. end if
  362. tlevel = 1
  363. z1 = 0.0
  364. pnow = pstart
  365. if (master .and. PDEBUGINFO > 0) then
  366. p = PDEBUGINFO
  367. print *, 'DEBUG dumping the info for profile #', p
  368. print *, 'DEBUG p =', p, ': lon =', obs(pstart(p)) % lon, ', lat =', obs(pstart(p)) % lat
  369. print *, 'DEBUG now dumping the layer depths:'
  370. end if
  371. ! mark all obs of this type as bad; unmask the best obs within a layer
  372. !
  373. do o = 1, nobs
  374. if (trim(obs(o) % id) == trim(obstag)) then
  375. obs(o) % status = .false.
  376. end if
  377. end do
  378. do k = 1, nk
  379. !call get_mod_fld_new(trim(fname), dz2d, 1, 'dp ', k, tlevel, ni, nj)
  380. ! [CKB,FM]
  381. call io_mod_fld(dz2d, 1, (/ 1 /), 'dp ', 2, &
  382. k, tlevel, ni, nj, 'get',FLOAT(obs(1)%date))
  383. do p = 1, nprof
  384. dz_cell(:, :) = dz2d(ipiv(p) : ipiv(p) + 1, jpiv(p) : jpiv(p) + 1)
  385. dz = dz_cell(1, 1) * a1(p) + dz_cell(2, 1) * a2(p)&
  386. + dz_cell(1, 2) * a3(p) + dz_cell(2, 2) * a4(p)
  387. dz = dz / ONEMETER
  388. z2(p) = z1(p) + dz
  389. zcentre = (z1(p) + z2(p)) / 2.0
  390. best = -1
  391. dmin = 1.0d+10
  392. dmax = -1.0d+10
  393. if (master .and. PDEBUGINFO > 0 .and. p == PDEBUGINFO) then
  394. print *, 'DEBUG p =', p, ', k =', k, ', z =', z1(p), '-', z2(p)
  395. end if
  396. do while (pnow(p) <= pend(p))
  397. o = pnow(p)
  398. ! check that the depth is within the layer
  399. !
  400. if (obs(o) % depth > z2(p)) then
  401. ! go to next profile; this obs will be dealt with when
  402. ! processing the next layer
  403. exit
  404. end if
  405. ! from this point on, the obs counter will be increased at the
  406. ! end of this loop
  407. ! store profile and layer number (overwrite the original profile
  408. ! id and vertical counter value)
  409. !
  410. obs(o) % i_orig_grid = p
  411. obs(o) % j_orig_grid = k
  412. obs(o) % h = z2(p) - z1(p)
  413. if (obs(o) % depth < z1(p)) then
  414. pnow(p) = pnow(p) + 1
  415. cycle ! next obs
  416. end if
  417. ! update `dmin' and `dmax'
  418. !
  419. dmin = min(dmin, obs(o) % d)
  420. dmax = max(dmax, obs(o) % d)
  421. if (best < 1) then
  422. best = o
  423. obs(best) % status = .true.
  424. else if (abs(obs(o) % depth - zcentre) < abs(obs(best) % depth - zcentre)) then
  425. obs(best) % status = .false. ! thrash the previous best obs
  426. best = o
  427. obs(best) % status = .true.
  428. end if
  429. pnow(p) = pnow(p) + 1
  430. end do ! o
  431. ! update the observation error variance if the difference between
  432. ! `dmin' and `dmax' is big enough
  433. !
  434. if (best < 1) then
  435. cycle
  436. end if
  437. if (.false.) then ! out for now; use the closest obs instead
  438. if (dmax - dmin > 0) then
  439. obs(best) % var = sqrt(obs(best) % var + ((dmax - dmin) / 2) ** 2)
  440. end if
  441. end if
  442. end do ! p
  443. z1 = z2
  444. end do ! k
  445. ! There are a number of ways the vertical variability can be
  446. ! used for updating the obs error variance.
  447. !
  448. ! Below, the following approach is used.
  449. !
  450. ! Calculate two estimates for vertical gradient using the closest data
  451. ! points (if available). Estimate the difference at (VARCOEFF1 * h)
  452. ! vertical distance from the current obs, where VARCOEFF1 is the portion
  453. ! of the layer thickness (typically around 0.1-0.3), and h is the layer
  454. ! thickness. Use the square of this difference as an estimate for the
  455. ! respresentation error variance. If the closest obs is in another layer
  456. ! -- decrease this estimate by a factor of VARCOEFF2 (typically around 2).
  457. ! Use the largest estimate between the two (when both are avalaible).
  458. !
  459. do p = 1, nprof
  460. do o = pstart(p), pend(p)
  461. k = obs(o) % j_orig_grid
  462. if (obs(o) % status) then
  463. var1 = -999.0
  464. var2 = -999.0
  465. if (o - 1 >= pstart(p)) then
  466. var1 = ((obs(o) % d - obs(o - 1) % d) /&
  467. (obs(o) % depth - obs(o - 1) % depth) * obs(o) % h * VARCOEFF1) ** 2
  468. if (obs(o - 1) % j_orig_grid /= k) then
  469. var1 = var1 / VARCOEFF2
  470. end if
  471. end if
  472. if (o + 1 <= pend(p)) then
  473. var2 = ((obs(o) % d - obs(o + 1) % d) /&
  474. (obs(o) % depth - obs(o + 1) % depth) * obs(o) % h * VARCOEFF1) ** 2
  475. if (obs(o + 1) % j_orig_grid /= k) then
  476. var2 = var2 / VARCOEFF2
  477. end if
  478. end if
  479. if (var1 < 0.0 .and. var2 < 0.0) then
  480. cycle
  481. end if
  482. obs(o) % var = obs(o) % var + max(var1, var2)
  483. end if
  484. end do
  485. end do
  486. if (master .and. PDEBUGINFO > 0) then
  487. p = PDEBUGINFO
  488. print *, 'DEBUG now dumping the obs info:'
  489. do o = pstart(p), pend(p)
  490. print *, 'DEBUG o =', o, ', status =', obs(o) % status, &
  491. ', d =', obs(o) % d, ', z =', obs(o) % depth,&
  492. ', k =', obs(o) % j_orig_grid, ', h =', obs(o) % h,&
  493. ', var =', obs(o) % var
  494. end do
  495. end if
  496. deallocate(dz2d)
  497. deallocate(pnow)
  498. deallocate(z2)
  499. deallocate(z1)
  500. deallocate(a4)
  501. deallocate(a3)
  502. deallocate(a2)
  503. deallocate(a1)
  504. deallocate(jpiv)
  505. deallocate(ipiv)
  506. ! now compact the obs array
  507. !
  508. nobsnew = 0
  509. nobs_thistype = 0
  510. nobs_othertype = 0
  511. do o = 1, nobs
  512. if (obs(o) % status) then
  513. nobsnew = nobsnew + 1
  514. obs(nobsnew) = obs(o)
  515. if (trim(obs(o) % id) == trim(obstag)) then
  516. nobs_thistype = nobs_thistype + 1
  517. else
  518. nobs_othertype = nobs_othertype + 1
  519. end if
  520. end if
  521. end do
  522. obs(nobsnew + 1 : nobs) % status = .false.
  523. nobs = nobsnew
  524. ! replace the original profiles by the thinned ones
  525. !
  526. call insitu_setprofiles(trim(obstag), nobs, obs)
  527. if (master) then
  528. print *, ' thinning completed:', nobs_thistype, ' "', trim(obstag), '" obs retained'
  529. if (nobs_othertype > 0) then
  530. print *, ' ', nobs_othertype, 'obs of other type(s) retained'
  531. end if
  532. end if
  533. end subroutine insitu_prepareobs
  534. ! Write profiles to a NetCDF file
  535. !
  536. subroutine insitu_writeprofiles(fname, obstag, nobs, obs)
  537. use nfw_mod
  538. character(*), intent(in) :: fname
  539. character(*), intent(in) :: obstag
  540. integer, intent(inout) :: nobs
  541. type(measurement), dimension(:), intent(inout) :: obs
  542. ! profiles
  543. !
  544. integer :: p
  545. integer :: npoints, npoints_max
  546. ! I/O
  547. !
  548. integer :: ncid
  549. integer :: nprof_id(1), nk_id(1), dids(2)
  550. integer :: lat_id, lon_id, ipiv_id, jpiv_id, npoints_id, depth_id, v_id, variance_id
  551. character(STRLEN) :: varname
  552. real(8), allocatable, dimension(:, :) :: v
  553. if (.not. allocated(pstart)) then
  554. call insitu_setprofiles(trim(obstag), nobs, obs)
  555. end if
  556. call nfw_create(fname, nf_write, ncid)
  557. call nfw_def_dim(fname, ncid, 'nprof', nprof, nprof_id(1))
  558. call nfw_def_var(fname, ncid, 'lat', nf_double, 1, nprof_id, lat_id)
  559. call nfw_def_var(fname, ncid, 'lon', nf_double, 1, nprof_id, lon_id)
  560. call nfw_def_var(fname, ncid, 'ipiv', nf_int, 1, nprof_id, ipiv_id)
  561. call nfw_def_var(fname, ncid, 'jpiv', nf_int, 1, nprof_id, jpiv_id)
  562. call nfw_def_var(fname, ncid, 'npoints', nf_int, 1, nprof_id, npoints_id)
  563. npoints_max = maxval(pend - pstart) + 1
  564. call nfw_def_dim(fname, ncid, 'nk', npoints_max, nk_id(1))
  565. dids(1) = nk_id(1)
  566. dids(2) = nprof_id(1)
  567. call nfw_def_var(fname, ncid, 'depth', nf_double, 2, dids, depth_id)
  568. if (trim(obstag) == 'SAL' .or. trim(obstag) == 'GSAL') then
  569. varname = 'salt'
  570. else if (trim(obstag) == 'TEM' .or. trim(obstag) == 'GTEM') then
  571. varname = 'temp'
  572. else
  573. varname = trim(obstag)
  574. end if
  575. call nfw_def_var(fname, ncid, trim(varname), nf_double, 2, dids, v_id)
  576. call nfw_def_var(fname, ncid, 'variance', nf_double, 2, dids, variance_id)
  577. call nfw_enddef(fname, ncid)
  578. call nfw_put_var_double(fname, ncid, lat_id, obs(pstart) % lat)
  579. call nfw_put_var_double(fname, ncid, lon_id, obs(pstart) % lon)
  580. call nfw_put_var_int(fname, ncid, ipiv_id, obs(pstart) % ipiv)
  581. call nfw_put_var_int(fname, ncid, jpiv_id, obs(pstart) % jpiv)
  582. call nfw_put_var_int(fname, ncid, npoints_id, pend - pstart + 1)
  583. ! depth
  584. !
  585. allocate(v(npoints_max, nprof))
  586. v = -999.0
  587. do p = 1, nprof
  588. npoints = pend(p) - pstart(p) + 1
  589. v(1 : npoints, p) = obs(pstart(p) : pend(p)) % depth
  590. end do
  591. call nfw_put_var_double(fname, ncid, depth_id, v)
  592. ! data
  593. !
  594. v = -999.0
  595. do p = 1, nprof
  596. npoints = pend(p) - pstart(p) + 1
  597. v(1 : npoints, p) = obs(pstart(p) : pend(p)) % d
  598. end do
  599. call nfw_put_var_double(fname, ncid, v_id, v)
  600. ! data error variance
  601. !
  602. v = -999.0
  603. do p = 1, nprof
  604. npoints = pend(p) - pstart(p) + 1
  605. v(1 : npoints, p) = obs(pstart(p) : pend(p)) % var
  606. end do
  607. call nfw_put_var_double(fname, ncid, variance_id, v)
  608. call nfw_close(fname, ncid)
  609. deallocate(v)
  610. deallocate(pstart)
  611. deallocate(pend)
  612. end subroutine insitu_writeprofiles
  613. ! This subroutine appends the interpolated ensemble mean and the ensemble
  614. ! error variance to the assimilated profile data SAL.nc or TEM.nc. It also
  615. ! overwrites the observation error variance with latest values.
  616. !
  617. subroutine insitu_writeforecast(obstag, nobs, nrens, S, obs)
  618. use nfw_mod
  619. implicit none
  620. character(*), intent(in) :: obstag
  621. integer, intent(in) :: nobs
  622. integer, intent(in) :: nrens
  623. real, dimension(nobs, nrens), intent(in) :: S
  624. type(measurement), dimension(nobs), intent(inout) :: obs
  625. character(STRLEN) :: fname
  626. real, dimension(nobs) :: Smean, Svar
  627. integer :: i, p
  628. integer :: ncid
  629. integer :: dids(2)
  630. integer :: v_id, variance_id
  631. integer :: npoints_max, npoints
  632. real(8), allocatable, dimension(:, :) :: v
  633. ! need to set profiles for the given observation type
  634. !
  635. call insitu_setprofiles(obstag, nobs, obs)
  636. write(fname, '(a, ".nc")') trim(obstag)
  637. print *, 'Appending interpolated forecast for "', trim(obstag),&
  638. '" to "', trim(fname), '"'
  639. Smean = sum(S, DIM = 2) / nrens
  640. Svar = 0.0
  641. do i = 1, nobs
  642. Svar(i) = sum((S(i, :) - Smean(i)) ** 2)
  643. end do
  644. Svar = Svar / real(nrens - 1)
  645. call nfw_open(fname, nf_write, ncid)
  646. call nfw_inq_dimid(fname, ncid, 'nk', dids(1))
  647. call nfw_inq_dimid(fname, ncid, 'nprof', dids(2))
  648. call nfw_redef(fname, ncid)
  649. call nfw_def_var(fname, ncid, 'forecast', nf_double, 2, dids, v_id)
  650. call nfw_def_var(fname, ncid, 'forecast_variance', nf_double, 2, dids, variance_id)
  651. call nfw_enddef(fname, ncid)
  652. npoints_max = maxval(pend - pstart) + 1
  653. allocate(v(npoints_max, nprof))
  654. v = -999.0
  655. do p = 1, nprof
  656. npoints = pend(p) - pstart(p) + 1
  657. v(1 : npoints, p) = Smean(pstart(p) : pend(p))
  658. end do
  659. call nfw_put_var_double(fname, ncid, v_id, v)
  660. v = -999.0
  661. do p = 1, nprof
  662. npoints = pend(p) - pstart(p) + 1
  663. v(1 : npoints, p) = Svar(pstart(p) : pend(p))
  664. end do
  665. call nfw_put_var_double(fname, ncid, variance_id, v)
  666. ! update observation error variance
  667. !
  668. call nfw_redef(fname, ncid)
  669. call nfw_rename_var(fname, ncid, 'variance', 'variance_orig')
  670. call nfw_def_var(fname, ncid, 'variance', nf_double, 2, dids, variance_id)
  671. call nfw_enddef(fname, ncid)
  672. v = -999.0
  673. do p = 1, nprof
  674. npoints = pend(p) - pstart(p) + 1
  675. v(1 : npoints, p) = obs(pstart(p) : pend(p)) % var
  676. end do
  677. call nfw_put_var_double(fname, ncid, variance_id, v)
  678. call nfw_close(fname, ncid)
  679. deallocate(v)
  680. end subroutine insitu_writeforecast
  681. end module m_insitu