EnKF.F90 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334
  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. #if defined(_G95_)
  71. integer, intrinsic :: iargc
  72. #else
  73. integer, external :: iargc
  74. #endif
  75. ! NFIELD is the number of fields (x N) passed for the update during a call to
  76. ! update_fields(). In TOPAZ4 NFIELD = 2 if there is 1 GB of RAM per node, and
  77. ! NFIELD = 4 if there are 2 GB of RAM. Higher value of NFIELD reduces the
  78. ! number of times X5tmp.uf is read from disk, which is the main bottleneck
  79. ! for the analysis time right now.
  80. !
  81. integer, parameter :: NFIELD = 8
  82. character(512) :: options
  83. integer :: nrens
  84. integer, allocatable, dimension(:) :: enslist ! [FM] List of existing
  85. ! ensemble members
  86. real, allocatable, dimension(:,:) :: modlon, modlat, depths, readfld
  87. real, allocatable, dimension(:,:) :: S ! ensemble observations HE
  88. real, allocatable, dimension(:) :: d ! d - Hx
  89. integer k, m
  90. ! "New" variables used in the parallelization
  91. integer, dimension(:,:), allocatable :: nlobs_array
  92. real(4), allocatable :: fld(:,:)
  93. real(8) rtc, time0, time1, time2
  94. ! Additional fields
  95. character(len=3) :: cmem
  96. character(len=80) :: memfile
  97. integer :: fieldcounter
  98. character(100) :: text_string
  99. real :: rdummy
  100. integer :: idm, jdm, kdm
  101. real :: mindx
  102. real :: meandx
  103. integer :: m1, m2, nfields
  104. real :: infls(NFIELD)
  105. #if defined(QMPI)
  106. call start_mpi()
  107. #endif
  108. ! Read the characteristics of the assimilation to be carried out.
  109. if (iargc() /= 1) then
  110. print *, 'Usage: EnKF <parameter file>'
  111. print *, ' EnKF -h'
  112. print *, 'Options:'
  113. print *, ' -h -- describe parameter fie format'
  114. call stop_mpi()
  115. else
  116. call getarg(1, options)
  117. if (trim(options) == "-h") then
  118. call prm_describe()
  119. call stop_mpi()
  120. end if
  121. end if
  122. if (master) then
  123. print *
  124. print '(a, a)', ' EnKF version ', ENKF_VERSION
  125. print *
  126. end if
  127. call prm_read()
  128. call prm_print()
  129. ! get model dimensions
  130. !
  131. ! Change FM May 2013. Goal is to avoid using parse_blkdat that requires a
  132. ! file with unknown format
  133. !call parse_blkdat('idm ', 'integer', rdummy, idm)
  134. !call parse_blkdat('jdm ', 'integer', rdummy, jdm)
  135. !call parse_blkdat('kdm ', 'integer', rdummy, kdm)
  136. CALL get_mod_xyz(idm, jdm, kdm)
  137. WRITE(*,*), 'The model dimensions are ', idm,jdm,kdm
  138. ! End Change FM May 2013.
  139. allocate(modlon(idm, jdm))
  140. allocate(readfld(idm, jdm))
  141. allocate(modlat(idm, jdm))
  142. allocate(depths(idm, jdm))
  143. allocate(nlobs_array(idm, jdm))
  144. ! get model grid
  145. !
  146. call get_mod_grid(modlon, modlat, depths, mindx, meandx, idm, jdm)
  147. ! set a variable random seed
  148. !
  149. call set_random_seed2
  150. ! initialise point output
  151. !
  152. call p2nc_init
  153. time0 = rtc()
  154. ! read measurements
  155. !
  156. if (master) then
  157. print *, 'EnKF: reading observations'
  158. end if
  159. call obs_readobs
  160. if (master) then
  161. print '(a, i6)', ' # of obs = ', nobs
  162. print '(a, a, a, e10.3, a, e10.3)', ' first obs = "', trim(obs(1) % id),&
  163. '", v = ', obs(1) % d, ', var = ', obs(1) % var
  164. print '(a, a, a, e10.3, a, e10.3)', ' last obs = "', trim(obs(nobs) % id),&
  165. '", v = ', obs(nobs) % d, ', var = ', obs(nobs) % var
  166. end if
  167. if (master) then
  168. print *
  169. end if
  170. ! read ensemble size and store in A
  171. !
  172. ! [CKB,FM] changed
  173. call get_mod_nrens(nrens)
  174. allocate( enslist(nrens) )
  175. call get_mod_nrens(nrens, enslist)
  176. ! end [CKB, FM]
  177. if (master) then
  178. print '(a, i4, a)', ' EnKF: ', nrens, ' ensemble members found'
  179. end if
  180. if (ENSSIZE > 0) then
  181. ENSSIZE = min(nrens, ENSSIZE)
  182. else
  183. ENSSIZE = nrens
  184. end if
  185. if (master) then
  186. print '(a, i4, a)', ' EnKF: ', ENSSIZE, ' ensemble members used'
  187. end if
  188. if (master) then
  189. print *
  190. end if
  191. ! PS - preprocess the obs using the information about the ensemble fields
  192. ! here (if necessary), before running prep_4_EnKF(). This is necessary e.g.
  193. ! for assimilating in-situ data because of the dynamic vertical geometry in
  194. ! HYCOM
  195. !
  196. call obs_prepareobs
  197. allocate(S(nobs, ENSSIZE), d(nobs))
  198. call prep_4_EnKF(ENSSIZE,enslist, d, S, depths, meandx / 1000.0, idm, jdm, kdm)
  199. if (master) then
  200. print *, 'EnKF: finished initialisation, time = ', rtc() - time0
  201. end if
  202. ! (no parallelization was required before this point)
  203. time1 = rtc()
  204. allocate(X5(ENSSIZE, ENSSIZE, idm))
  205. allocate(X5check(ENSSIZE, ENSSIZE, idm))
  206. call calc_X5(ENSSIZE, modlon, modlat, depths, mindx, meandx, d, S,&
  207. LOCRAD, RFACTOR2, nlobs_array, idm, jdm)
  208. deallocate(d, S, X5check)
  209. if (master) then
  210. print *, 'EnKF: finished calculation of X5, time = ', rtc() - time0
  211. end if
  212. allocate(fld(idm * jdm, ENSSIZE * NFIELD))
  213. #if defined(QMPI)
  214. call barrier()
  215. #endif
  216. ! get fieldnames and fieldlevels
  217. !
  218. call get_analysisfields()
  219. call distribute_iterations(numfields)
  220. #if defined(QMPI)
  221. call barrier() !KAL - just for "niceness" of output
  222. #endif
  223. time2 = rtc()
  224. do m1 = my_first_iteration, my_last_iteration, NFIELD
  225. m2 = min(my_last_iteration, m1 + NFIELD - 1)
  226. nfields = m2 - m1 + 1
  227. do m = m1, m2
  228. print '(a, i2, a, i3, a, a6, a, i3, a, f11.0)',&
  229. "I am ", qmpi_proc_num, ', m = ', m, ", field = ",&
  230. fieldnames(m), ", k = ", fieldlevel(m), ", time = ",&
  231. rtc() - time2
  232. do k = 1, ENSSIZE
  233. write(cmem, '(i3.3)') k
  234. memfile = 'forecast' // cmem
  235. !call get_mod_fld_new(trim(memfile), readfld, k, fieldnames(m),&
  236. ! fieldlevel(m), 1, idm, jdm)
  237. ! [CKB,FM]
  238. call io_mod_fld(readfld, k, enslist,fieldnames(m),fieldtype(m), &
  239. fieldlevel(m), 1, idm, jdm, 'get',FLOAT(obs(1)%date))
  240. ! end CKB,FM
  241. ! reshaping and conversion to real(4)
  242. fld(:, ENSSIZE * (m - m1) + k) = reshape(readfld, (/idm * jdm/))
  243. end do
  244. call p2nc_storeforecast(idm, jdm, ENSSIZE, numfields, m, fld(:, ENSSIZE * (m - m1) + 1 : ENSSIZE * (m + 1 - m1)))
  245. infls(m - m1 + 1) = prm_getinfl(trim(fieldnames(m)));
  246. end do
  247. call update_fields(idm, jdm, ENSSIZE, nfields, nlobs_array, depths,&
  248. fld(1,1), infls)
  249. do m = m1, m2
  250. fieldcounter = (m - my_first_iteration) + 1
  251. do k = 1, ENSSIZE
  252. write(cmem,'(i3.3)') k
  253. memfile = 'analysis' // cmem
  254. ! reshaping and conversion to real(8)
  255. readfld = reshape(fld(:, ENSSIZE * (m - m1) + k), (/idm, jdm/))
  256. write(text_string, '(a, i3.3)') '_proc', qmpi_proc_num
  257. !call put_mod_fld(trim(memfile) // trim(text_string), readfld, k,&
  258. ! fieldnames(m), fieldlevel(m), 1, fieldcounter, idm, jdm)
  259. ! [FM,CKB]
  260. call io_mod_fld(readfld, k, enslist, fieldnames(m), fieldtype(m), &
  261. fieldlevel(m), 1, idm, jdm, 'put',FLOAT(obs(1)%date))
  262. ! end FM,CKB
  263. end do
  264. end do
  265. end do
  266. deallocate(X5)
  267. deallocate(fld)
  268. call p2nc_writeforecast
  269. ! Barrier only necessary for timings
  270. #if defined(QMPI)
  271. call barrier()
  272. #endif
  273. if (master) then
  274. print *, 'EnKF: time for initialization = ', time1 - time0
  275. print *, 'EnKF: time for X5 calculation = ', time2 - time1
  276. print *, 'EnKF: time for ensemble update = ', rtc() - time2
  277. print *, 'EnKF: total time = ', rtc() - time0
  278. end if
  279. print *, 'EnKF: Finished'
  280. call stop_mpi()
  281. end program EnKF
  282. #if defined(_G95_)
  283. ! not tested! - PS
  284. !
  285. real function rtc()
  286. integer :: c
  287. call system_clock(count=c)
  288. rtc = dfloat(c)
  289. end function rtc
  290. #endif