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.5 ! 2013-07 (J. Chanut) Compliant with time splitting changes |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | !! ssh_wzv : after ssh & now vertical velocity |
---|
15 | !! ssh_nxt : filter ans swap the ssh arrays |
---|
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 iom ! I/O library |
---|
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 obc_par ! open boundary cond. parameter |
---|
29 | USE obc_oce |
---|
30 | USE bdy_oce |
---|
31 | USE bdy_par |
---|
32 | USE bdydyn2d ! bdy_ssh routine |
---|
33 | USE diaar5, ONLY: lk_diaar5 |
---|
34 | USE iom |
---|
35 | USE sbcrnf ! River runoff |
---|
36 | #if defined key_agrif |
---|
37 | USE agrif_opa_update |
---|
38 | USE agrif_opa_interp |
---|
39 | #endif |
---|
40 | #if defined key_asminc |
---|
41 | USE asminc ! Assimilation increment |
---|
42 | #endif |
---|
43 | USE wrk_nemo ! Memory Allocation |
---|
44 | USE timing ! Timing |
---|
45 | |
---|
46 | #if defined key_dynspg_ts |
---|
47 | ! jchanut: Needed modules for dynamics update: |
---|
48 | USE eosbn2 ! equation of state (eos_bn2 routine) |
---|
49 | USE zpshde ! partial step: hor. derivative (zps_hde routine) |
---|
50 | USE dynadv ! advection (dyn_adv routine) |
---|
51 | USE dynvor ! vorticity term (dyn_vor routine) |
---|
52 | USE dynhpg ! hydrostatic pressure grad. (dyn_hpg routine) |
---|
53 | USE dynldf ! lateral momentum diffusion (dyn_ldf routine) |
---|
54 | USE dynspg_oce ! surface pressure gradient (dyn_spg routine) |
---|
55 | USE dynspg ! surface pressure gradient (dyn_spg routine) |
---|
56 | USE dynnept ! simp. form of Neptune effect(dyn_nept_cor routine) |
---|
57 | USE asminc ! assimilation increments (dyn_asm_inc routine) |
---|
58 | USE bdydyn3d ! (bdy_dyn3d_dmp routine) |
---|
59 | USE dynspg_ts ! for advective velocities issued from time splitting |
---|
60 | USE zdfbfr ! Bottom friction |
---|
61 | #if defined key_agrif |
---|
62 | USE agrif_opa_sponge ! Momemtum and tracers sponges |
---|
63 | #endif |
---|
64 | #endif |
---|
65 | |
---|
66 | IMPLICIT NONE |
---|
67 | PRIVATE |
---|
68 | |
---|
69 | PUBLIC ssh_wzv ! called by step.F90 |
---|
70 | PUBLIC ssh_nxt ! called by step.F90 |
---|
71 | |
---|
72 | !! * Substitutions |
---|
73 | # include "domzgr_substitute.h90" |
---|
74 | # include "vectopt_loop_substitute.h90" |
---|
75 | !!---------------------------------------------------------------------- |
---|
76 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
77 | !! $Id$ |
---|
78 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
79 | !!---------------------------------------------------------------------- |
---|
80 | CONTAINS |
---|
81 | |
---|
82 | SUBROUTINE ssh_wzv( kt ) |
---|
83 | !!---------------------------------------------------------------------- |
---|
84 | !! *** ROUTINE ssh_wzv *** |
---|
85 | !! |
---|
86 | !! ** Purpose : compute the after ssh (ssha), the now vertical velocity |
---|
87 | !! and update the now vertical coordinate (lk_vvl=T). |
---|
88 | !! |
---|
89 | !! ** Method : - Using the incompressibility hypothesis, the vertical |
---|
90 | !! velocity is computed by integrating the horizontal divergence |
---|
91 | !! from the bottom to the surface minus the scale factor evolution. |
---|
92 | !! The boundary conditions are w=0 at the bottom (no flux) and. |
---|
93 | !! |
---|
94 | !! ** action : ssha : after sea surface height |
---|
95 | !! wn : now vertical velocity |
---|
96 | !! sshu_a, sshv_a, sshf_a : after sea surface height (lk_vvl=T) |
---|
97 | !! hu, hv, hur, hvr : ocean depth and its inverse at u-,v-points |
---|
98 | !! |
---|
99 | !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling. |
---|
100 | !!---------------------------------------------------------------------- |
---|
101 | INTEGER, INTENT(in) :: kt ! time step |
---|
102 | ! |
---|
103 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
104 | INTEGER :: indic |
---|
105 | REAL(wp) :: zcoefu, zcoefv, zcoeff, z2dt, z1_2dt, z1_rau0 ! local scalars |
---|
106 | REAL(wp), POINTER, DIMENSION(:,: ) :: z2d, zhdiv |
---|
107 | REAL(wp), POINTER, DIMENSION(:,:,:) :: z3d |
---|
108 | !!---------------------------------------------------------------------- |
---|
109 | ! |
---|
110 | IF( nn_timing == 1 ) CALL timing_start('ssh_wzv') |
---|
111 | ! |
---|
112 | CALL wrk_alloc( jpi, jpj, z2d, zhdiv ) |
---|
113 | ! |
---|
114 | IF( kt == nit000 ) THEN |
---|
115 | ! |
---|
116 | IF(lwp) WRITE(numout,*) |
---|
117 | IF(lwp) WRITE(numout,*) 'ssh_wzv : after sea surface height and now vertical velocity ' |
---|
118 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
119 | ! |
---|
120 | wn(:,:,jpk) = 0._wp ! bottom boundary condition: w=0 (set once for all) |
---|
121 | ! |
---|
122 | IF( lk_vvl ) THEN ! before and now Sea SSH at u-, v-, f-points (vvl case only) |
---|
123 | DO jj = 1, jpjm1 |
---|
124 | DO ji = 1, jpim1 ! caution: use of Vector Opt. not possible |
---|
125 | zcoefu = 0.5 * umask(ji,jj,1) / ( e1u(ji,jj) * e2u(ji,jj) ) |
---|
126 | zcoefv = 0.5 * vmask(ji,jj,1) / ( e1v(ji,jj) * e2v(ji,jj) ) |
---|
127 | zcoeff = 0.25 * umask(ji,jj,1) * umask(ji,jj+1,1) |
---|
128 | sshu_b(ji,jj) = zcoefu * ( e1t(ji ,jj) * e2t(ji ,jj) * sshb(ji ,jj) & |
---|
129 | & + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) ) |
---|
130 | sshv_b(ji,jj) = zcoefv * ( e1t(ji,jj ) * e2t(ji,jj ) * sshb(ji,jj ) & |
---|
131 | & + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) ) |
---|
132 | sshu_n(ji,jj) = zcoefu * ( e1t(ji ,jj) * e2t(ji ,jj) * sshn(ji ,jj) & |
---|
133 | & + e1t(ji+1,jj) * e2t(ji+1,jj) * sshn(ji+1,jj) ) |
---|
134 | sshv_n(ji,jj) = zcoefv * ( e1t(ji,jj ) * e2t(ji,jj ) * sshn(ji,jj ) & |
---|
135 | & + e1t(ji,jj+1) * e2t(ji,jj+1) * sshn(ji,jj+1) ) |
---|
136 | END DO |
---|
137 | END DO |
---|
138 | CALL lbc_lnk( sshu_b, 'U', 1. ) ; CALL lbc_lnk( sshu_n, 'U', 1. ) |
---|
139 | CALL lbc_lnk( sshv_b, 'V', 1. ) ; CALL lbc_lnk( sshv_n, 'V', 1. ) |
---|
140 | DO jj = 1, jpjm1 |
---|
141 | DO ji = 1, jpim1 ! NO Vector Opt. |
---|
142 | sshf_n(ji,jj) = 0.5 * umask(ji,jj,1) * umask(ji,jj+1,1) & |
---|
143 | & / ( e1f(ji,jj ) * e2f(ji,jj ) ) & |
---|
144 | & * ( e1u(ji,jj ) * e2u(ji,jj ) * sshu_n(ji,jj ) & |
---|
145 | & + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) ) |
---|
146 | END DO |
---|
147 | END DO |
---|
148 | CALL lbc_lnk( sshf_n, 'F', 1. ) |
---|
149 | ENDIF |
---|
150 | ! |
---|
151 | ENDIF |
---|
152 | |
---|
153 | ! !------------------------------------------! |
---|
154 | IF( lk_vvl ) THEN ! Regridding: Update Now Vertical coord. ! (only in vvl case) |
---|
155 | ! !------------------------------------------! |
---|
156 | DO jk = 1, jpkm1 |
---|
157 | fsdept(:,:,jk) = fsdept_n(:,:,jk) ! now local depths stored in fsdep. arrays |
---|
158 | fsdepw(:,:,jk) = fsdepw_n(:,:,jk) |
---|
159 | fsde3w(:,:,jk) = fsde3w_n(:,:,jk) |
---|
160 | ! |
---|
161 | fse3t (:,:,jk) = fse3t_n (:,:,jk) ! vertical scale factors stored in fse3. arrays |
---|
162 | fse3u (:,:,jk) = fse3u_n (:,:,jk) |
---|
163 | fse3v (:,:,jk) = fse3v_n (:,:,jk) |
---|
164 | fse3f (:,:,jk) = fse3f_n (:,:,jk) |
---|
165 | fse3w (:,:,jk) = fse3w_n (:,:,jk) |
---|
166 | fse3uw(:,:,jk) = fse3uw_n(:,:,jk) |
---|
167 | fse3vw(:,:,jk) = fse3vw_n(:,:,jk) |
---|
168 | END DO |
---|
169 | ! |
---|
170 | hu(:,:) = hu_0(:,:) + sshu_n(:,:) ! now ocean depth (at u- and v-points) |
---|
171 | hv(:,:) = hv_0(:,:) + sshv_n(:,:) |
---|
172 | ! bg jchanut tschanges |
---|
173 | #if defined key_dynspg_ts |
---|
174 | ht(:,:) = ht_0(:,:) + sshn(:,:) ! now ocean depth (at t- and f-points) |
---|
175 | hf(:,:) = hf_0(:,:) + sshf_n(:,:) |
---|
176 | #endif |
---|
177 | ! end jchanut tschanges |
---|
178 | ! ! now masked inverse of the ocean depth (at u- and v-points) |
---|
179 | hur(:,:) = umask(:,:,1) / ( hu(:,:) + 1._wp - umask(:,:,1) ) |
---|
180 | hvr(:,:) = vmask(:,:,1) / ( hv(:,:) + 1._wp - vmask(:,:,1) ) |
---|
181 | ! |
---|
182 | ENDIF |
---|
183 | ! |
---|
184 | CALL div_cur( kt ) ! Horizontal divergence & Relative vorticity |
---|
185 | ! |
---|
186 | z2dt = 2._wp * rdt ! set time step size (Euler/Leapfrog) |
---|
187 | IF( neuler == 0 .AND. kt == nit000 ) z2dt = rdt |
---|
188 | |
---|
189 | ! !------------------------------! |
---|
190 | ! ! After Sea Surface Height ! |
---|
191 | ! !------------------------------! |
---|
192 | zhdiv(:,:) = 0._wp |
---|
193 | DO jk = 1, jpkm1 ! Horizontal divergence of barotropic transports |
---|
194 | zhdiv(:,:) = zhdiv(:,:) + fse3t(:,:,jk) * hdivn(:,:,jk) |
---|
195 | END DO |
---|
196 | ! ! Sea surface elevation time stepping |
---|
197 | ! In forward Euler time stepping case, the same formulation as in the leap-frog case can be used |
---|
198 | ! because emp_b field is initialized with the vlaues of emp field. Hence, 0.5 * ( emp + emp_b ) = emp |
---|
199 | z1_rau0 = 0.5 / rau0 |
---|
200 | ssha(:,:) = ( sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) ) ) * tmask(:,:,1) |
---|
201 | |
---|
202 | #if defined key_agrif |
---|
203 | CALL agrif_ssh( kt ) |
---|
204 | #endif |
---|
205 | #if defined key_obc |
---|
206 | IF( Agrif_Root() ) THEN |
---|
207 | ssha(:,:) = ssha(:,:) * obctmsk(:,:) |
---|
208 | CALL lbc_lnk( ssha, 'T', 1. ) ! absolutly compulsory !! (jmm) |
---|
209 | ENDIF |
---|
210 | #endif |
---|
211 | #if defined key_bdy |
---|
212 | ! bg jchanut tschanges |
---|
213 | ! These lines are not necessary with time splitting since |
---|
214 | ! boundary condition on sea level is set during ts loop |
---|
215 | IF (lk_bdy) THEN |
---|
216 | CALL lbc_lnk( ssha, 'T', 1. ) ! Not sure that's necessary |
---|
217 | CALL bdy_ssh( ssha ) ! Duplicate sea level across open boundaries |
---|
218 | ENDIF |
---|
219 | #endif |
---|
220 | ! end jchanut tschanges |
---|
221 | #if defined key_asminc |
---|
222 | ! ! Include the IAU weighted SSH increment |
---|
223 | IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN |
---|
224 | CALL ssh_asm_inc( kt ) |
---|
225 | ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:) |
---|
226 | ENDIF |
---|
227 | #endif |
---|
228 | ! !------------------------------! |
---|
229 | ! ! Now Vertical Velocity ! |
---|
230 | ! !------------------------------! |
---|
231 | z1_2dt = 1.e0 / z2dt |
---|
232 | DO jk = jpkm1, 1, -1 ! integrate from the bottom the hor. divergence |
---|
233 | ! - ML - need 3 lines here because replacement of fse3t by its expression yields too long lines otherwise |
---|
234 | wn(:,:,jk) = wn(:,:,jk+1) - fse3t_n(:,:,jk) * hdivn(:,:,jk) & |
---|
235 | & - ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) & |
---|
236 | & * tmask(:,:,jk) * z1_2dt |
---|
237 | #if defined key_bdy |
---|
238 | wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:) |
---|
239 | #endif |
---|
240 | END DO |
---|
241 | |
---|
242 | ! bg jchanut tschanges |
---|
243 | #if defined key_dynspg_ts |
---|
244 | ! In case the time splitting case, update most of all momentum trends here: |
---|
245 | ! Note that the computation of vertical velocity above, hence "after" sea level is necessary |
---|
246 | ! to compute momentum advection for the rhs of barotropic loop: |
---|
247 | CALL eos ( tsn, rhd, rhop ) ! now in situ density for hpg computation |
---|
248 | |
---|
249 | IF( ln_zps ) CALL zps_hde( kt, jpts, tsn, gtsu, gtsv, & ! zps: now hor. derivative |
---|
250 | & rhd, gru , grv ) ! of t, s, rd at the last ocean level |
---|
251 | CALL zdf_bfr( kt ) ! bottom friction (if quadratic) |
---|
252 | |
---|
253 | ua(:,:,:) = 0.e0 ! set dynamics trends to zero |
---|
254 | va(:,:,:) = 0.e0 |
---|
255 | IF( ln_asmiau .AND. & |
---|
256 | & ln_dyninc ) CALL dyn_asm_inc( kt ) ! apply dynamics assimilation increment |
---|
257 | IF( ln_neptsimp ) CALL dyn_nept_cor( kt ) ! subtract Neptune velocities (simplified) |
---|
258 | IF( lk_bdy ) CALL bdy_dyn3d_dmp( kt ) ! bdy damping trends |
---|
259 | CALL dyn_adv( kt ) ! advection (vector or flux form) |
---|
260 | CALL dyn_vor( kt ) ! vorticity term including Coriolis |
---|
261 | CALL dyn_ldf( kt ) ! lateral mixing |
---|
262 | IF( ln_neptsimp ) CALL dyn_nept_cor( kt ) ! add Neptune velocities (simplified) |
---|
263 | #if defined key_agrif |
---|
264 | IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_dyn ! momentum sponge |
---|
265 | #endif |
---|
266 | CALL dyn_hpg( kt ) ! horizontal gradient of Hydrostatic pressure |
---|
267 | CALL dyn_spg( kt, indic ) ! surface pressure gradient |
---|
268 | |
---|
269 | ua_bak(:,:,:) = ua(:,:,:) ! save next velocities (not trends !) |
---|
270 | va_bak(:,:,:) = va(:,:,:) |
---|
271 | |
---|
272 | ! At this stage: |
---|
273 | ! 1) ssha, sshu_a, sshv_a have been updated. |
---|
274 | ! 2) Time averaged barotropic velocities around step kt+1 (ua_b, va_b) as well as |
---|
275 | ! "advective" barotropic velocities (un_adv, vn_adv) at step kt in agreement |
---|
276 | ! with continuity equation are available. |
---|
277 | ! 3) 3d velocities (ua, va) hold ua_b, va_b issued from time splitting as barotropic components. |
---|
278 | ! Below => Update "now" velocities, divergence, then vertical velocity |
---|
279 | ! NB: hdivn is NOT updated such that hdivb is not updated also. This means that horizontal |
---|
280 | ! momentum diffusion is still performed on Instantaneous barotropic arrays. I experencied |
---|
281 | ! some issues with UBS with the current method. Uncomment "divup" below to update the divergence. |
---|
282 | ! |
---|
283 | CALL wrk_alloc( jpi,jpj,jpk, z3d ) |
---|
284 | ! |
---|
285 | DO jk = 1, jpkm1 |
---|
286 | ! Correct velocities: |
---|
287 | un(:,:,jk) = ( un(:,:,jk) + un_adv(:,:) - un_b(:,:) )*umask(:,:,jk) |
---|
288 | vn(:,:,jk) = ( vn(:,:,jk) + vn_adv(:,:) - vn_b(:,:) )*vmask(:,:,jk) |
---|
289 | ! |
---|
290 | ! Compute divergence: |
---|
291 | DO jj = 2, jpjm1 |
---|
292 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
293 | z3d(ji,jj,jk) = & |
---|
294 | ( e2u(ji,jj)*fse3u(ji,jj,jk) * un(ji,jj,jk) - e2u(ji-1,jj)*fse3u(ji-1,jj,jk) * un(ji-1,jj,jk) & |
---|
295 | + e1v(ji,jj)*fse3v(ji,jj,jk) * vn(ji,jj,jk) - e1v(ji,jj-1)*fse3v(ji,jj-1,jk) * vn(ji,jj-1,jk) ) & |
---|
296 | / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) |
---|
297 | END DO |
---|
298 | END DO |
---|
299 | END DO |
---|
300 | ! |
---|
301 | IF( ln_rnf ) CALL sbc_rnf_div( z3d ) ! runoffs (update divergence) |
---|
302 | ! |
---|
303 | ! |
---|
304 | ! Set new vertical velocities: |
---|
305 | DO jk = jpkm1, 1, -1 ! integrate from the bottom the hor. divergence |
---|
306 | ! - ML - need 3 lines here because replacement of fse3t by its expression yields too long lines otherwise |
---|
307 | wn(:,:,jk) = wn(:,:,jk+1) - fse3t_n(:,:,jk) * z3d(:,:,jk) & |
---|
308 | & - (fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) & |
---|
309 | & * tmask(:,:,jk) * z1_2dt |
---|
310 | #if defined key_bdy |
---|
311 | ! JC: line below is purely cosmetic I guess: |
---|
312 | wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:) |
---|
313 | #endif |
---|
314 | END DO |
---|
315 | |
---|
316 | #if defined key_agrif |
---|
317 | ! Set verticaly velocity to zero along open boundaries (cosmetic) |
---|
318 | IF( .NOT. AGRIF_Root() ) THEN |
---|
319 | DO jk = 1, jpkm1 |
---|
320 | IF ((nbondi == 1).OR.(nbondi == 2)) wn(nlci-1 , : ,jk) = 0.e0 ! east |
---|
321 | IF ((nbondi == -1).OR.(nbondi == 2)) wn(2 , : ,jk) = 0.e0 ! west |
---|
322 | IF ((nbondj == 1).OR.(nbondj == 2)) wn(: ,nlcj-1 ,jk) = 0.e0 ! north |
---|
323 | IF ((nbondj == -1).OR.(nbondj == 2)) wn(: ,2 ,jk) = 0.e0 ! south |
---|
324 | END DO |
---|
325 | ENDIF |
---|
326 | #endif |
---|
327 | ! |
---|
328 | CALL wrk_dealloc( jpi,jpj,jpk, z3d ) |
---|
329 | |
---|
330 | !! bg divup |
---|
331 | !! ! |
---|
332 | !! DO jk = 1, jpkm1 |
---|
333 | !! ! Correct velocities: |
---|
334 | !! un(:,:,jk) = ( un(:,:,jk) + un_adv(:,:) - un_b(:,:) )*umask(:,:,jk) |
---|
335 | !! vn(:,:,jk) = ( vn(:,:,jk) + vn_adv(:,:) - vn_b(:,:) )*vmask(:,:,jk) |
---|
336 | !! ! |
---|
337 | !! END DO |
---|
338 | !! ! |
---|
339 | !! ! |
---|
340 | !! CALL div_cur( kt ) ! Horizontal divergence & Relative vorticity |
---|
341 | !! ! |
---|
342 | !! ! Set new vertical velocities: |
---|
343 | !! DO jk = jpkm1, 1, -1 ! integrate from the bottom the hor. divergence |
---|
344 | !! ! - ML - need 3 lines here because replacement of fse3t by its expression yields too long lines otherwise |
---|
345 | !! wn(:,:,jk) = wn(:,:,jk+1) - fse3t_n(:,:,jk) * hdivn(:,:,jk) & |
---|
346 | !! & - (fse3t_a(:,:,jk) - fse3t_b(:,:,jk) ) & |
---|
347 | !! & * tmask(:,:,jk) * z1_2dt |
---|
348 | !!#if defined key_bdy |
---|
349 | !! ! JC: line below is purely cosmetic I guess: |
---|
350 | !! wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:) |
---|
351 | !!#endif |
---|
352 | !! END DO |
---|
353 | !! end divup |
---|
354 | ! |
---|
355 | ! |
---|
356 | ! End of time splitting case |
---|
357 | #else |
---|
358 | ! ! Sea Surface Height at u-,v- and f-points (vvl case only) |
---|
359 | IF( lk_vvl ) THEN ! (required only in key_vvl case) |
---|
360 | DO jj = 1, jpjm1 |
---|
361 | DO ji = 1, jpim1 ! NO Vector Opt. |
---|
362 | sshu_a(ji,jj) = 0.5 * umask(ji,jj,1) / ( e1u(ji ,jj) * e2u(ji ,jj) ) & |
---|
363 | & * ( e1t(ji ,jj) * e2t(ji ,jj) * ssha(ji ,jj) & |
---|
364 | & + e1t(ji+1,jj) * e2t(ji+1,jj) * ssha(ji+1,jj) ) |
---|
365 | sshv_a(ji,jj) = 0.5 * vmask(ji,jj,1) / ( e1v(ji,jj ) * e2v(ji,jj ) ) & |
---|
366 | & * ( e1t(ji,jj ) * e2t(ji,jj ) * ssha(ji,jj ) & |
---|
367 | & + e1t(ji,jj+1) * e2t(ji,jj+1) * ssha(ji,jj+1) ) |
---|
368 | END DO |
---|
369 | END DO |
---|
370 | CALL lbc_lnk( sshu_a, 'U', 1. ) ; CALL lbc_lnk( sshv_a, 'V', 1. ) ! Boundaries conditions |
---|
371 | ENDIF |
---|
372 | #endif |
---|
373 | ! !------------------------------! |
---|
374 | ! ! outputs ! |
---|
375 | ! !------------------------------! |
---|
376 | CALL iom_put( "woce", wn ) ! vertical velocity |
---|
377 | CALL iom_put( "ssh" , sshn ) ! sea surface height |
---|
378 | CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) ) ! square of sea surface height |
---|
379 | IF( lk_diaar5 ) THEN ! vertical mass transport & its square value |
---|
380 | ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. |
---|
381 | CALL wrk_alloc( jpi,jpj,jpk, z3d ) |
---|
382 | z2d(:,:) = rau0 * e1t(:,:) * e2t(:,:) |
---|
383 | DO jk = 1, jpk |
---|
384 | z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) |
---|
385 | END DO |
---|
386 | CALL iom_put( "w_masstr" , z3d ) |
---|
387 | CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) |
---|
388 | CALL wrk_dealloc( jpi,jpj,jpk, z3d ) |
---|
389 | ENDIF |
---|
390 | ! |
---|
391 | IF(ln_ctl) CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha - : ', mask1=tmask, ovlap=1 ) |
---|
392 | ! |
---|
393 | CALL wrk_dealloc( jpi, jpj, z2d, zhdiv ) |
---|
394 | ! |
---|
395 | IF( nn_timing == 1 ) CALL timing_stop('ssh_wzv') |
---|
396 | ! |
---|
397 | END SUBROUTINE ssh_wzv |
---|
398 | |
---|
399 | |
---|
400 | SUBROUTINE ssh_nxt( kt ) |
---|
401 | !!---------------------------------------------------------------------- |
---|
402 | !! *** ROUTINE ssh_nxt *** |
---|
403 | !! |
---|
404 | !! ** Purpose : achieve the sea surface height time stepping by |
---|
405 | !! applying Asselin time filter and swapping the arrays |
---|
406 | !! ssha already computed in ssh_wzv |
---|
407 | !! |
---|
408 | !! ** Method : - apply Asselin time fiter to now ssh (excluding the forcing |
---|
409 | !! from the filter, see Leclair and Madec 2010) and swap : |
---|
410 | !! sshn = ssha + atfp * ( sshb -2 sshn + ssha ) |
---|
411 | !! - atfp * rdt * ( emp_b - emp ) / rau0 |
---|
412 | !! sshn = ssha |
---|
413 | !! |
---|
414 | !! ** action : - sshb, sshn : before & now sea surface height |
---|
415 | !! ready for the next time step |
---|
416 | !! |
---|
417 | !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling. |
---|
418 | !!---------------------------------------------------------------------- |
---|
419 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
420 | !! |
---|
421 | INTEGER :: ji, jj ! dummy loop indices |
---|
422 | REAL(wp) :: zec ! temporary scalar |
---|
423 | !!---------------------------------------------------------------------- |
---|
424 | ! |
---|
425 | IF( nn_timing == 1 ) CALL timing_start('ssh_nxt') |
---|
426 | ! |
---|
427 | IF( kt == nit000 ) THEN |
---|
428 | IF(lwp) WRITE(numout,*) |
---|
429 | IF(lwp) WRITE(numout,*) 'ssh_nxt : next sea surface height (Asselin time filter + swap)' |
---|
430 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
431 | ENDIF |
---|
432 | |
---|
433 | ! !--------------------------! |
---|
434 | IF( lk_vvl ) THEN ! Variable volume levels ! (ssh at t-, u-, v, f-points) |
---|
435 | ! !--------------------------! |
---|
436 | ! |
---|
437 | #if defined key_dynspg_ts |
---|
438 | IF(( neuler == 0 .AND. kt == nit000 ).OR.(lk_dynspg_ts.AND.ln_bt_fw)) THEN !** Euler time-stepping: no filter |
---|
439 | #else |
---|
440 | IF ( neuler == 0 .AND. kt == nit000 ) THEN |
---|
441 | #endif |
---|
442 | sshb (:,:) = sshn (:,:) |
---|
443 | sshn (:,:) = ssha (:,:) ! now <-- after (before already = now) |
---|
444 | sshu_b(:,:) = sshu_n(:,:) |
---|
445 | sshv_b(:,:) = sshv_n(:,:) |
---|
446 | sshu_n(:,:) = sshu_a(:,:) |
---|
447 | sshv_n(:,:) = sshv_a(:,:) |
---|
448 | DO jj = 1, jpjm1 ! ssh now at f-point |
---|
449 | DO ji = 1, jpim1 ! NO Vector Opt. |
---|
450 | sshf_n(ji,jj) = 0.5 * umask(ji,jj,1) * umask(ji,jj+1,1) & |
---|
451 | & / ( e1f(ji,jj ) * e2f(ji,jj ) ) & |
---|
452 | & * ( e1u(ji,jj ) * e2u(ji,jj ) * sshu_n(ji,jj ) & |
---|
453 | & + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) ) |
---|
454 | END DO |
---|
455 | END DO |
---|
456 | CALL lbc_lnk( sshf_n, 'F', 1. ) ! Boundaries conditions |
---|
457 | ! |
---|
458 | ELSE !** Leap-Frog time-stepping: Asselin filter + swap |
---|
459 | zec = atfp * rdt / rau0 |
---|
460 | DO jj = 1, jpj |
---|
461 | DO ji = 1, jpi ! before <-- now filtered |
---|
462 | sshb (ji,jj) = sshn (ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) ) & |
---|
463 | & - zec * ( emp_b(ji,jj) - emp(ji,jj) ) * tmask(ji,jj,1) |
---|
464 | sshn (ji,jj) = ssha (ji,jj) ! now <-- after |
---|
465 | sshu_n(ji,jj) = sshu_a(ji,jj) |
---|
466 | sshv_n(ji,jj) = sshv_a(ji,jj) |
---|
467 | END DO |
---|
468 | END DO |
---|
469 | DO jj = 1, jpjm1 ! ssh now at f-point |
---|
470 | DO ji = 1, jpim1 ! NO Vector Opt. |
---|
471 | sshf_n(ji,jj) = 0.5 * umask(ji,jj,1) * umask(ji,jj+1,1) & |
---|
472 | & / ( e1f(ji,jj ) * e2f(ji,jj ) ) & |
---|
473 | & * ( e1u(ji,jj ) * e2u(ji,jj ) * sshu_n(ji,jj ) & |
---|
474 | & + e1u(ji,jj+1) * e2u(ji,jj+1) * sshu_n(ji,jj+1) ) |
---|
475 | END DO |
---|
476 | END DO |
---|
477 | CALL lbc_lnk( sshf_n, 'F', 1. ) ! Boundaries conditions |
---|
478 | ! |
---|
479 | DO jj = 1, jpjm1 ! ssh before at u- & v-points |
---|
480 | DO ji = 1, jpim1 ! NO Vector Opt. |
---|
481 | sshu_b(ji,jj) = 0.5 * umask(ji,jj,1) / ( e1u(ji ,jj) * e2u(ji ,jj) ) & |
---|
482 | & * ( e1t(ji ,jj) * e2t(ji ,jj) * sshb(ji ,jj) & |
---|
483 | & + e1t(ji+1,jj) * e2t(ji+1,jj) * sshb(ji+1,jj) ) |
---|
484 | sshv_b(ji,jj) = 0.5 * vmask(ji,jj,1) / ( e1v(ji,jj ) * e2v(ji,jj ) ) & |
---|
485 | & * ( e1t(ji,jj ) * e2t(ji,jj ) * sshb(ji,jj ) & |
---|
486 | & + e1t(ji,jj+1) * e2t(ji,jj+1) * sshb(ji,jj+1) ) |
---|
487 | END DO |
---|
488 | END DO |
---|
489 | CALL lbc_lnk( sshu_b, 'U', 1. ) |
---|
490 | CALL lbc_lnk( sshv_b, 'V', 1. ) ! Boundaries conditions |
---|
491 | ! |
---|
492 | ENDIF |
---|
493 | ! !--------------------------! |
---|
494 | ELSE ! fixed levels ! (ssh at t-point only) |
---|
495 | ! !--------------------------! |
---|
496 | ! |
---|
497 | #if defined key_dynspg_ts |
---|
498 | IF(( neuler == 0 .AND. kt == nit000 ).OR.(lk_dynspg_ts.AND.ln_bt_fw)) THEN !** Euler time-stepping: no filter |
---|
499 | #else |
---|
500 | IF( neuler == 0 .AND. kt == nit000 ) THEN |
---|
501 | #endif |
---|
502 | sshb(:,:) = sshn(:,:) |
---|
503 | sshn(:,:) = ssha(:,:) ! now <-- after (before already = now) |
---|
504 | ! |
---|
505 | ELSE ! Leap-Frog time-stepping: Asselin filter + swap |
---|
506 | DO jj = 1, jpj |
---|
507 | DO ji = 1, jpi ! before <-- now filtered |
---|
508 | sshb(ji,jj) = sshn(ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) ) |
---|
509 | sshn(ji,jj) = ssha(ji,jj) ! now <-- after |
---|
510 | END DO |
---|
511 | END DO |
---|
512 | ENDIF |
---|
513 | ! |
---|
514 | ENDIF |
---|
515 | ! |
---|
516 | ! Update velocity at AGRIF zoom boundaries |
---|
517 | #if defined key_agrif |
---|
518 | IF ( .NOT.Agrif_Root() ) CALL Agrif_Update_Dyn( kt ) |
---|
519 | #endif |
---|
520 | ! |
---|
521 | IF(ln_ctl) CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb - : ', mask1=tmask, ovlap=1 ) |
---|
522 | ! |
---|
523 | IF( nn_timing == 1 ) CALL timing_stop('ssh_nxt') |
---|
524 | ! |
---|
525 | END SUBROUTINE ssh_nxt |
---|
526 | |
---|
527 | !!====================================================================== |
---|
528 | END MODULE sshwzv |
---|