EnKF.f90 10 KB

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