sshwzv.F90 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282
  1. MODULE sshwzv
  2. !!==============================================================================
  3. !! *** MODULE sshwzv ***
  4. !! Ocean dynamics : sea surface height and vertical velocity
  5. !!==============================================================================
  6. !! History : 3.1 ! 2009-02 (G. Madec, M. Leclair) Original code
  7. !! 3.3 ! 2010-04 (M. Leclair, G. Madec) modified LF-RA
  8. !! - ! 2010-05 (K. Mogensen, A. Weaver, M. Martin, D. Lea) Assimilation interface
  9. !! - ! 2010-09 (D.Storkey and E.O'Dea) bug fixes for BDY module
  10. !! 3.3 ! 2011-10 (M. Leclair) split former ssh_wzv routine and remove all vvl related work
  11. !!----------------------------------------------------------------------
  12. !!----------------------------------------------------------------------
  13. !! ssh_nxt : after ssh
  14. !! ssh_swp : filter ans swap the ssh arrays
  15. !! wzv : compute now vertical velocity
  16. !!----------------------------------------------------------------------
  17. USE oce ! ocean dynamics and tracers variables
  18. USE dom_oce ! ocean space and time domain variables
  19. USE sbc_oce ! surface boundary condition: ocean
  20. USE domvvl ! Variable volume
  21. USE divcur ! hor. divergence and curl (div & cur routines)
  22. USE restart ! only for lrst_oce
  23. USE in_out_manager ! I/O manager
  24. USE prtctl ! Print control
  25. USE phycst
  26. USE lbclnk ! ocean lateral boundary condition (or mpp link)
  27. USE lib_mpp ! MPP library
  28. USE bdy_oce
  29. USE bdy_par
  30. USE bdydyn2d ! bdy_ssh routine
  31. #if defined key_agrif
  32. USE agrif_opa_interp
  33. #endif
  34. #if defined key_asminc
  35. USE asminc ! Assimilation increment
  36. #endif
  37. USE wrk_nemo ! Memory Allocation
  38. USE timing ! Timing
  39. IMPLICIT NONE
  40. PRIVATE
  41. PUBLIC ssh_nxt ! called by step.F90
  42. PUBLIC wzv ! called by step.F90
  43. PUBLIC ssh_swp ! called by step.F90
  44. !! * Substitutions
  45. # include "domzgr_substitute.h90"
  46. # include "vectopt_loop_substitute.h90"
  47. !!----------------------------------------------------------------------
  48. !! NEMO/OPA 3.3 , NEMO Consortium (2010)
  49. !! $Id: sshwzv.F90 5628 2015-07-22 20:26:35Z mathiot $
  50. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  51. !!----------------------------------------------------------------------
  52. CONTAINS
  53. SUBROUTINE ssh_nxt( kt )
  54. !!----------------------------------------------------------------------
  55. !! *** ROUTINE ssh_nxt ***
  56. !!
  57. !! ** Purpose : compute the after ssh (ssha)
  58. !!
  59. !! ** Method : - Using the incompressibility hypothesis, the ssh increment
  60. !! is computed by integrating the horizontal divergence and multiply by
  61. !! by the time step.
  62. !!
  63. !! ** action : ssha : after sea surface height
  64. !!
  65. !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling.
  66. !!----------------------------------------------------------------------
  67. !
  68. REAL(wp), POINTER, DIMENSION(:,: ) :: zhdiv
  69. INTEGER, INTENT(in) :: kt ! time step
  70. !
  71. INTEGER :: jk ! dummy loop indice
  72. REAL(wp) :: z2dt, z1_rau0 ! local scalars
  73. !!----------------------------------------------------------------------
  74. !
  75. IF( nn_timing == 1 ) CALL timing_start('ssh_nxt')
  76. !
  77. CALL wrk_alloc( jpi, jpj, zhdiv )
  78. !
  79. IF( kt == nit000 ) THEN
  80. !
  81. IF(lwp) WRITE(numout,*)
  82. IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
  83. IF(lwp) WRITE(numout,*) '~~~~~~~ '
  84. !
  85. ENDIF
  86. !
  87. CALL div_cur( kt ) ! Horizontal divergence & Relative vorticity
  88. !
  89. z2dt = 2._wp * rdt ! set time step size (Euler/Leapfrog)
  90. IF( neuler == 0 .AND. kt == nit000 ) z2dt = rdt
  91. ! !------------------------------!
  92. ! ! After Sea Surface Height !
  93. ! !------------------------------!
  94. zhdiv(:,:) = 0._wp
  95. DO jk = 1, jpkm1 ! Horizontal divergence of barotropic transports
  96. zhdiv(:,:) = zhdiv(:,:) + fse3t_n(:,:,jk) * hdivn(:,:,jk)
  97. END DO
  98. ! ! Sea surface elevation time stepping
  99. ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
  100. ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
  101. !
  102. z1_rau0 = 0.5_wp * r1_rau0
  103. ssha(:,:) = ( sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) ) ) * ssmask(:,:)
  104. #if ! defined key_dynspg_ts
  105. ! These lines are not necessary with time splitting since
  106. ! boundary condition on sea level is set during ts loop
  107. #if defined key_agrif
  108. CALL agrif_ssh( kt )
  109. #endif
  110. #if defined key_bdy
  111. IF (lk_bdy) THEN
  112. CALL lbc_lnk( ssha, 'T', 1. ) ! Not sure that's necessary
  113. CALL bdy_ssh( ssha ) ! Duplicate sea level across open boundaries
  114. ENDIF
  115. #endif
  116. #endif
  117. #if defined key_asminc
  118. ! ! Include the IAU weighted SSH increment
  119. IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN
  120. CALL ssh_asm_inc( kt )
  121. ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
  122. ENDIF
  123. #endif
  124. ! !------------------------------!
  125. ! ! outputs !
  126. ! !------------------------------!
  127. !
  128. IF(ln_ctl) CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha - : ', mask1=tmask, ovlap=1 )
  129. !
  130. CALL wrk_dealloc( jpi, jpj, zhdiv )
  131. !
  132. IF( nn_timing == 1 ) CALL timing_stop('ssh_nxt')
  133. !
  134. END SUBROUTINE ssh_nxt
  135. SUBROUTINE wzv( kt )
  136. !!----------------------------------------------------------------------
  137. !! *** ROUTINE wzv ***
  138. !!
  139. !! ** Purpose : compute the now vertical velocity
  140. !!
  141. !! ** Method : - Using the incompressibility hypothesis, the vertical
  142. !! velocity is computed by integrating the horizontal divergence
  143. !! from the bottom to the surface minus the scale factor evolution.
  144. !! The boundary conditions are w=0 at the bottom (no flux) and.
  145. !!
  146. !! ** action : wn : now vertical velocity
  147. !!
  148. !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling.
  149. !!----------------------------------------------------------------------
  150. !
  151. INTEGER, INTENT(in) :: kt ! time step
  152. REAL(wp), POINTER, DIMENSION(:,: ) :: z2d
  153. REAL(wp), POINTER, DIMENSION(:,:,:) :: z3d, zhdiv
  154. !
  155. INTEGER :: ji, jj, jk ! dummy loop indices
  156. REAL(wp) :: z1_2dt ! local scalars
  157. !!----------------------------------------------------------------------
  158. IF( nn_timing == 1 ) CALL timing_start('wzv')
  159. !
  160. IF( kt == nit000 ) THEN
  161. !
  162. IF(lwp) WRITE(numout,*)
  163. IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
  164. IF(lwp) WRITE(numout,*) '~~~~~ '
  165. !
  166. wn(:,:,jpk) = 0._wp ! bottom boundary condition: w=0 (set once for all)
  167. !
  168. ENDIF
  169. ! !------------------------------!
  170. ! ! Now Vertical Velocity !
  171. ! !------------------------------!
  172. z1_2dt = 1. / ( 2. * rdt ) ! set time step size (Euler/Leapfrog)
  173. IF( neuler == 0 .AND. kt == nit000 ) z1_2dt = 1. / rdt
  174. !
  175. IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN ! z_tilde and layer cases
  176. CALL wrk_alloc( jpi, jpj, jpk, zhdiv )
  177. !
  178. DO jk = 1, jpkm1
  179. ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
  180. ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
  181. DO jj = 2, jpjm1
  182. DO ji = fs_2, fs_jpim1 ! vector opt.
  183. zhdiv(ji,jj,jk) = r1_e12t(ji,jj) * ( un_td(ji,jj,jk) - un_td(ji-1,jj,jk) + vn_td(ji,jj,jk) - vn_td(ji,jj-1,jk) )
  184. END DO
  185. END DO
  186. END DO
  187. CALL lbc_lnk(zhdiv, 'T', 1.) ! - ML - Perhaps not necessary: not used for horizontal "connexions"
  188. ! ! Is it problematic to have a wrong vertical velocity in boundary cells?
  189. ! ! Same question holds for hdivn. Perhaps just for security
  190. DO jk = jpkm1, 1, -1 ! integrate from the bottom the hor. divergence
  191. ! computation of w
  192. wn(:,:,jk) = wn(:,:,jk+1) - ( fse3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk) &
  193. & + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
  194. END DO
  195. ! IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
  196. CALL wrk_dealloc( jpi, jpj, jpk, zhdiv )
  197. ELSE ! z_star and linear free surface cases
  198. DO jk = jpkm1, 1, -1 ! integrate from the bottom the hor. divergence
  199. ! computation of w
  200. wn(:,:,jk) = wn(:,:,jk+1) - ( fse3t_n(:,:,jk) * hdivn(:,:,jk) &
  201. & + z1_2dt * ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) ) * tmask(:,:,jk)
  202. END DO
  203. ENDIF
  204. #if defined key_bdy
  205. IF (lk_bdy) THEN
  206. DO jk = 1, jpkm1
  207. wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
  208. END DO
  209. ENDIF
  210. #endif
  211. !
  212. IF( nn_timing == 1 ) CALL timing_stop('wzv')
  213. END SUBROUTINE wzv
  214. SUBROUTINE ssh_swp( kt )
  215. !!----------------------------------------------------------------------
  216. !! *** ROUTINE ssh_nxt ***
  217. !!
  218. !! ** Purpose : achieve the sea surface height time stepping by
  219. !! applying Asselin time filter and swapping the arrays
  220. !! ssha already computed in ssh_nxt
  221. !!
  222. !! ** Method : - apply Asselin time fiter to now ssh (excluding the forcing
  223. !! from the filter, see Leclair and Madec 2010) and swap :
  224. !! sshn = ssha + atfp * ( sshb -2 sshn + ssha )
  225. !! - atfp * rdt * ( emp_b - emp ) / rau0
  226. !! sshn = ssha
  227. !!
  228. !! ** action : - sshb, sshn : before & now sea surface height
  229. !! ready for the next time step
  230. !!
  231. !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling.
  232. !!----------------------------------------------------------------------
  233. INTEGER, INTENT(in) :: kt ! ocean time-step index
  234. !!----------------------------------------------------------------------
  235. !
  236. IF( nn_timing == 1 ) CALL timing_start('ssh_swp')
  237. !
  238. IF( kt == nit000 ) THEN
  239. IF(lwp) WRITE(numout,*)
  240. IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
  241. IF(lwp) WRITE(numout,*) '~~~~~~~ '
  242. ENDIF
  243. # if defined key_dynspg_ts
  244. IF( ( neuler == 0 .AND. kt == nit000 ) .OR. ln_bt_fw ) THEN !** Euler time-stepping: no filter
  245. # else
  246. IF ( neuler == 0 .AND. kt == nit000 ) THEN !** Euler time-stepping at first time-step : no filter
  247. #endif
  248. sshb(:,:) = sshn(:,:) ! before <-- now
  249. sshn(:,:) = ssha(:,:) ! now <-- after (before already = now)
  250. ELSE !** Leap-Frog time-stepping: Asselin filter + swap
  251. sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) ) ! before <-- now filtered
  252. IF( lk_vvl ) sshb(:,:) = sshb(:,:) - atfp * rdt / rau0 * ( emp_b(:,:) - emp(:,:) &
  253. & - rnf_b(:,:) + rnf(:,:) &
  254. & + fwfisf_b(:,:) - fwfisf(:,:) ) * ssmask(:,:)
  255. sshn(:,:) = ssha(:,:) ! now <-- after
  256. ENDIF
  257. !
  258. IF(ln_ctl) CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb - : ', mask1=tmask, ovlap=1 )
  259. !
  260. IF( nn_timing == 1 ) CALL timing_stop('ssh_swp')
  261. !
  262. END SUBROUTINE ssh_swp
  263. !!======================================================================
  264. END MODULE sshwzv