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