dtadyn.F90 36 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717
  1. MODULE dtadyn
  2. !!======================================================================
  3. !! *** MODULE dtadyn ***
  4. !! Off-line : interpolation of the physical fields
  5. !!======================================================================
  6. !! History : OPA ! 1992-01 (M. Imbard) Original code
  7. !! 8.0 ! 1998-04 (L.Bopp MA Foujols) slopes for isopyc.
  8. !! - ! 1998-05 (L. Bopp) read output of coupled run
  9. !! 8.2 ! 2001-01 (M. Levy et M. Benjelloul) add netcdf FORMAT
  10. !! NEMO 1.0 ! 2005-03 (O. Aumont and A. El Moussaoui) F90
  11. !! - ! 2005-12 (C. Ethe) Adapted for DEGINT
  12. !! 3.0 ! 2007-06 (C. Ethe) use of iom module
  13. !! 3.3 ! 2010-11 (C. Ethe) Full reorganization of the off-line: phasing with the on-line
  14. !! 3.4 ! 2011-05 (C. Ethe) Use of fldread
  15. !!----------------------------------------------------------------------
  16. !!----------------------------------------------------------------------
  17. !! dta_dyn_init : initialization, namelist read, and SAVEs control
  18. !! dta_dyn : Interpolation of the fields
  19. !!----------------------------------------------------------------------
  20. USE oce ! ocean dynamics and tracers variables
  21. USE c1d ! 1D configuration: lk_c1d
  22. USE dom_oce ! ocean domain: variables
  23. USE domvvl ! variable volume
  24. USE zdf_oce ! ocean vertical physics: variables
  25. USE sbc_oce ! surface module: variables
  26. USE trc_oce ! share ocean/biogeo variables
  27. USE phycst ! physical constants
  28. USE trabbl ! active tracer: bottom boundary layer
  29. USE ldfslp ! lateral diffusion: iso-neutral slopes
  30. USE sbcrnf ! river runoffs
  31. USE ldfeiv ! eddy induced velocity coef.
  32. USE ldftra_oce ! ocean tracer lateral physics
  33. USE zdfmxl ! vertical physics: mixed layer depth
  34. USE eosbn2 ! equation of state - Brunt Vaisala frequency
  35. USE lbclnk ! ocean lateral boundary conditions (or mpp link)
  36. USE zpshde ! z-coord. with partial steps: horizontal derivatives
  37. USE in_out_manager ! I/O manager
  38. USE iom ! I/O library
  39. USE lib_mpp ! distributed memory computing library
  40. USE prtctl ! print control
  41. USE fldread ! read input fields
  42. USE wrk_nemo ! Memory allocation
  43. USE timing ! Timing
  44. USE trc, ONLY : ln_rsttr, numrtr, numrtw, lrst_trc
  45. IMPLICIT NONE
  46. PRIVATE
  47. PUBLIC dta_dyn_init ! called by opa.F90
  48. PUBLIC dta_dyn ! called by step.F90
  49. PUBLIC dta_dyn_swp ! called by step.F90
  50. CHARACTER(len=100) :: cn_dir !: Root directory for location of ssr files
  51. LOGICAL :: ln_ssh_ini !: initial ssh from dyn file (T) or not (F) - ssh is then read from passive tracer restart
  52. LOGICAL :: ln_dynrnf !: read runoff data in file (T) or set to zero (F)
  53. LOGICAL :: ln_dynrnf_depth !: read runoff data in file (T) or set to zero (F)
  54. REAL(wp) :: fwbcorr
  55. INTEGER , PARAMETER :: jpfld = 20 ! maximum number of fields to read
  56. INTEGER , SAVE :: jf_tem ! index of temperature
  57. INTEGER , SAVE :: jf_sal ! index of salinity
  58. INTEGER , SAVE :: jf_uwd ! index of u-transport
  59. INTEGER , SAVE :: jf_vwd ! index of v-transport
  60. INTEGER , SAVE :: jf_wwd ! index of v-transport
  61. INTEGER , SAVE :: jf_avt ! index of Kz
  62. INTEGER , SAVE :: jf_mld ! index of mixed layer deptht
  63. INTEGER , SAVE :: jf_emp ! index of water flux
  64. INTEGER , SAVE :: jf_empb ! index of water flux
  65. INTEGER , SAVE :: jf_qsr ! index of solar radiation
  66. INTEGER , SAVE :: jf_wnd ! index of wind speed
  67. INTEGER , SAVE :: jf_ice ! index of sea ice cover
  68. INTEGER , SAVE :: jf_rnf ! index of river runoff
  69. INTEGER , SAVE :: jf_fmf ! index of downward salt flux
  70. INTEGER , SAVE :: jf_ubl ! index of u-bbl coef
  71. INTEGER , SAVE :: jf_vbl ! index of v-bbl coef
  72. INTEGER , SAVE :: jf_div ! index of e3t
  73. TYPE(FLD), ALLOCATABLE, SAVE, DIMENSION(:) :: sf_dyn ! structure of input fields (file informations, fields read)
  74. ! !
  75. REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: uslpdta ! zonal isopycnal slopes
  76. REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: vslpdta ! meridional isopycnal slopes
  77. REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpidta ! zonal diapycnal slopes
  78. REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpjdta ! meridional diapycnal slopes
  79. INTEGER, SAVE :: nprevrec, nsecdyn
  80. !! * Substitutions
  81. # include "domzgr_substitute.h90"
  82. # include "vectopt_loop_substitute.h90"
  83. !!----------------------------------------------------------------------
  84. !! NEMO/OFF 3.3 , NEMO Consortium (2010)
  85. !! $Id: dtadyn.F90 4990 2014-12-15 16:42:49Z timgraham $
  86. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  87. !!----------------------------------------------------------------------
  88. CONTAINS
  89. SUBROUTINE dta_dyn( kt )
  90. !!----------------------------------------------------------------------
  91. !! *** ROUTINE dta_dyn ***
  92. !!
  93. !! ** Purpose : Prepares dynamics and physics fields from a NEMO run
  94. !! for an off-line simulation of passive tracers
  95. !!
  96. !! ** Method : calculates the position of data
  97. !! - computes slopes if needed
  98. !! - interpolates data if needed
  99. !!----------------------------------------------------------------------
  100. !
  101. USE oce, ONLY: zhdivtr => ua
  102. INTEGER, INTENT(in) :: kt ! ocean time-step index
  103. INTEGER :: ji, jj, jk
  104. REAL(wp), POINTER, DIMENSION(:,:) :: zemp
  105. !
  106. !!----------------------------------------------------------------------
  107. !
  108. IF( nn_timing == 1 ) CALL timing_start( 'dta_dyn')
  109. !
  110. !
  111. nsecdyn = nsec_year + nsec1jan000 ! number of seconds between Jan. 1st 00h of nit000 year and the middle of time step
  112. !
  113. IF( kt == nit000 ) THEN ; nprevrec = 0
  114. ELSE ; nprevrec = sf_dyn(jf_tem)%nrec_a(2)
  115. ENDIF
  116. !
  117. CALL fld_read( kt, 1, sf_dyn ) != read data at kt time step ==!
  118. !
  119. IF( lk_ldfslp .AND. .NOT.lk_c1d ) CALL dta_dyn_slp( kt ) ! Computation of slopes
  120. !
  121. tsn(:,:,:,jp_tem) = sf_dyn(jf_tem)%fnow(:,:,:) * tmask(:,:,:) ! temperature
  122. tsn(:,:,:,jp_sal) = sf_dyn(jf_sal)%fnow(:,:,:) * tmask(:,:,:) ! salinity
  123. wndm(:,:) = sf_dyn(jf_wnd)%fnow(:,:,1) * tmask(:,:,1) ! wind speed - needed for gas exchange
  124. fmmflx(:,:) = sf_dyn(jf_fmf)%fnow(:,:,1) * tmask(:,:,1) ! downward salt flux (v3.5+)
  125. fr_i(:,:) = sf_dyn(jf_ice)%fnow(:,:,1) * tmask(:,:,1) ! Sea-ice fraction
  126. qsr (:,:) = sf_dyn(jf_qsr)%fnow(:,:,1) * tmask(:,:,1) ! solar radiation
  127. emp (:,:) = sf_dyn(jf_emp)%fnow(:,:,1) * tmask(:,:,1) ! E-P
  128. IF( ln_dynrnf ) THEN
  129. rnf (:,:) = sf_dyn(jf_rnf)%fnow(:,:,1) * tmask(:,:,1) ! E-P
  130. IF( ln_dynrnf_depth .AND. lk_vvl ) CALL dta_dyn_hrnf
  131. ENDIF
  132. !
  133. un(:,:,:) = sf_dyn(jf_uwd)%fnow(:,:,:) * umask(:,:,:) ! effective u-transport
  134. vn(:,:,:) = sf_dyn(jf_vwd)%fnow(:,:,:) * vmask(:,:,:) ! effective v-transport
  135. wn(:,:,:) = sf_dyn(jf_wwd)%fnow(:,:,:) * tmask(:,:,:) ! effective v-transport
  136. !
  137. IF( lk_vvl ) THEN
  138. CALL wrk_alloc(jpi, jpj, zemp )
  139. zhdivtr(:,:,:) = sf_dyn(jf_div)%fnow(:,:,:) * tmask(:,:,:) ! effective u-transport
  140. emp_b (:,:) = sf_dyn(jf_empb)%fnow(:,:,1) * tmask(:,:,1) ! E-P
  141. zemp(:,:) = 0.5_wp * ( emp(:,:) + emp_b(:,:) ) + rnf(:,:) + fwbcorr * tmask(:,:,1)
  142. CALL dta_dyn_ssh( kt, zhdivtr, sshb, zemp, ssha, fse3t_a(:,:,:) ) != ssh, vertical scale factor & vertical transport
  143. CALL wrk_dealloc(jpi, jpj, zemp )
  144. ! Write in the tracer restart file
  145. ! *******************************
  146. IF( lrst_trc ) THEN
  147. IF(lwp) WRITE(numout,*)
  148. IF(lwp) WRITE(numout,*) 'dta_dyn_ssh : ssh field written in tracer restart file ', &
  149. & 'at it= ', kt,' date= ', ndastp
  150. IF(lwp) WRITE(numout,*) '~~~~'
  151. CALL iom_rstput( kt, nitrst, numrtw, 'sshn', ssha )
  152. CALL iom_rstput( kt, nitrst, numrtw, 'sshb', sshn )
  153. ENDIF
  154. ENDIF
  155. !
  156. CALL eos ( tsn, rhd, rhop, gdept_0(:,:,:) ) ! In any case, we need rhop
  157. CALL eos_rab( tsn, rab_n ) ! now local thermal/haline expension ratio at T-points
  158. CALL bn2 ( tsn, rab_n, rn2 ) ! before Brunt-Vaisala frequency need for zdfmxl
  159. rn2b(:,:,:) = rn2(:,:,:) ! need for zdfmxl
  160. CALL zdf_mxl( kt ) ! In any case, we need mxl
  161. !
  162. hmld(:,:) = sf_dyn(jf_mld)%fnow(:,:,1) * tmask(:,:,1) ! mixed layer depht
  163. avt(:,:,:) = sf_dyn(jf_avt)%fnow(:,:,:) * tmask(:,:,:) ! vertical diffusive coefficient
  164. !
  165. #if defined key_trabbl && ! defined key_c1d
  166. ahu_bbl(:,:) = sf_dyn(jf_ubl)%fnow(:,:,1) * umask(:,:,1) ! bbl diffusive coef
  167. ahv_bbl(:,:) = sf_dyn(jf_vbl)%fnow(:,:,1) * vmask(:,:,1)
  168. #endif
  169. !
  170. !
  171. CALL eos( tsn, rhd, rhop, gdept_0(:,:,:) ) ! In any case, we need rhop
  172. !
  173. IF(ln_ctl) THEN ! print control
  174. CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' tn - : ', mask1=tmask, ovlap=1, kdim=jpk )
  175. CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_sal), clinfo1=' sn - : ', mask1=tmask, ovlap=1, kdim=jpk )
  176. CALL prt_ctl(tab3d_1=un , clinfo1=' un - : ', mask1=umask, ovlap=1, kdim=jpk )
  177. CALL prt_ctl(tab3d_1=vn , clinfo1=' vn - : ', mask1=vmask, ovlap=1, kdim=jpk )
  178. CALL prt_ctl(tab3d_1=wn , clinfo1=' wn - : ', mask1=tmask, ovlap=1, kdim=jpk )
  179. CALL prt_ctl(tab3d_1=avt , clinfo1=' kz - : ', mask1=tmask, ovlap=1, kdim=jpk )
  180. ! CALL prt_ctl(tab2d_1=fr_i , clinfo1=' fr_i - : ', mask1=tmask, ovlap=1 )
  181. ! CALL prt_ctl(tab2d_1=hmld , clinfo1=' hmld - : ', mask1=tmask, ovlap=1 )
  182. ! CALL prt_ctl(tab2d_1=fmmflx , clinfo1=' fmmflx - : ', mask1=tmask, ovlap=1 )
  183. ! CALL prt_ctl(tab2d_1=emp , clinfo1=' emp - : ', mask1=tmask, ovlap=1 )
  184. ! CALL prt_ctl(tab2d_1=wndm , clinfo1=' wspd - : ', mask1=tmask, ovlap=1 )
  185. ! CALL prt_ctl(tab2d_1=qsr , clinfo1=' qsr - : ', mask1=tmask, ovlap=1 )
  186. ENDIF
  187. !
  188. IF( nn_timing == 1 ) CALL timing_stop( 'dta_dyn')
  189. !
  190. END SUBROUTINE dta_dyn
  191. SUBROUTINE dta_dyn_init
  192. !!----------------------------------------------------------------------
  193. !! *** ROUTINE dta_dyn_init ***
  194. !!
  195. !! ** Purpose : Initialisation of the dynamical data
  196. !! ** Method : - read the data namdta_dyn namelist
  197. !!
  198. !! ** Action : - read parameters
  199. !!----------------------------------------------------------------------
  200. INTEGER :: ierr, ierr0, ierr1, ierr2, ierr3 ! return error code
  201. INTEGER :: ifpr ! dummy loop indice
  202. INTEGER :: jfld ! dummy loop arguments
  203. INTEGER :: inum, idv, idimv ! local integer
  204. INTEGER :: ios ! Local integer output status for namelist read
  205. INTEGER :: ji, jj, jk
  206. REAL(wp) :: zcoef
  207. INTEGER :: nkrnf_max
  208. REAL(wp) :: hrnf_max
  209. !!
  210. CHARACTER(len=100) :: cn_dir ! Root directory for location of core files
  211. TYPE(FLD_N), DIMENSION(jpfld) :: slf_d ! array of namelist informations on the fields to read
  212. TYPE(FLD_N) :: sn_uwd, sn_vwd, sn_wwd, sn_empb, sn_emp ! informations about the fields to be read
  213. TYPE(FLD_N) :: sn_tem , sn_sal , sn_avt ! " "
  214. TYPE(FLD_N) :: sn_mld, sn_qsr, sn_wnd , sn_ice , sn_fmf ! " "
  215. TYPE(FLD_N) :: sn_ubl, sn_vbl, sn_rnf ! " "
  216. TYPE(FLD_N) :: sn_div ! informations about the fields to be read
  217. !!----------------------------------------------------------------------
  218. NAMELIST/namdta_dyn/cn_dir, ln_dynrnf, ln_dynrnf_depth, ln_ssh_ini, fwbcorr, &
  219. & sn_uwd, sn_vwd, sn_wwd, sn_emp, &
  220. & sn_avt, sn_tem, sn_sal, sn_mld , sn_qsr , &
  221. & sn_wnd, sn_ice, sn_fmf, &
  222. & sn_ubl, sn_vbl, sn_rnf, &
  223. & sn_empb, sn_div
  224. !
  225. REWIND( numnam_ref ) ! Namelist namdta_dyn in reference namelist : Offline: init. of dynamical data
  226. READ ( numnam_ref, namdta_dyn, IOSTAT = ios, ERR = 901)
  227. 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdta_dyn in reference namelist', lwp )
  228. REWIND( numnam_cfg ) ! Namelist namdta_dyn in configuration namelist : Offline: init. of dynamical data
  229. READ ( numnam_cfg, namdta_dyn, IOSTAT = ios, ERR = 902 )
  230. 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdta_dyn in configuration namelist', lwp )
  231. IF(lwm) WRITE ( numond, namdta_dyn )
  232. ! ! store namelist information in an array
  233. ! ! Control print
  234. IF(lwp) THEN
  235. WRITE(numout,*)
  236. WRITE(numout,*) 'dta_dyn : offline dynamics '
  237. WRITE(numout,*) '~~~~~~~ '
  238. WRITE(numout,*) ' Namelist namdta_dyn'
  239. WRITE(numout,*) ' ssh initialised from dyn file (T) or not (F) ln_ssh_ini = ', ln_ssh_ini
  240. WRITE(numout,*) ' runoffs option enabled (T) or not (F) ln_dynrnf = ', ln_dynrnf
  241. WRITE(numout,*) ' runoffs is spread in vertical ln_dynrnf_depth = ', ln_dynrnf_depth
  242. WRITE(numout,*) ' annual global mean of empmr for ssh correction fwbcorr = ', fwbcorr
  243. WRITE(numout,*)
  244. ENDIF
  245. !
  246. jf_uwd = 1 ; jf_vwd = 2 ; jf_wwd = 3 ; jf_emp = 4 ; jf_avt = 5
  247. jf_tem = 6 ; jf_sal = 7 ; jf_mld = 8 ; jf_qsr = 9
  248. jf_wnd = 10 ; jf_ice = 11 ; jf_fmf = 12 ; jfld = jf_fmf
  249. !
  250. slf_d(jf_uwd) = sn_uwd ; slf_d(jf_vwd) = sn_vwd ; slf_d(jf_wwd) = sn_wwd
  251. slf_d(jf_emp) = sn_emp ; slf_d(jf_avt) = sn_avt
  252. slf_d(jf_tem) = sn_tem ; slf_d(jf_sal) = sn_sal ; slf_d(jf_mld) = sn_mld
  253. slf_d(jf_qsr) = sn_qsr ; slf_d(jf_wnd) = sn_wnd ; slf_d(jf_ice) = sn_ice
  254. slf_d(jf_fmf) = sn_fmf
  255. !
  256. IF( lk_vvl ) THEN
  257. jf_div = jfld + 1 ; jf_empb = jfld + 2 ; jfld = jf_empb
  258. slf_d(jf_div) = sn_div ; slf_d(jf_empb) = sn_empb
  259. ENDIF
  260. !
  261. IF( lk_trabbl ) THEN
  262. jf_ubl = jfld + 1 ; jf_vbl = jfld + 2 ; jfld = jf_vbl
  263. slf_d(jf_ubl) = sn_ubl ; slf_d(jf_vbl) = sn_vbl
  264. ENDIF
  265. !
  266. IF( ln_dynrnf ) THEN
  267. jf_rnf = jfld + 1 ; jfld = jf_rnf
  268. slf_d(jf_rnf) = sn_rnf
  269. ELSE
  270. rnf(:,:) = 0._wp
  271. ENDIF
  272. ALLOCATE( sf_dyn(jfld), STAT=ierr ) ! set sf structure
  273. IF( ierr > 0 ) THEN
  274. CALL ctl_stop( 'dta_dyn: unable to allocate sf structure' ) ; RETURN
  275. ENDIF
  276. ! ! fill sf with slf_i and control print
  277. CALL fld_fill( sf_dyn, slf_d, cn_dir, 'dta_dyn_init', 'Data in file', 'namdta_dyn' )
  278. !
  279. ! Open file for each variable to get his number of dimension
  280. DO ifpr = 1, jfld
  281. CALL fld_clopn( sf_dyn(ifpr), nyear, nmonth, nday )
  282. idv = iom_varid( sf_dyn(ifpr)%num , slf_d(ifpr)%clvar ) ! id of the variable sdjf%clvar
  283. idimv = iom_file ( sf_dyn(ifpr)%num )%ndims(idv) ! number of dimension for variable sdjf%clvar
  284. IF( sf_dyn(ifpr)%num /= 0 ) CALL iom_close( sf_dyn(ifpr)%num ) ! close file if already open
  285. ierr1=0
  286. IF( idimv == 3 ) THEN ! 2D variable
  287. ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,1) , STAT=ierr0 )
  288. IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,1,2) , STAT=ierr1 )
  289. ELSE ! 3D variable
  290. ALLOCATE( sf_dyn(ifpr)%fnow(jpi,jpj,jpk) , STAT=ierr0 )
  291. IF( slf_d(ifpr)%ln_tint ) ALLOCATE( sf_dyn(ifpr)%fdta(jpi,jpj,jpk,2), STAT=ierr1 )
  292. ENDIF
  293. IF( ierr0 + ierr1 > 0 ) THEN
  294. CALL ctl_stop( 'dta_dyn_init : unable to allocate sf_dyn array structure' ) ; RETURN
  295. ENDIF
  296. END DO
  297. !
  298. IF( lk_ldfslp .AND. .NOT.lk_c1d ) THEN ! slopes
  299. IF( sf_dyn(jf_tem)%ln_tint ) THEN ! time interpolation
  300. ALLOCATE( uslpdta (jpi,jpj,jpk,2), vslpdta (jpi,jpj,jpk,2), &
  301. & wslpidta(jpi,jpj,jpk,2), wslpjdta(jpi,jpj,jpk,2), STAT=ierr2 )
  302. !
  303. IF( ierr2 > 0 ) THEN
  304. CALL ctl_stop( 'dta_dyn_init : unable to allocate slope arrays' ) ; RETURN
  305. ENDIF
  306. ENDIF
  307. ENDIF
  308. !
  309. IF( lk_vvl ) THEN
  310. IF( ln_ssh_ini ) THEN ! Restart: read in restart file
  311. IF(lwp) WRITE(numout,*) ' sshn forcing fields read in the dynamics restart file for initialisation'
  312. CALL iom_open( 'restart', inum )
  313. CALL iom_get( inum, jpdom_autoglo, 'sshn', sshn(:,:) )
  314. CALL iom_get( inum, jpdom_autoglo, 'sshb', sshb(:,:) )
  315. CALL iom_close( inum ) ! close file
  316. ELSE
  317. IF(lwp) WRITE(numout,*) ' sshn forcing fields read in passive tracers restart file for initialisation'
  318. CALL iom_get( numrtr, jpdom_autoglo, 'sshn', sshn(:,:) )
  319. CALL iom_get( numrtr, jpdom_autoglo, 'sshb', sshb(:,:) )
  320. ENDIF
  321. !
  322. DO jk = 1, jpkm1
  323. fse3t_n(:,:,jk) = e3t_0(:,:,jk) * ( 1._wp + sshn(:,:) * tmask(:,:,1) / ( ht_0(:,:) + 1.0 - tmask(:,:,1) ) )
  324. ENDDO
  325. fse3t_a(:,:,jpk) = e3t_0(:,:,jpk)
  326. ! Horizontal scale factor interpolations
  327. ! --------------------------------------
  328. CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3u_n(:,:,:), 'U' )
  329. CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3v_n(:,:,:), 'V' )
  330. ! Vertical scale factor interpolations
  331. ! ------------------------------------
  332. CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n(:,:,:), 'W' )
  333. fse3t_b(:,:,:) = fse3t_n(:,:,:)
  334. fse3u_b(:,:,:) = fse3u_n(:,:,:)
  335. fse3v_b(:,:,:) = fse3v_n(:,:,:)
  336. ! t- and w- points depth
  337. ! ----------------------
  338. fsdept_n(:,:,1) = 0.5_wp * fse3w_n(:,:,1)
  339. fsdepw_n(:,:,1) = 0.0_wp
  340. DO jk = 2, jpk
  341. DO jj = 1,jpj
  342. DO ji = 1,jpi
  343. ! zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk)) ! 0 everywhere
  344. ! tmask = wmask, ie everywhere expect at jk = mikt
  345. ! 1 for jk =
  346. ! mikt
  347. zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
  348. fsdepw_n(ji,jj,jk) = fsdepw_n(ji,jj,jk-1) + fse3t_n(ji,jj,jk-1)
  349. fsdept_n(ji,jj,jk) = zcoef * ( fsdepw_n(ji,jj,jk ) + 0.5 * fse3w_n(ji,jj,jk)) &
  350. & + (1-zcoef) * ( fsdept_n(ji,jj,jk-1) + fse3w_n(ji,jj,jk))
  351. END DO
  352. END DO
  353. END DO
  354. fsdept_b(:,:,:) = fsdept_n(:,:,:)
  355. fsdepw_b(:,:,:) = fsdepw_n(:,:,:)
  356. !
  357. ENDIF
  358. !
  359. IF( ln_dynrnf .AND. ln_dynrnf_depth ) THEN ! read depht over which runoffs are distributed
  360. IF(lwp) WRITE(numout,*)
  361. IF(lwp) WRITE(numout,*) ' read in the file depht over which runoffs are distributed'
  362. CALL iom_open ( "runoffs", inum ) ! open file
  363. CALL iom_get ( inum, jpdom_data, 'rodepth', h_rnf ) ! read the river mouth array
  364. CALL iom_close( inum ) ! close file
  365. !
  366. nk_rnf(:,:) = 0 ! set the number of level over which river runoffs are applied
  367. DO jj = 1, jpj
  368. DO ji = 1, jpi
  369. IF( h_rnf(ji,jj) > 0._wp ) THEN
  370. jk = 2
  371. DO WHILE ( jk /= mbkt(ji,jj) .AND. gdept_0(ji,jj,jk) < h_rnf(ji,jj) ) ; jk = jk + 1
  372. END DO
  373. nk_rnf(ji,jj) = jk
  374. ELSEIF( h_rnf(ji,jj) == -1._wp ) THEN ; nk_rnf(ji,jj) = 1
  375. ELSEIF( h_rnf(ji,jj) == -999._wp ) THEN ; nk_rnf(ji,jj) = mbkt(ji,jj)
  376. ELSE
  377. CALL ctl_stop( 'sbc_rnf_init: runoff depth not positive, and not -999 or -1, rnf value in file fort.999' )
  378. WRITE(999,*) 'ji, jj, h_rnf(ji,jj) :', ji, jj, h_rnf(ji,jj)
  379. ENDIF
  380. END DO
  381. END DO
  382. DO jj = 1, jpj ! set the associated depth
  383. DO ji = 1, jpi
  384. h_rnf(ji,jj) = 0._wp
  385. DO jk = 1, nk_rnf(ji,jj)
  386. h_rnf(ji,jj) = h_rnf(ji,jj) + fse3t(ji,jj,jk)
  387. END DO
  388. END DO
  389. END DO
  390. ELSE ! runoffs applied at the surface
  391. nk_rnf(:,:) = 1
  392. h_rnf (:,:) = fse3t(:,:,1)
  393. ENDIF
  394. nkrnf_max = MAXVAL( nk_rnf(:,:) )
  395. hrnf_max = MAXVAL( h_rnf(:,:) )
  396. IF( lk_mpp ) THEN
  397. CALL mpp_max( nkrnf_max ) ! max over the global domain
  398. CALL mpp_max( hrnf_max ) ! max over the global domain
  399. ENDIF
  400. IF(lwp) WRITE(numout,*) ' '
  401. IF(lwp) WRITE(numout,*) ' max depht of runoff : ', hrnf_max,' max level : ', nkrnf_max
  402. IF(lwp) WRITE(numout,*) ' '
  403. !
  404. CALL dta_dyn( nit000 )
  405. !
  406. END SUBROUTINE dta_dyn_init
  407. SUBROUTINE dta_dyn_swp( kt )
  408. !!---------------------------------------------------------------------
  409. !! *** ROUTINE dta_dyn_swp ***
  410. !!
  411. !! ** Purpose : Swap and the data and compute the vertical scale factor at U/V/W point
  412. !! and the depht
  413. !!
  414. !!---------------------------------------------------------------------
  415. INTEGER, INTENT(in) :: kt ! time step
  416. INTEGER :: ji, jj, jk
  417. REAL(wp) :: zcoef
  418. !
  419. !!---------------------------------------------------------------------
  420. IF( kt == nit000 ) THEN
  421. IF(lwp) WRITE(numout,*)
  422. IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
  423. IF(lwp) WRITE(numout,*) '~~~~~~~ '
  424. ENDIF
  425. sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:)) ! before <-- now filtered
  426. sshn(:,:) = ssha(:,:)
  427. fse3t_n(:,:,:) = fse3t_a(:,:,:)
  428. ! Reconstruction of all vertical scale factors at now and before time steps
  429. ! =============================================================================
  430. ! Horizontal scale factor interpolations
  431. ! --------------------------------------
  432. CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3u_n(:,:,:), 'U' )
  433. CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3v_n(:,:,:), 'V' )
  434. ! Vertical scale factor interpolations
  435. ! ------------------------------------
  436. CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n (:,:,:), 'W' )
  437. fse3t_b(:,:,:) = fse3t_n(:,:,:)
  438. fse3u_b(:,:,:) = fse3u_n(:,:,:)
  439. fse3v_b(:,:,:) = fse3v_n(:,:,:)
  440. ! t- and w- points depth
  441. ! ----------------------
  442. fsdept_n(:,:,1) = 0.5_wp * fse3w_n(:,:,1)
  443. fsdepw_n(:,:,1) = 0.0_wp
  444. DO jk = 2, jpk
  445. DO jj = 1,jpj
  446. DO ji = 1,jpi
  447. zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
  448. fsdepw_n(ji,jj,jk) = fsdepw_n(ji,jj,jk-1) + fse3t_n(ji,jj,jk-1)
  449. fsdept_n(ji,jj,jk) = zcoef * ( fsdepw_n(ji,jj,jk ) + 0.5 * fse3w_n(ji,jj,jk)) &
  450. & + (1-zcoef) * ( fsdept_n(ji,jj,jk-1) + fse3w_n(ji,jj,jk))
  451. END DO
  452. END DO
  453. END DO
  454. fsdept_b(:,:,:) = fsdept_n(:,:,:)
  455. fsdepw_b(:,:,:) = fsdepw_n(:,:,:)
  456. !
  457. END SUBROUTINE dta_dyn_swp
  458. SUBROUTINE dta_dyn_ssh( kt, phdivtr, psshb, pemp, pssha, pe3ta )
  459. !!----------------------------------------------------------------------
  460. !! *** ROUTINE dta_dyn_wzv ***
  461. !!
  462. !! ** Purpose : compute the after ssh (ssha) and the now vertical velocity
  463. !!
  464. !! ** Method : Using the incompressibility hypothesis,
  465. !! - the ssh increment is computed by integrating the horizontal divergence
  466. !! and multiply by the time step.
  467. !!
  468. !! - compute the after scale factor : repartition of ssh INCREMENT proportionnaly
  469. !! to the level thickness ( z-star case )
  470. !!
  471. !! - the vertical velocity is computed by integrating the horizontal divergence
  472. !! from the bottom to the surface minus the scale factor evolution.
  473. !! The boundary conditions are w=0 at the bottom (no flux)
  474. !!
  475. !! ** action : ssha / e3t_a / wn
  476. !!
  477. !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling.
  478. !!----------------------------------------------------------------------
  479. !! * Arguments
  480. INTEGER, INTENT(in ) :: kt ! time-step
  481. REAL(wp), DIMENSION(jpi,jpj,jpk) , INTENT(in ) :: phdivtr ! horizontal divergence transport
  482. REAL(wp), DIMENSION(jpi,jpj) , OPTIONAL, INTENT(in ) :: psshb ! now ssh
  483. REAL(wp), DIMENSION(jpi,jpj) , OPTIONAL, INTENT(in ) :: pemp ! evaporation minus precipitation
  484. REAL(wp), DIMENSION(jpi,jpj) , OPTIONAL, INTENT(inout) :: pssha ! after ssh
  485. REAL(wp), DIMENSION(jpi,jpj,jpk), OPTIONAL, INTENT(out) :: pe3ta ! after vertical scale factor
  486. !! * Local declarations
  487. INTEGER :: jk
  488. REAL(wp), DIMENSION(jpi,jpj) :: zhdiv
  489. REAL(wp) :: z2dt
  490. !!----------------------------------------------------------------------
  491. !
  492. z2dt = 2._wp * rdt
  493. !
  494. zhdiv(:,:) = 0._wp
  495. DO jk = 1, jpkm1
  496. zhdiv(:,:) = zhdiv(:,:) + phdivtr(:,:,jk) * tmask(:,:,jk)
  497. END DO
  498. ! ! Sea surface elevation time-stepping
  499. pssha(:,:) = ( psshb(:,:) - z2dt * ( r1_rau0 * pemp(:,:) + zhdiv(:,:) ) ) * ssmask(:,:)
  500. ! !
  501. ! ! After acale factors at t-points ( z_star coordinate )
  502. DO jk = 1, jpkm1
  503. pe3ta(:,:,jk) = e3t_0(:,:,jk) * ( 1._wp + pssha(:,:) * tmask(:,:,1) / ( ht_0(:,:) + 1.0 - tmask(:,:,1) ) )
  504. END DO
  505. !
  506. END SUBROUTINE dta_dyn_ssh
  507. SUBROUTINE dta_dyn_hrnf
  508. !!----------------------------------------------------------------------
  509. !! *** ROUTINE sbc_rnf ***
  510. !!
  511. !! ** Purpose : update the horizontal divergence with the runoff inflow
  512. !!
  513. !! ** Method :
  514. !! CAUTION : rnf is positive (inflow) decreasing the
  515. !! divergence and expressed in m/s
  516. !!
  517. !! ** Action : phdivn decreased by the runoff inflow
  518. !!----------------------------------------------------------------------
  519. !!
  520. INTEGER :: ji, jj, jk ! dummy loop indices
  521. !!----------------------------------------------------------------------
  522. !
  523. DO jj = 1, jpj ! update the depth over which runoffs are distributed
  524. DO ji = 1, jpi
  525. h_rnf(ji,jj) = 0._wp
  526. DO jk = 1, nk_rnf(ji,jj) ! recalculates h_rnf to be the depth in metres
  527. h_rnf(ji,jj) = h_rnf(ji,jj) + fse3t(ji,jj,jk) ! to the bottom of the relevant grid box
  528. END DO
  529. END DO
  530. END DO
  531. !
  532. END SUBROUTINE dta_dyn_hrnf
  533. SUBROUTINE dta_dyn_slp( kt )
  534. !!---------------------------------------------------------------------
  535. !! *** ROUTINE dta_dyn_slp ***
  536. !!
  537. !! ** Purpose : Computation of slope
  538. !!
  539. !!---------------------------------------------------------------------
  540. USE oce, ONLY: zts => tsa
  541. !
  542. INTEGER, INTENT(in) :: kt ! time step
  543. !
  544. INTEGER :: ji, jj ! dummy loop indices
  545. REAL(wp) :: ztinta ! ratio applied to after records when doing time interpolation
  546. REAL(wp) :: ztintb ! ratio applied to before records when doing time interpolation
  547. INTEGER :: iswap
  548. REAL(wp), POINTER, DIMENSION(:,:,:) :: zuslp, zvslp, zwslpi, zwslpj
  549. !!---------------------------------------------------------------------
  550. !
  551. CALL wrk_alloc(jpi, jpj, jpk, zuslp, zvslp, zwslpi, zwslpj )
  552. !
  553. IF( sf_dyn(jf_tem)%ln_tint ) THEN ! Computes slopes (here avt is used as workspace)
  554. IF( kt == nit000 ) THEN
  555. zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,1) * tmask(:,:,:) ! temperature
  556. zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,1) * tmask(:,:,:) ! salinity
  557. avt(:,:,:) = sf_dyn(jf_avt)%fdta(:,:,:,1) * tmask(:,:,:) ! vertical diffusive coef.
  558. CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj )
  559. uslpdta (:,:,:,1) = zuslp (:,:,:)
  560. vslpdta (:,:,:,1) = zvslp (:,:,:)
  561. wslpidta(:,:,:,1) = zwslpi(:,:,:)
  562. wslpjdta(:,:,:,1) = zwslpj(:,:,:)
  563. !
  564. zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,2) * tmask(:,:,:) ! temperature
  565. zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,2) * tmask(:,:,:) ! salinity
  566. avt(:,:,:) = sf_dyn(jf_avt)%fdta(:,:,:,2) * tmask(:,:,:) ! vertical diffusive coef.
  567. CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj )
  568. uslpdta (:,:,:,2) = zuslp (:,:,:)
  569. vslpdta (:,:,:,2) = zvslp (:,:,:)
  570. wslpidta(:,:,:,2) = zwslpi(:,:,:)
  571. wslpjdta(:,:,:,2) = zwslpj(:,:,:)
  572. ELSE
  573. !
  574. iswap = 0
  575. IF( sf_dyn(jf_tem)%nrec_a(2) - nprevrec /= 0 ) iswap = 1
  576. IF( nsecdyn > sf_dyn(jf_tem)%nrec_b(2) .AND. iswap == 1 ) THEN ! read/update the after data
  577. IF(lwp) WRITE(numout,*) ' Compute new slopes at kt = ', kt
  578. uslpdta (:,:,:,1) = uslpdta (:,:,:,2) ! swap the data
  579. vslpdta (:,:,:,1) = vslpdta (:,:,:,2)
  580. wslpidta(:,:,:,1) = wslpidta(:,:,:,2)
  581. wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2)
  582. !
  583. zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fdta(:,:,:,2) * tmask(:,:,:) ! temperature
  584. zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fdta(:,:,:,2) * tmask(:,:,:) ! salinity
  585. avt(:,:,:) = sf_dyn(jf_avt)%fdta(:,:,:,2) * tmask(:,:,:) ! vertical diffusive coef.
  586. CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj )
  587. !
  588. uslpdta (:,:,:,2) = zuslp (:,:,:)
  589. vslpdta (:,:,:,2) = zvslp (:,:,:)
  590. wslpidta(:,:,:,2) = zwslpi(:,:,:)
  591. wslpjdta(:,:,:,2) = zwslpj(:,:,:)
  592. ENDIF
  593. ENDIF
  594. ENDIF
  595. !
  596. IF( sf_dyn(jf_tem)%ln_tint ) THEN
  597. ztinta = REAL( nsecdyn - sf_dyn(jf_tem)%nrec_b(2), wp ) &
  598. & / REAL( sf_dyn(jf_tem)%nrec_a(2) - sf_dyn(jf_tem)%nrec_b(2), wp )
  599. ztintb = 1. - ztinta
  600. #if defined key_ldfslp && ! defined key_c1d
  601. uslp (:,:,:) = ztintb * uslpdta (:,:,:,1) + ztinta * uslpdta (:,:,:,2)
  602. vslp (:,:,:) = ztintb * vslpdta (:,:,:,1) + ztinta * vslpdta (:,:,:,2)
  603. wslpi(:,:,:) = ztintb * wslpidta(:,:,:,1) + ztinta * wslpidta(:,:,:,2)
  604. wslpj(:,:,:) = ztintb * wslpjdta(:,:,:,1) + ztinta * wslpjdta(:,:,:,2)
  605. #endif
  606. ELSE
  607. zts(:,:,:,jp_tem) = sf_dyn(jf_tem)%fnow(:,:,:) * tmask(:,:,:) ! temperature
  608. zts(:,:,:,jp_sal) = sf_dyn(jf_sal)%fnow(:,:,:) * tmask(:,:,:) ! salinity
  609. avt(:,:,:) = sf_dyn(jf_avt)%fnow(:,:,:) * tmask(:,:,:) ! vertical diffusive coef.
  610. CALL compute_slopes( kt, zts, zuslp, zvslp, zwslpi, zwslpj )
  611. !
  612. #if defined key_ldfslp && ! defined key_c1d
  613. uslp (:,:,:) = zuslp (:,:,:)
  614. vslp (:,:,:) = zvslp (:,:,:)
  615. wslpi(:,:,:) = zwslpi(:,:,:)
  616. wslpj(:,:,:) = zwslpj(:,:,:)
  617. #endif
  618. ENDIF
  619. !
  620. CALL wrk_dealloc(jpi, jpj, jpk, zuslp, zvslp, zwslpi, zwslpj )
  621. !
  622. END SUBROUTINE dta_dyn_slp
  623. SUBROUTINE compute_slopes( kt, pts, puslp, pvslp, pwslpi, pwslpj )
  624. !!---------------------------------------------------------------------
  625. !! *** ROUTINE dta_dyn_slp ***
  626. !!
  627. !! ** Purpose : Computation of slope
  628. !!
  629. !!---------------------------------------------------------------------
  630. INTEGER , INTENT(in ) :: kt ! time step
  631. REAL(wp), DIMENSION(jpi,jpj,jpk,jpts), INTENT(in ) :: pts ! temperature/salinity
  632. REAL(wp), DIMENSION(jpi,jpj,jpk) , INTENT(out) :: puslp ! zonal isopycnal slopes
  633. REAL(wp), DIMENSION(jpi,jpj,jpk) , INTENT(out) :: pvslp ! meridional isopycnal slopes
  634. REAL(wp), DIMENSION(jpi,jpj,jpk) , INTENT(out) :: pwslpi ! zonal diapycnal slopes
  635. REAL(wp), DIMENSION(jpi,jpj,jpk) , INTENT(out) :: pwslpj ! meridional diapycnal slopes
  636. !!---------------------------------------------------------------------
  637. #if defined key_ldfslp && ! defined key_c1d
  638. CALL eos ( pts, rhd, rhop, gdept_0(:,:,:) )
  639. CALL eos_rab( pts, rab_n ) ! now local thermal/haline expension ratio at T-points
  640. CALL bn2 ( pts, rab_n, rn2 ) ! now Brunt-Vaisala
  641. ! Partial steps: before Horizontal DErivative
  642. IF( ln_zps .AND. .NOT. ln_isfcav) &
  643. & CALL zps_hde ( kt, jpts, pts, gtsu, gtsv, & ! Partial steps: before horizontal gradient
  644. & rhd, gru , grv ) ! of t, s, rd at the last ocean level
  645. IF( ln_zps .AND. ln_isfcav) &
  646. & CALL zps_hde_isf( kt, jpts, pts, gtsu, gtsv, & ! Partial steps for top cell (ISF)
  647. & rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv , &
  648. & gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi ) ! of t, s, rd at the first ocean level
  649. rn2b(:,:,:) = rn2(:,:,:) ! need for zdfmxl
  650. CALL zdf_mxl( kt ) ! mixed layer depth
  651. CALL ldf_slp( kt, rhd, rn2 ) ! slopes
  652. puslp (:,:,:) = uslp (:,:,:)
  653. pvslp (:,:,:) = vslp (:,:,:)
  654. pwslpi(:,:,:) = wslpi(:,:,:)
  655. pwslpj(:,:,:) = wslpj(:,:,:)
  656. #else
  657. puslp (:,:,:) = 0. ! to avoid warning when compiling
  658. pvslp (:,:,:) = 0.
  659. pwslpi(:,:,:) = 0.
  660. pwslpj(:,:,:) = 0.
  661. #endif
  662. !
  663. END SUBROUTINE compute_slopes
  664. !!======================================================================
  665. END MODULE dtadyn