EnKF.F90 10.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330
  1. ! File: EnKF.F90
  2. !
  3. ! Created: ???
  4. !
  5. ! Last modified: 20/04/2010
  6. !
  7. ! Purpose: Main program for EnKF analysis
  8. !
  9. ! Description: The workflow is as follows:
  10. ! -- read model parameters
  11. ! -- read obs
  12. ! -- conduct necessary pre-processing of obs (superobing)
  13. ! -- calculate ensemble observations
  14. ! -- calculate X5
  15. ! -- update the ensemble
  16. !
  17. ! Modifications:
  18. ! 20/9/2011 PS:
  19. ! Modified code to allow individual inflations for each of
  20. ! `NFIELD' fields updated in a batch - thanks to Ehouarn Simon
  21. ! for spotting this inconsistency
  22. ! 6/8/2010 PS:
  23. ! Small changes in calls to calc_X5() and update_fields() to
  24. ! reflect changes in interfaces.
  25. ! 6/7/2010 PS:
  26. ! Moved point output to a separate module m_point2nc.F90
  27. ! 25/5/2010 PS:
  28. ! Added inflation as a 4th command line argument
  29. ! 20/5/2010 PS:
  30. ! Set NFIELD = 4. This requires 4 GB per node in TOPAZ and
  31. ! "medium" memory model on Hexagon (a single allocation for a
  32. ! variable over 2GB)
  33. ! 20/4/2010 PS:
  34. ! Set NFIELD = 4. This will require 2 GB per node in TOPAZ.
  35. ! Thanks to Alok Gupta for hinting this possibility.
  36. ! 10/4/2010 PS:
  37. ! Moved variable `field' from real(8) to real(4);
  38. ! set NFIELD = 2.
  39. ! Prior history:
  40. ! Not documented.
  41. ! 15/4/2016 Francois Massonnet (FM): Make changes to be
  42. ! NEMO-compliant. Targeted for NEMO3.6 at BSC,
  43. ! Barcelona, but based on previous experience
  44. ! at UCL and on work from Chris Konig-Beaty [CKB]
  45. program EnKF
  46. #if defined(QMPI)
  47. use qmpi
  48. #else
  49. use qmpi_fake
  50. #endif
  51. use m_parameters
  52. use distribute
  53. use mod_measurement
  54. use m_get_mod_grid
  55. use m_get_mod_nrens
  56. use m_get_mod_xyz ! Added by Francois Massonnet [FM] May 2013 and Apr 2016 !
  57. use m_obs
  58. use m_local_analysis
  59. use m_prep_4_EnKF
  60. use m_set_random_seed2
  61. !use m_get_mod_fld ! Taken out and simplified by m_io_mod_fld
  62. !use m_put_mod_fld
  63. use m_io_mod_fld ![CKB,FM] was: m_get_mod_fld and m_put_mod_fld
  64. use mod_analysisfields
  65. use m_parse_blkdat
  66. use m_random
  67. use m_point2nc
  68. implicit none
  69. character(*), parameter :: ENKF_VERSION = "2.11"
  70. integer, external :: iargc
  71. ! NFIELD is the number of fields (x N) passed for the update during a call to
  72. ! update_fields(). In TOPAZ4 NFIELD = 2 if there is 1 GB of RAM per node, and
  73. ! NFIELD = 4 if there are 2 GB of RAM. Higher value of NFIELD reduces the
  74. ! number of times X5tmp.uf is read from disk, which is the main bottleneck
  75. ! for the analysis time right now.
  76. !
  77. integer, parameter :: NFIELD = 8
  78. character(512) :: options
  79. integer :: nrens
  80. integer, allocatable, dimension(:) :: enslist ! [FM] List of existing
  81. ! ensemble members
  82. real, allocatable, dimension(:,:) :: modlon, modlat, depths, readfld
  83. real, allocatable, dimension(:,:) :: S ! ensemble observations HE
  84. real, allocatable, dimension(:) :: d ! d - Hx
  85. integer k, m
  86. ! "New" variables used in the parallelization
  87. integer, dimension(:,:), allocatable :: nlobs_array
  88. real(4), allocatable :: fld(:,:)
  89. real(8) rtc, time0, time1, time2
  90. ! Additional fields
  91. character(len=3) :: cmem
  92. character(len=80) :: memfile
  93. integer :: fieldcounter
  94. character(100) :: text_string
  95. real :: rdummy
  96. integer :: idm, jdm, kdm
  97. real :: mindx
  98. real :: meandx
  99. integer :: m1, m2, nfields
  100. real :: infls(NFIELD)
  101. #if defined(QMPI)
  102. call start_mpi()
  103. #endif
  104. ! Read the characteristics of the assimilation to be carried out.
  105. if (iargc() /= 1) then
  106. print *, 'Usage: EnKF <parameter file>'
  107. print *, ' EnKF -h'
  108. print *, 'Options:'
  109. print *, ' -h -- describe parameter fie format'
  110. call stop_mpi()
  111. else
  112. call getarg(1, options)
  113. if (trim(options) == "-h") then
  114. call prm_describe()
  115. call stop_mpi()
  116. end if
  117. end if
  118. if (master) then
  119. print *
  120. print '(a, a)', ' EnKF version ', ENKF_VERSION
  121. print *
  122. end if
  123. call prm_read()
  124. call prm_print()
  125. ! get model dimensions
  126. !
  127. ! Change FM May 2013. Goal is to avoid using parse_blkdat that requires a
  128. ! file with unknown format
  129. !call parse_blkdat('idm ', 'integer', rdummy, idm)
  130. !call parse_blkdat('jdm ', 'integer', rdummy, jdm)
  131. !call parse_blkdat('kdm ', 'integer', rdummy, kdm)
  132. CALL get_mod_xyz(idm, jdm, kdm)
  133. WRITE(*,*), 'The model dimensions are ', idm,jdm,kdm
  134. ! End Change FM May 2013.
  135. allocate(modlon(idm, jdm))
  136. allocate(readfld(idm, jdm))
  137. allocate(modlat(idm, jdm))
  138. allocate(depths(idm, jdm))
  139. allocate(nlobs_array(idm, jdm))
  140. ! get model grid
  141. !
  142. call get_mod_grid(modlon, modlat, depths, mindx, meandx, idm, jdm)
  143. ! set a variable random seed
  144. !
  145. call set_random_seed2
  146. ! initialise point output
  147. !
  148. call p2nc_init
  149. time0 = rtc()
  150. ! read measurements
  151. !
  152. if (master) then
  153. print *, 'EnKF: reading observations'
  154. end if
  155. call obs_readobs
  156. if (master) then
  157. print '(a, i6)', ' # of obs = ', nobs
  158. print '(a, a, a, e10.3, a, e10.3)', ' first obs = "', trim(obs(1) % id),&
  159. '", v = ', obs(1) % d, ', var = ', obs(1) % var
  160. print '(a, a, a, e10.3, a, e10.3)', ' last obs = "', trim(obs(nobs) % id),&
  161. '", v = ', obs(nobs) % d, ', var = ', obs(nobs) % var
  162. end if
  163. if (master) then
  164. print *
  165. end if
  166. ! read ensemble size and store in A
  167. !
  168. ! [CKB,FM] changed
  169. call get_mod_nrens(nrens)
  170. allocate( enslist(nrens) )
  171. call get_mod_nrens(nrens, enslist)
  172. ! end [CKB, FM]
  173. if (master) then
  174. print '(a, i4, a)', ' EnKF: ', nrens, ' ensemble members found'
  175. end if
  176. if (ENSSIZE > 0) then
  177. ENSSIZE = min(nrens, ENSSIZE)
  178. else
  179. ENSSIZE = nrens
  180. end if
  181. if (master) then
  182. print '(a, i4, a)', ' EnKF: ', ENSSIZE, ' ensemble members used'
  183. end if
  184. if (master) then
  185. print *
  186. end if
  187. ! PS - preprocess the obs using the information about the ensemble fields
  188. ! here (if necessary), before running prep_4_EnKF(). This is necessary e.g.
  189. ! for assimilating in-situ data because of the dynamic vertical geometry in
  190. ! HYCOM
  191. !
  192. call obs_prepareobs
  193. allocate(S(nobs, ENSSIZE), d(nobs))
  194. call prep_4_EnKF(ENSSIZE,enslist, d, S, depths, meandx / 1000.0, idm, jdm, kdm)
  195. if (master) then
  196. print *, 'EnKF: finished initialisation, time = ', rtc() - time0
  197. end if
  198. ! (no parallelization was required before this point)
  199. time1 = rtc()
  200. allocate(X5(ENSSIZE, ENSSIZE, idm))
  201. allocate(X5check(ENSSIZE, ENSSIZE, idm))
  202. call calc_X5(ENSSIZE, modlon, modlat, depths, mindx, meandx, d, S,&
  203. LOCRAD, RFACTOR2, nlobs_array, idm, jdm)
  204. deallocate(d, S, X5check)
  205. if (master) then
  206. print *, 'EnKF: finished calculation of X5, time = ', rtc() - time0
  207. end if
  208. allocate(fld(idm * jdm, ENSSIZE * NFIELD))
  209. #if defined(QMPI)
  210. call barrier()
  211. #endif
  212. ! get fieldnames and fieldlevels
  213. !
  214. call get_analysisfields()
  215. call distribute_iterations(numfields)
  216. #if defined(QMPI)
  217. call barrier() !KAL - just for "niceness" of output
  218. #endif
  219. time2 = rtc()
  220. do m1 = my_first_iteration, my_last_iteration, NFIELD
  221. m2 = min(my_last_iteration, m1 + NFIELD - 1)
  222. nfields = m2 - m1 + 1
  223. do m = m1, m2
  224. print '(a, i2, a, i3, a, a6, a, i3, a, f11.0)',&
  225. "I am ", qmpi_proc_num, ', m = ', m, ", field = ",&
  226. fieldnames(m), ", k = ", fieldlevel(m), ", time = ",&
  227. rtc() - time2
  228. do k = 1, ENSSIZE
  229. write(cmem, '(i3.3)') k
  230. memfile = 'forecast' // cmem
  231. !call get_mod_fld_new(trim(memfile), readfld, k, fieldnames(m),&
  232. ! fieldlevel(m), 1, idm, jdm)
  233. ! [CKB,FM]
  234. call io_mod_fld(readfld, k, enslist,fieldnames(m),fieldtype(m), &
  235. fieldlevel(m), 1, idm, jdm, 'get',FLOAT(obs(1)%date))
  236. ! end CKB,FM
  237. ! reshaping and conversion to real(4)
  238. fld(:, ENSSIZE * (m - m1) + k) = reshape(readfld, (/idm * jdm/))
  239. end do
  240. call p2nc_storeforecast(idm, jdm, ENSSIZE, numfields, m, fld(:, ENSSIZE * (m - m1) + 1 : ENSSIZE * (m + 1 - m1)))
  241. infls(m - m1 + 1) = prm_getinfl(trim(fieldnames(m)));
  242. end do
  243. call update_fields(idm, jdm, ENSSIZE, nfields, nlobs_array, depths,&
  244. fld(1,1), infls)
  245. do m = m1, m2
  246. fieldcounter = (m - my_first_iteration) + 1
  247. do k = 1, ENSSIZE
  248. write(cmem,'(i3.3)') k
  249. memfile = 'analysis' // cmem
  250. ! reshaping and conversion to real(8)
  251. readfld = reshape(fld(:, ENSSIZE * (m - m1) + k), (/idm, jdm/))
  252. write(text_string, '(a, i3.3)') '_proc', qmpi_proc_num
  253. !call put_mod_fld(trim(memfile) // trim(text_string), readfld, k,&
  254. ! fieldnames(m), fieldlevel(m), 1, fieldcounter, idm, jdm)
  255. ! [FM,CKB]
  256. call io_mod_fld(readfld, k, enslist, fieldnames(m), fieldtype(m), &
  257. fieldlevel(m), 1, idm, jdm, 'put',FLOAT(obs(1)%date))
  258. ! end FM,CKB
  259. end do
  260. end do
  261. end do
  262. deallocate(X5)
  263. deallocate(fld)
  264. call p2nc_writeforecast
  265. ! Barrier only necessary for timings
  266. #if defined(QMPI)
  267. call barrier()
  268. #endif
  269. if (master) then
  270. print *, 'EnKF: time for initialization = ', time1 - time0
  271. print *, 'EnKF: time for X5 calculation = ', time2 - time1
  272. print *, 'EnKF: time for ensemble update = ', rtc() - time2
  273. print *, 'EnKF: total time = ', rtc() - time0
  274. end if
  275. print *, 'EnKF: Finished'
  276. call stop_mpi()
  277. end program EnKF
  278. #if defined(_G95_)
  279. ! not tested! - PS
  280. !
  281. real function rtc()
  282. integer :: c
  283. call system_clock(count=c)
  284. rtc = dfloat(c)
  285. end function rtc
  286. #endif