m_superobs.F90 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312
  1. ! File: m_superobs.F90
  2. !
  3. ! Created: 02 Sep 2008
  4. !
  5. ! Author: Pavel Sakov
  6. ! NERSC
  7. !
  8. ! Purpose: Superobing observations to model grid
  9. !
  10. ! Description: Conducts the following operations:
  11. ! - determine the number of observations with this tag
  12. ! - sort observations according to the pivot values
  13. ! - calculate superobs
  14. !
  15. ! Modifications: 14.10.2009 PS: added cycle over the data age, so that it only
  16. ! superobs the data of the same age. Do not set
  17. ! age if you assiilate data of different age in
  18. ! one go.
  19. ! 15.11.2009 PS: fixed a defect at l.102: it should be
  20. ! "thisob = obs_now(sorted(o))", not
  21. ! "thisob = obs_now(o)"
  22. ! 17.11.2009 PS: extended to handle the 3D case
  23. module m_superobs
  24. use mod_measurement
  25. use m_bilincoeff
  26. implicit none
  27. integer, parameter, private :: STRLEN = 512
  28. logical, parameter, private :: TEST = .false.
  29. contains
  30. subroutine superob(obstag, nobs, obs, ni, nj, modlon, modlat, nnewobs, newobs, is3d)
  31. character(*), intent(in) :: obstag
  32. integer, intent(in) :: nobs
  33. type(measurement), intent(inout), dimension(:) :: obs
  34. integer, intent(in) :: ni, nj
  35. real, dimension(:,:), intent(in) :: modlon, modlat
  36. integer, intent(inout) :: nnewobs
  37. type(measurement), intent(inout), dimension(:) :: newobs
  38. logical, intent(in), optional :: is3d
  39. integer :: age_min, age_max, nobs_total, nobs_now, age
  40. integer :: o, iprev, jprev, kprev, ii, ii_now
  41. logical, dimension(nobs) :: mask
  42. integer, dimension(nobs) :: sorted
  43. type(measurement), dimension(nobs) :: obs_now
  44. real(8), dimension(1) :: nobs_real
  45. type(measurement) :: thisob
  46. real :: n, nmax, valsum, valsqsum, varinvsum, lonsum, latsum, depthsum, valmax, valmin
  47. real :: a1sum, a2sum, a3sum, a4sum
  48. integer :: nlon_pos, nlon_neg
  49. real :: lonsum_abs
  50. integer, dimension(nobs) :: kpiv ! vertical index for 3D case
  51. integer, dimension(nobs) :: ids ! ids of obs contributing to this superob
  52. integer :: fid
  53. kpiv(1:nobs)=0
  54. ! find the range of the data age
  55. !
  56. age_min = minval(obs % date)
  57. age_max = maxval(obs % date)
  58. print *, 'min age =', age_min
  59. print *, 'max age =', age_max
  60. ! get the total number of observations to process
  61. !
  62. mask = .false.
  63. do o = 1, nobs
  64. if (trim(obs(o) % id) == trim(obstag)) then
  65. mask(o) = .true.
  66. end if
  67. obs(o) % orig_id = o
  68. end do
  69. nobs_total = count(mask)
  70. print *, 'total # of obs of all types =', nobs
  71. print *, 'total # of obs of type "', trim(obstag), '" =', nobs_total
  72. if (TEST) then
  73. open(101, file = 'superobs.txt', access = 'sequential', status = 'replace')
  74. end if
  75. ii = 0
  76. do age = age_min, age_max
  77. ! trim() prevents vectorising below
  78. mask = .false.
  79. do o = 1, nobs
  80. if (trim(obs(o) % id) == trim(obstag) .and. obs(o) % date == age .and. obs(o) % status) then
  81. mask(o) = .true.
  82. end if
  83. end do
  84. nobs_now = count(mask)
  85. print *, 'age =', age
  86. print *, ' nobs =', nobs_now
  87. if (nobs_now == 0) then
  88. cycle
  89. end if
  90. obs_now(1 : nobs_now) = pack(obs(1 : nobs), mask)
  91. nobs_real(1) = nobs_now
  92. if (.not. present(is3d) .or. .not. is3d) then
  93. call sortgriddedobs(nobs_real, obs_now % ipiv, obs_now % jpiv, sorted)
  94. else
  95. kpiv = z2k(obs_now % depth)
  96. call sortgriddedobs3d(nobs_real, obs_now % ipiv, obs_now % jpiv,&
  97. kpiv, sorted)
  98. end if
  99. iprev = 0
  100. jprev = 0
  101. kprev = 0
  102. nmax = 0
  103. ii_now = 0
  104. do o = 1, nobs_now + 1
  105. if (o <= nobs_now) then
  106. thisob = obs_now(sorted(o))
  107. else
  108. thisob % ipiv = -1 ! to force write of the previous measurement
  109. end if
  110. if (thisob % ipiv /= iprev .or. thisob % jpiv /= jprev .or. kpiv(sorted(o)) /= kprev) then
  111. if (ii_now > 0) then ! write the previous measurement
  112. newobs(ii) % d = valsum / n
  113. newobs(ii) % var = 1.0d0 / varinvsum
  114. newobs(ii) % id = obstag
  115. if (nlon_pos == 0 .or. nlon_neg == 0 .or. lonsum_abs / n < 90.0d0) then
  116. newobs(ii) % lon = lonsum / n
  117. else
  118. lonsum = lonsum + real(nlon_neg) * 360.0d0;
  119. newobs(ii) % lon = lonsum / n
  120. if (newobs(ii) % lon > 180.0d0) then
  121. newobs(ii) % lon = newobs(ii) % lon - 360.0d0
  122. end if
  123. end if
  124. newobs(ii) % lat = latsum / n
  125. newobs(ii) % depth = depthsum / n
  126. newobs(ii) % ipiv = iprev
  127. newobs(ii) % jpiv = jprev
  128. newobs(ii) % ns = 0 ! not 100% sure
  129. newobs(ii) % a1 = a1sum / n
  130. newobs(ii) % a2 = a2sum / n
  131. newobs(ii) % a3 = a3sum / n
  132. newobs(ii) % a4 = a4sum / n
  133. newobs(ii) % status = .true.
  134. newobs(ii) % i_orig_grid = -1
  135. newobs(ii) % j_orig_grid = -1
  136. newobs(ii) % h = n
  137. newobs(ii) % date = age
  138. newobs(ii) % orig_id = ids(1) ! ID of the first ob
  139. nmax = max(n, nmax)
  140. if (TEST) then
  141. write(101, '(a, g10.3)') 'total # of obs = ', n
  142. write(101, '(a, i6)') ' index = ', ii
  143. write(101, '(a, g10.3)') ' d = ', newobs(ii) % d
  144. write(101, '(a, g10.3)') ' var = ', newobs(ii) % var
  145. write(101, '(a, g10.3)') ' lon = ', newobs(ii) % lon
  146. write(101, '(a, g10.3)') ' lat = ', newobs(ii) % lat
  147. write(101, '(a, i4)') ' ipiv = ', newobs(ii) % ipiv
  148. write(101, '(a, i4)') ' jpiv = ', newobs(ii) % jpiv
  149. write(101, '(a, g10.3)') ' depth = ', newobs(ii) % depth
  150. write(101, '(a, g10.3)') ' a1 = ', newobs(ii) % a1
  151. write(101, '(a, g10.3)') ' a2 = ', newobs(ii) % a2
  152. write(101, '(a, g10.3)') ' a3 = ', newobs(ii) % a3
  153. write(101, '(a, g10.3)') ' a4 = ', newobs(ii) % a4
  154. write(101, '(a)') '---'
  155. call superobs_dump(trim(obstag), ii, ids, int(n))
  156. end if
  157. end if
  158. if (o > nobs_now) then
  159. exit
  160. end if
  161. ii = ii + 1
  162. ii_now = ii_now + 1
  163. if (TEST) then
  164. write(101, '(a, i6)') 'new superob, index = ', ii
  165. end if
  166. n = 0.0
  167. valsum = 0.0d0
  168. valsqsum = 0.0d0
  169. varinvsum = 0.0d0
  170. lonsum = 0.0d0
  171. latsum = 0.0d0
  172. depthsum = 0.0
  173. a1sum = 0.0d0
  174. a2sum = 0.0d0
  175. a3sum = 0.0d0
  176. a4sum = 0.0d0
  177. valmax = -1.0d+20
  178. valmin = 1.0d+20
  179. iprev = thisob % ipiv
  180. jprev = thisob % jpiv
  181. kprev = kpiv(sorted(o))
  182. nlon_pos = 0
  183. nlon_neg = 0
  184. lonsum_abs = 0.0d0
  185. end if
  186. n = n + 1.0
  187. valsum = valsum + thisob % d
  188. valsqsum = valsqsum + (thisob % d) ** 2
  189. varinvsum = varinvsum + 1.0 / thisob % var
  190. lonsum = lonsum + thisob % lon
  191. lonsum_abs = lonsum_abs + abs(thisob % lon)
  192. if (thisob % lon >= 0.0) then
  193. nlon_pos = nlon_pos + 1
  194. else
  195. nlon_neg = nlon_neg + 1
  196. end if
  197. latsum = latsum + thisob % lat
  198. depthsum = depthsum + thisob % depth
  199. a1sum = a1sum + thisob % a1
  200. a2sum = a2sum + thisob % a2
  201. a3sum = a3sum + thisob % a3
  202. a4sum = a4sum + thisob % a4
  203. valmin = min(valmin, thisob % d)
  204. valmax = max(valmax, thisob % d)
  205. ids(int(n)) = thisob % orig_id;
  206. if (TEST) then
  207. write(101, '(a, i6)') ' obs index = ', sorted(o)
  208. write(101, '(a, g10.3)') ' d = ', thisob % d
  209. write(101, '(a, g10.3)') ' var = ', thisob % var
  210. write(101, '(a, g10.3)') ' lon = ', thisob % lon
  211. write(101, '(a, g10.3)') ' lat = ', thisob % lat
  212. write(101, '(a, i4)') ' ipiv = ', thisob % ipiv
  213. write(101, '(a, i4)') ' jpiv = ', thisob % jpiv
  214. write(101, '(a, g10.3)') ' depth = ', thisob % depth
  215. write(101, '(a, g10.3)') ' a1 = ', thisob % a1
  216. write(101, '(a, g10.3)') ' a2 = ', thisob % a2
  217. write(101, '(a, g10.3)') ' a3 = ', thisob % a3
  218. write(101, '(a, g10.3)') ' a4 = ', thisob % a4
  219. end if
  220. end do ! obs for this age
  221. print *, ' nsuperobs =', ii_now
  222. end do ! age
  223. if (TEST) then
  224. close(101)
  225. end if
  226. nnewobs = ii
  227. print *, 'Superobing("', trim(obstag), '"):'
  228. print *, ' ', nobs, 'observations ->', nnewobs, 'observations'
  229. print *, ' max # of obs found in a grid cell =', int(nmax)
  230. end subroutine superob
  231. function z2k(z)
  232. real, intent(in), dimension(:) :: z
  233. integer, dimension(size(z)) :: z2k
  234. integer :: i, nz
  235. nz = size(z)
  236. do i = 1, nz
  237. if (z(i) < 3.0d0) then
  238. z2k(i) = 1
  239. elseif (z(i) < 6.0d0) then
  240. z2k(i) = 2
  241. elseif (z(i) < 10.0d0) then
  242. z2k(i) = 3
  243. elseif (z(i) < 100.0d0) then
  244. z2k(i) = int(z(i) / 10.0d0) + 3
  245. elseif (z(i) < 1000.0d0) then
  246. z2k(i) = int(z(i) / 25.0d0) + 9
  247. else
  248. z2k(i) = int(z(i) / 50.0d0) + 29
  249. end if
  250. end do
  251. end function z2k
  252. subroutine superobs_dump(tag, id, ids, n)
  253. use nfw_mod
  254. character(*) :: tag
  255. integer, intent(in) :: id
  256. integer, intent(in) :: ids(n)
  257. integer, intent(in) :: n
  258. character(STRLEN) :: fname
  259. character(64) :: dname
  260. character(64) :: vname
  261. integer :: ncid, did(1), vid
  262. if (id > NF_MAX_DIMS) then
  263. return
  264. end if
  265. write(fname, '(a, a, a)') 'superobs-', trim(tag), '.nc'
  266. if (id == 1) then
  267. print *, 'dumping obs ids for each superob to "', trim(fname), '"'
  268. call nfw_create(fname, nf_clobber, ncid)
  269. else
  270. call nfw_open(fname, nf_write, ncid)
  271. call nfw_redef(fname, ncid)
  272. end if
  273. write(dname, '(a,i0)') 'd', id
  274. call nfw_def_dim(fname, ncid, trim(dname), n, did(1))
  275. write(vname, '(a,i0)') 'v', id
  276. call nfw_def_var(fname, ncid, trim(vname), nf_int, 1, did(1), vid)
  277. call nfw_enddef(fname, ncid)
  278. call nfw_put_var_int(fname, ncid, vid, ids)
  279. call nfw_close(fname, ncid)
  280. end subroutine superobs_dump
  281. end module m_superobs