m_point2nc.F90 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339
  1. ! File: m_point2nc.F90
  2. !
  3. ! Created: 6 July 2010
  4. !
  5. ! Last modified: 6/7/2010
  6. !
  7. ! Author: Pavel Sakov
  8. ! NERSC
  9. !
  10. ! Purpose: Output of assimilation related information for selected points
  11. ! to files in NetCDF format, 1 file per point.
  12. !
  13. ! Description: This module reads a list of points from a file "point2nc.txt"
  14. ! in the working NetCDF directory. It then dumps the
  15. ! assimilation related information for these points in NetCDF
  16. ! format to files named enkf_III,JJJ.nc, where III and JJJ - i
  17. ! and j grid coordinates.
  18. !
  19. ! Modifications: PS 4/8/2010 "point2nc.txt" now allows comments etc. E.g. put
  20. ! "#" in front of an entry to comment it out.
  21. module m_point2nc
  22. use m_parameters
  23. implicit none
  24. integer, private :: FID = 31
  25. integer, parameter, private :: STRLEN = 512
  26. public p2nc_init
  27. public p2nc_testthiscell
  28. public p2nc_writeobs
  29. public p2nc_storeforecast
  30. public p2nc_writeforecast
  31. integer, private :: npoints
  32. integer, allocatable, dimension(:), private :: icoords, jcoords
  33. real(4), allocatable, dimension(:, :, :) :: forecast
  34. contains
  35. ! Initialise the point output.
  36. !
  37. subroutine p2nc_init()
  38. #if defined (QMPI)
  39. use qmpi
  40. #else
  41. use qmpi_fake
  42. #endif
  43. character(STRLEN) :: line
  44. integer :: iostatus
  45. integer :: i, j, n
  46. npoints = 0
  47. open(FID, file = trim(POINTFNAME), status = 'old', iostat = iostatus)
  48. if (iostatus /= 0) then
  49. if (master) then
  50. print *, 'WARNING: could not open "', trim(POINTFNAME), '" for reading'
  51. print *, ' no point output will be performed'
  52. end if
  53. return
  54. end if
  55. do while (.true.)
  56. read(FID, '(a)', iostat = iostatus) line
  57. if (iostatus == 0) then
  58. read(line, *, iostat = iostatus) i, j
  59. if (iostatus == 0) then
  60. npoints = npoints + 1
  61. end if
  62. else
  63. exit
  64. end if
  65. end do
  66. close(FID)
  67. if (master) then
  68. print '(a, i3, a)', ' p2nc: ', npoints, ' points specified'
  69. end if
  70. allocate(icoords(npoints), jcoords(npoints))
  71. open(FID, file = trim(POINTFNAME), status = 'old', iostat = iostatus)
  72. if (iostatus /= 0) then
  73. print *, 'ERROR: point2nc: I/O problem'
  74. stop
  75. end if
  76. n = 0
  77. do while (n < npoints)
  78. read(FID, '(a)', iostat = iostatus) line
  79. if (iostatus == 0) then
  80. read(line, *, iostat = iostatus) i, j
  81. if (iostatus == 0) then
  82. n = n + 1
  83. icoords(n) = i
  84. jcoords(n) = j
  85. if (master) then
  86. print '(a, i3, a, i4, a, i4)', ' point', n, ': i =', i, ', j =', j
  87. end if
  88. end if
  89. end if
  90. end do
  91. close(FID)
  92. if (master) then
  93. print *
  94. end if
  95. end subroutine p2nc_init
  96. ! Test if the output is requested for the point (i, j)
  97. !
  98. function p2nc_testthiscell(i, j)
  99. logical :: p2nc_testthiscell
  100. integer, intent(in) :: i, j
  101. integer :: p
  102. p2nc_testthiscell = .false.
  103. do p = 1, npoints
  104. if (i == icoords(p) .and. j == jcoords(p)) then
  105. p2nc_testthiscell = .true.
  106. return
  107. end if
  108. end do
  109. end function p2nc_testthiscell
  110. ! Write the assimilation parameters (local observations and the X5 matrices)
  111. ! to the point output files.
  112. !
  113. subroutine p2nc_writeobs(i, j, nlobs, nens, X5, lon, lat, depth, rfactor,&
  114. ids, lobs, Hx, S, ss, lfactors)
  115. use mod_measurement
  116. use m_obs
  117. use nfw_mod
  118. integer, intent(in) :: i, j, nlobs, nens
  119. real, intent(in) :: X5(nens, nens)
  120. real, intent(in) :: lon(1), lat(1), depth(1)
  121. real, intent(in), optional :: rfactor(1)
  122. integer, intent(in), optional :: ids(nlobs)
  123. type(measurement), intent(in), optional :: lobs(nlobs)
  124. real, intent(in), optional :: Hx(nlobs)
  125. real, intent(in), optional :: S(nlobs, nens)
  126. real, intent(in), optional :: ss(nlobs), lfactors(nlobs)
  127. character(STRLEN) :: fname
  128. character(STRLEN) :: typename
  129. integer :: ncid
  130. integer :: dids(2)
  131. integer :: vid_ids, vid_lon, vid_lat, vid_val, vid_var, vid_hx, vid_s, vid_x5
  132. integer :: vid_a1, vid_a2, vid_a3, vid_a4, vid_otype, vid_ss, vid_lfactors
  133. integer :: otype(nlobs)
  134. integer :: o, ot
  135. write(fname, '(a, i0.3, ",", i0.3, ".nc")') 'enkf_', i, j
  136. call nfw_create(fname, nf_write, ncid)
  137. call nfw_def_dim(fname, ncid, 'p', nlobs, dids(2))
  138. call nfw_def_dim(fname, ncid, 'm', nens, dids(1))
  139. if (nlobs > 0) then
  140. call nfw_def_var(fname, ncid, 'obs_ids', nf_int, 1, dids(2), vid_ids)
  141. call nfw_def_var(fname, ncid, 'Hx', nf_double, 1, dids(2), vid_hx)
  142. call nfw_def_var(fname, ncid, 'lon', nf_double, 1, dids(2), vid_lon)
  143. call nfw_def_var(fname, ncid, 'lat', nf_double, 1, dids(2), vid_lat)
  144. call nfw_def_var(fname, ncid, 'obs_val', nf_double, 1, dids(2), vid_val)
  145. call nfw_def_var(fname, ncid, 'obs_var', nf_double, 1, dids(2), vid_var)
  146. call nfw_def_var(fname, ncid, 'a1', nf_double, 1, dids(2), vid_a1)
  147. call nfw_def_var(fname, ncid, 'a2', nf_double, 1, dids(2), vid_a2)
  148. call nfw_def_var(fname, ncid, 'a3', nf_double, 1, dids(2), vid_a3)
  149. call nfw_def_var(fname, ncid, 'a4', nf_double, 1, dids(2), vid_a4)
  150. call nfw_def_var(fname, ncid, 'obs_type', nf_int, 1, dids(2), vid_otype)
  151. call nfw_def_var(fname, ncid, 'S', nf_double, 2, dids, vid_s)
  152. call nfw_def_var(fname, ncid, 's', nf_double, 1, dids(2), vid_ss)
  153. call nfw_def_var(fname, ncid, 'lfactors', nf_double, 1, dids(2), vid_lfactors)
  154. end if
  155. dids(2) = dids(1)
  156. call nfw_def_var(fname, ncid, 'X5', nf_double, 2, dids, vid_x5)
  157. call nfw_put_att_double(fname, ncid, nf_global, 'lon', nf_double, 1, lon)
  158. call nfw_put_att_double(fname, ncid, nf_global, 'lat', nf_double, 1, lat)
  159. call nfw_put_att_double(fname, ncid, nf_global, 'depth', nf_double, 1, depth)
  160. call nfw_put_att_double(fname, ncid, nf_global, 'rfactor', nf_double, 1, rfactor)
  161. do ot = 1, nuobs
  162. write(typename, '(a, i1)') 'obstype', ot
  163. call nfw_put_att_text(fname, ncid, nf_global, typename, len_trim(unique_obs(ot)), trim(unique_obs(ot)))
  164. end do
  165. call nfw_enddef(fname, ncid)
  166. if (nlobs > 0) then
  167. call nfw_put_var_double(fname, ncid, vid_hx, Hx)
  168. call nfw_put_var_int(fname, ncid, vid_ids, ids)
  169. call nfw_put_var_double(fname, ncid, vid_lon, lobs % lon)
  170. call nfw_put_var_double(fname, ncid, vid_lat, lobs % lat)
  171. call nfw_put_var_double(fname, ncid, vid_val, lobs % d)
  172. call nfw_put_var_double(fname, ncid, vid_var, lobs % var)
  173. call nfw_put_var_double(fname, ncid, vid_a1, lobs % a1)
  174. call nfw_put_var_double(fname, ncid, vid_a2, lobs % a2)
  175. call nfw_put_var_double(fname, ncid, vid_a3, lobs % a3)
  176. call nfw_put_var_double(fname, ncid, vid_a4, lobs % a4)
  177. otype = 0
  178. do o = 1, nlobs
  179. do ot = 1, nuobs
  180. if (trim(lobs(o) % id) == trim(unique_obs(ot))) then
  181. otype(o) = ot
  182. end if
  183. end do
  184. end do
  185. call nfw_put_var_int(fname, ncid, vid_otype, otype)
  186. call nfw_put_var_double(fname, ncid, vid_s, transpose(S))
  187. call nfw_put_var_double(fname, ncid, vid_ss, ss)
  188. call nfw_put_var_double(fname, ncid, vid_lfactors, lfactors)
  189. end if
  190. call nfw_put_var_double(fname, ncid, vid_x5, transpose(X5))
  191. call nfw_close(fname, ncid)
  192. end subroutine p2nc_writeobs
  193. ! Store the values of the forecast field No. `fid' in each output point to
  194. ! the variable `forecast'.
  195. !
  196. subroutine p2nc_storeforecast(ni, nj, nrens, nfields, fid, field)
  197. integer, intent(in) :: ni, nj ! size of grid
  198. integer, intent(in) :: nrens
  199. integer, intent(in) :: nfields
  200. integer, intent(in) :: fid
  201. real(4), dimension(ni * nj, nrens), intent(in) :: field
  202. integer :: n
  203. if (npoints == 0) then
  204. return
  205. end if
  206. if (.not. allocated(forecast)) then
  207. allocate(forecast(nrens, npoints, nfields))
  208. end if
  209. do n = 1, npoints
  210. forecast(:, n, fid) = field((jcoords(n) - 1) * ni + icoords(n), :)
  211. end do
  212. end subroutine p2nc_storeforecast
  213. ! This procedure consolidates all forecast fields for each output point
  214. ! together in the variable `forecast' on the master node, and then writes
  215. ! them to the appropriate NetCDF files.
  216. !
  217. subroutine p2nc_writeforecast
  218. #if defined (QMPI)
  219. use qmpi
  220. #else
  221. use qmpi_fake
  222. #endif
  223. use distribute
  224. use nfw_mod
  225. use mod_analysisfields
  226. implicit none
  227. character(STRLEN) :: fname
  228. integer :: p, k, nf
  229. character(8) :: varname
  230. integer kstart
  231. integer ncid, dids(2), varid, nf2
  232. #if defined(QMPI)
  233. if (.not. master) then
  234. call send(forecast(:, :, my_first_iteration : my_last_iteration), 0, 0)
  235. return ! leave writing to master
  236. else
  237. do p = 2, qmpi_num_proc ! here p is the MPI node ID
  238. call receive(forecast(:, :, first_iteration(p) : last_iteration(p)), p - 1, 0)
  239. end do
  240. end if
  241. #endif
  242. ! only master keeps working here
  243. !
  244. do p = 1, npoints
  245. write(fname, '(a, i0.3, ",", i0.3, ".nc")') 'enkf_', icoords(p), jcoords(p)
  246. call nfw_open(fname, nf_write, ncid)
  247. call nfw_redef(fname, ncid)
  248. call nfw_inq_dimid(fname, ncid, 'm', dids(1))
  249. call nfw_enddef(fname, ncid)
  250. kstart = -1
  251. do k = 1, numfields
  252. if (kstart == -1) then
  253. kstart = k
  254. varname = fieldnames(k)
  255. end if
  256. ! check if there are more fields for this variable
  257. !
  258. if (k < numfields .and. fieldnames(k + 1) == varname) then
  259. cycle
  260. end if
  261. ! this is the last field for this variable - write the variable
  262. !
  263. nf = k - kstart + 1
  264. call nfw_redef(fname, ncid)
  265. if (nf == 1) then
  266. call nfw_def_var(fname, ncid, trim(varname), nf_float, 1, dids(1), varid)
  267. else
  268. if (.not. nfw_dim_exists(ncid, 'k')) then
  269. call nfw_def_dim(fname, ncid, 'k', nf, dids(2))
  270. else
  271. call nfw_inq_dimid(fname, ncid, 'k', dids(2))
  272. call nfw_inq_dimlen(fname, ncid, dids(2), nf2)
  273. if (nf /= nf2) then
  274. print *, 'ERROR: p2nc_writeforecast(): varname = "', trim(varname),&
  275. '", # levels = ', nf, '# levels in "', trim(fname), '" =', nf2
  276. print *, 'ERROR: p2nc_writeforecast(): returning'
  277. end if
  278. end if
  279. call nfw_def_var(fname, ncid, trim(varname), nf_float, 2, dids, varid)
  280. end if
  281. call nfw_enddef(fname, ncid)
  282. call nfw_put_var_real(fname, ncid, varid, forecast(:, p, kstart : kstart + nf - 1))
  283. kstart = -1
  284. end do
  285. call nfw_close(fname, ncid)
  286. end do
  287. end subroutine p2nc_writeforecast
  288. end module m_point2nc