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 | !! 4.0 ! 2018-12 (A. Coward) add mixed implicit/explicit advection |
---|
12 | !! 4.1 ! 2019-08 (A. Coward, D. Storkey) Rename ssh_nxt -> ssh_atf. Now only does time filtering. |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | !! ssh_nxt : after ssh |
---|
17 | !! ssh_atf : time filter the ssh arrays |
---|
18 | !! wzv : compute now vertical velocity |
---|
19 | !!---------------------------------------------------------------------- |
---|
20 | USE oce ! ocean dynamics and tracers variables |
---|
21 | USE isf_oce ! ice shelf |
---|
22 | USE dom_oce ! ocean space and time domain variables |
---|
23 | USE sbc_oce ! surface boundary condition: ocean |
---|
24 | USE domvvl ! Variable volume |
---|
25 | USE divhor ! horizontal divergence |
---|
26 | USE phycst ! physical constants |
---|
27 | USE bdy_oce , ONLY : ln_bdy, bdytmask ! Open BounDarY |
---|
28 | USE bdydyn2d ! bdy_ssh routine |
---|
29 | #if defined key_agrif |
---|
30 | USE agrif_oce_interp |
---|
31 | #endif |
---|
32 | ! |
---|
33 | USE iom |
---|
34 | USE in_out_manager ! I/O manager |
---|
35 | USE restart ! only for lrst_oce |
---|
36 | USE prtctl ! Print control |
---|
37 | USE lbclnk ! ocean lateral boundary condition (or mpp link) |
---|
38 | USE lib_mpp ! MPP library |
---|
39 | USE timing ! Timing |
---|
40 | USE wet_dry ! Wetting/Drying flux limiting |
---|
41 | |
---|
42 | IMPLICIT NONE |
---|
43 | PRIVATE |
---|
44 | |
---|
45 | PUBLIC ssh_nxt ! called by step.F90 |
---|
46 | PUBLIC wzv ! called by step.F90 |
---|
47 | PUBLIC wAimp ! called by step.F90 |
---|
48 | PUBLIC ssh_atf ! called by step.F90 |
---|
49 | |
---|
50 | !! * Substitutions |
---|
51 | # include "do_loop_substitute.h90" |
---|
52 | !!---------------------------------------------------------------------- |
---|
53 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
54 | !! $Id$ |
---|
55 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
56 | !!---------------------------------------------------------------------- |
---|
57 | CONTAINS |
---|
58 | |
---|
59 | SUBROUTINE ssh_nxt( kt, Kbb, Kmm, pssh, Kaa ) |
---|
60 | !!---------------------------------------------------------------------- |
---|
61 | !! *** ROUTINE ssh_nxt *** |
---|
62 | !! |
---|
63 | !! ** Purpose : compute the after ssh (ssh(Kaa)) |
---|
64 | !! |
---|
65 | !! ** Method : - Using the incompressibility hypothesis, the ssh increment |
---|
66 | !! is computed by integrating the horizontal divergence and multiply by |
---|
67 | !! by the time step. |
---|
68 | !! |
---|
69 | !! ** action : ssh(:,:,Kaa), after sea surface height |
---|
70 | !! |
---|
71 | !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling. |
---|
72 | !!---------------------------------------------------------------------- |
---|
73 | INTEGER , INTENT(in ) :: kt ! time step |
---|
74 | INTEGER , INTENT(in ) :: Kbb, Kmm, Kaa ! time level index |
---|
75 | REAL(wp), DIMENSION(jpi,jpj,jpt), INTENT(inout) :: pssh ! sea-surface height |
---|
76 | ! |
---|
77 | INTEGER :: jk ! dummy loop index |
---|
78 | REAL(wp) :: zcoef ! local scalar |
---|
79 | REAL(wp), DIMENSION(jpi,jpj) :: zhdiv ! 2D workspace |
---|
80 | !!---------------------------------------------------------------------- |
---|
81 | ! |
---|
82 | IF( ln_timing ) CALL timing_start('ssh_nxt') |
---|
83 | ! |
---|
84 | IF( kt == nit000 ) THEN |
---|
85 | IF(lwp) WRITE(numout,*) |
---|
86 | IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height' |
---|
87 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
88 | ENDIF |
---|
89 | ! |
---|
90 | zcoef = 0.5_wp * r1_rho0 |
---|
91 | |
---|
92 | ! !------------------------------! |
---|
93 | ! ! After Sea Surface Height ! |
---|
94 | ! !------------------------------! |
---|
95 | IF(ln_wd_il) THEN |
---|
96 | CALL wad_lmt(pssh(:,:,Kbb), zcoef * (emp_b(:,:) + emp(:,:)), rDt, Kmm, uu, vv ) |
---|
97 | ENDIF |
---|
98 | |
---|
99 | CALL div_hor( kt, Kbb, Kmm ) ! Horizontal divergence |
---|
100 | ! |
---|
101 | zhdiv(:,:) = 0._wp |
---|
102 | DO jk = 1, jpkm1 ! Horizontal divergence of barotropic transports |
---|
103 | zhdiv(:,:) = zhdiv(:,:) + e3t(:,:,jk,Kmm) * hdiv(:,:,jk) |
---|
104 | END DO |
---|
105 | ! ! Sea surface elevation time stepping |
---|
106 | ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to |
---|
107 | ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations. |
---|
108 | ! |
---|
109 | pssh(:,:,Kaa) = ( pssh(:,:,Kbb) - rDt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) ) ) * ssmask(:,:) |
---|
110 | ! |
---|
111 | #if defined key_agrif |
---|
112 | Kbb_a = Kbb; Kmm_a = Kmm; Krhs_a = Kaa; CALL agrif_ssh( kt ) |
---|
113 | #endif |
---|
114 | ! |
---|
115 | IF ( .NOT.ln_dynspg_ts ) THEN |
---|
116 | IF( ln_bdy ) THEN |
---|
117 | CALL lbc_lnk( 'sshwzv', pssh(:,:,Kaa), 'T', 1. ) ! Not sure that's necessary |
---|
118 | CALL bdy_ssh( pssh(:,:,Kaa) ) ! Duplicate sea level across open boundaries |
---|
119 | ENDIF |
---|
120 | ENDIF |
---|
121 | ! !------------------------------! |
---|
122 | ! ! outputs ! |
---|
123 | ! !------------------------------! |
---|
124 | ! |
---|
125 | IF(sn_cfctl%l_prtctl) CALL prt_ctl( tab2d_1=pssh(:,:,Kaa), clinfo1=' pssh(:,:,Kaa) - : ', mask1=tmask ) |
---|
126 | ! |
---|
127 | IF( ln_timing ) CALL timing_stop('ssh_nxt') |
---|
128 | ! |
---|
129 | END SUBROUTINE ssh_nxt |
---|
130 | |
---|
131 | |
---|
132 | SUBROUTINE wzv( kt, Kbb, Kmm, pww, Kaa ) |
---|
133 | !!---------------------------------------------------------------------- |
---|
134 | !! *** ROUTINE wzv *** |
---|
135 | !! |
---|
136 | !! ** Purpose : compute the now vertical velocity |
---|
137 | !! |
---|
138 | !! ** Method : - Using the incompressibility hypothesis, the vertical |
---|
139 | !! velocity is computed by integrating the horizontal divergence |
---|
140 | !! from the bottom to the surface minus the scale factor evolution. |
---|
141 | !! The boundary conditions are w=0 at the bottom (no flux) and. |
---|
142 | !! |
---|
143 | !! ** action : pww : now vertical velocity |
---|
144 | !! |
---|
145 | !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling. |
---|
146 | !!---------------------------------------------------------------------- |
---|
147 | INTEGER , INTENT(in) :: kt ! time step |
---|
148 | INTEGER , INTENT(in) :: Kbb, Kmm, Kaa ! time level indices |
---|
149 | REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) :: pww ! now vertical velocity |
---|
150 | ! |
---|
151 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
152 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: zhdiv |
---|
153 | !!---------------------------------------------------------------------- |
---|
154 | ! |
---|
155 | IF( ln_timing ) CALL timing_start('wzv') |
---|
156 | ! |
---|
157 | IF( kt == nit000 ) THEN |
---|
158 | IF(lwp) WRITE(numout,*) |
---|
159 | IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity ' |
---|
160 | IF(lwp) WRITE(numout,*) '~~~~~ ' |
---|
161 | ! |
---|
162 | pww(:,:,jpk) = 0._wp ! bottom boundary condition: w=0 (set once for all) |
---|
163 | ENDIF |
---|
164 | ! !------------------------------! |
---|
165 | ! ! Now Vertical Velocity ! |
---|
166 | ! !------------------------------! |
---|
167 | ! |
---|
168 | IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN ! z_tilde and layer cases |
---|
169 | ALLOCATE( zhdiv(jpi,jpj,jpk) ) |
---|
170 | ! |
---|
171 | DO jk = 1, jpkm1 |
---|
172 | ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t) |
---|
173 | ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way) |
---|
174 | DO_2D_00_00 |
---|
175 | zhdiv(ji,jj,jk) = r1_e1e2t(ji,jj) * ( un_td(ji,jj,jk) - un_td(ji-1,jj,jk) + vn_td(ji,jj,jk) - vn_td(ji,jj-1,jk) ) |
---|
176 | END_2D |
---|
177 | END DO |
---|
178 | CALL lbc_lnk('sshwzv', zhdiv, 'T', 1.) ! - ML - Perhaps not necessary: not used for horizontal "connexions" |
---|
179 | ! ! Is it problematic to have a wrong vertical velocity in boundary cells? |
---|
180 | ! ! Same question holds for hdiv. Perhaps just for security |
---|
181 | DO jk = jpkm1, 1, -1 ! integrate from the bottom the hor. divergence |
---|
182 | ! computation of w |
---|
183 | pww(:,:,jk) = pww(:,:,jk+1) - ( e3t(:,:,jk,Kmm) * hdiv(:,:,jk) + zhdiv(:,:,jk) & |
---|
184 | & + r1_Dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) ) ) * tmask(:,:,jk) |
---|
185 | END DO |
---|
186 | ! IF( ln_vvl_layer ) pww(:,:,:) = 0.e0 |
---|
187 | DEALLOCATE( zhdiv ) |
---|
188 | ELSE ! z_star and linear free surface cases |
---|
189 | DO jk = jpkm1, 1, -1 ! integrate from the bottom the hor. divergence |
---|
190 | ! computation of w |
---|
191 | pww(:,:,jk) = pww(:,:,jk+1) - ( e3t(:,:,jk,Kmm) * hdiv(:,:,jk) & |
---|
192 | & + r1_Dt * ( e3t(:,:,jk,Kaa) - e3t(:,:,jk,Kbb) ) ) * tmask(:,:,jk) |
---|
193 | END DO |
---|
194 | ENDIF |
---|
195 | |
---|
196 | IF( ln_bdy ) THEN |
---|
197 | DO jk = 1, jpkm1 |
---|
198 | pww(:,:,jk) = pww(:,:,jk) * bdytmask(:,:) |
---|
199 | END DO |
---|
200 | ENDIF |
---|
201 | ! |
---|
202 | IF( .NOT. AGRIF_Root() ) THEN |
---|
203 | ! |
---|
204 | ! Mask vertical velocity at first/last columns/row |
---|
205 | ! inside computational domain (cosmetic) |
---|
206 | DO jk = 1, jpkm1 |
---|
207 | ! --- West --- ! |
---|
208 | DO ji = mi0(2+nn_hls), mi1(2+nn_hls) |
---|
209 | DO jj = 1, jpj |
---|
210 | pww(ji,jj,jk) = 0._wp |
---|
211 | END DO |
---|
212 | END DO |
---|
213 | ! |
---|
214 | ! --- East --- ! |
---|
215 | DO ji = mi0(jpiglo-1-nn_hls), mi1(jpiglo-1-nn_hls) |
---|
216 | DO jj = 1, jpj |
---|
217 | pww(ji,jj,jk) = 0._wp |
---|
218 | END DO |
---|
219 | END DO |
---|
220 | ! |
---|
221 | ! --- South --- ! |
---|
222 | DO jj = mj0(2+nn_hls), mj1(2+nn_hls) |
---|
223 | DO ji = 1, jpi |
---|
224 | pww(ji,jj,jk) = 0._wp |
---|
225 | END DO |
---|
226 | END DO |
---|
227 | ! |
---|
228 | ! --- North --- ! |
---|
229 | DO jj = mj0(jpjglo-1-nn_hls), mj1(jpjglo-1-nn_hls) |
---|
230 | DO ji = 1, jpi |
---|
231 | pww(ji,jj,jk) = 0._wp |
---|
232 | END DO |
---|
233 | END DO |
---|
234 | END DO |
---|
235 | ! |
---|
236 | ENDIF |
---|
237 | ! |
---|
238 | IF( ln_timing ) CALL timing_stop('wzv') |
---|
239 | ! |
---|
240 | END SUBROUTINE wzv |
---|
241 | |
---|
242 | |
---|
243 | SUBROUTINE ssh_atf( kt, Kbb, Kmm, Kaa, pssh ) |
---|
244 | !!---------------------------------------------------------------------- |
---|
245 | !! *** ROUTINE ssh_atf *** |
---|
246 | !! |
---|
247 | !! ** Purpose : Apply Asselin time filter to now SSH. |
---|
248 | !! |
---|
249 | !! ** Method : - apply Asselin time fiter to now ssh (excluding the forcing |
---|
250 | !! from the filter, see Leclair and Madec 2010) and swap : |
---|
251 | !! pssh(:,:,Kmm) = pssh(:,:,Kaa) + rn_atfp * ( pssh(:,:,Kbb) -2 pssh(:,:,Kmm) + pssh(:,:,Kaa) ) |
---|
252 | !! - rn_atfp * rn_Dt * ( emp_b - emp ) / rho0 |
---|
253 | !! |
---|
254 | !! ** action : - pssh(:,:,Kmm) time filtered |
---|
255 | !! |
---|
256 | !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling. |
---|
257 | !!---------------------------------------------------------------------- |
---|
258 | INTEGER , INTENT(in ) :: kt ! ocean time-step index |
---|
259 | INTEGER , INTENT(in ) :: Kbb, Kmm, Kaa ! ocean time level indices |
---|
260 | REAL(wp), DIMENSION(jpi,jpj,jpt), INTENT(inout) :: pssh ! SSH field |
---|
261 | ! |
---|
262 | REAL(wp) :: zcoef ! local scalar |
---|
263 | !!---------------------------------------------------------------------- |
---|
264 | ! |
---|
265 | IF( ln_timing ) CALL timing_start('ssh_atf') |
---|
266 | ! |
---|
267 | IF( kt == nit000 ) THEN |
---|
268 | IF(lwp) WRITE(numout,*) |
---|
269 | IF(lwp) WRITE(numout,*) 'ssh_atf : Asselin time filter of sea surface height' |
---|
270 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
271 | ENDIF |
---|
272 | ! !== Euler time-stepping: no filter, just swap ==! |
---|
273 | IF ( .NOT.( l_1st_euler ) ) THEN ! Only do time filtering for leapfrog timesteps |
---|
274 | ! ! filtered "now" field |
---|
275 | pssh(:,:,Kmm) = pssh(:,:,Kmm) + rn_atfp * ( pssh(:,:,Kbb) - 2 * pssh(:,:,Kmm) + pssh(:,:,Kaa) ) |
---|
276 | IF( .NOT.ln_linssh ) THEN ! "now" <-- with forcing removed |
---|
277 | zcoef = rn_atfp * rn_Dt * r1_rho0 |
---|
278 | pssh(:,:,Kmm) = pssh(:,:,Kmm) - zcoef * ( emp_b(:,:) - emp (:,:) & |
---|
279 | & - rnf_b(:,:) + rnf (:,:) & |
---|
280 | & + fwfisf_cav_b(:,:) - fwfisf_cav(:,:) & |
---|
281 | & + fwfisf_par_b(:,:) - fwfisf_par(:,:) ) * ssmask(:,:) |
---|
282 | |
---|
283 | ! ice sheet coupling |
---|
284 | IF ( ln_isf .AND. ln_isfcpl .AND. kt == nit000+1) pssh(:,:,Kbb) = pssh(:,:,Kbb) - rn_atfp * rn_Dt * ( risfcpl_ssh(:,:) - 0.0 ) * ssmask(:,:) |
---|
285 | |
---|
286 | ENDIF |
---|
287 | ENDIF |
---|
288 | ! |
---|
289 | IF(sn_cfctl%l_prtctl) CALL prt_ctl( tab2d_1=pssh(:,:,Kmm), clinfo1=' pssh(:,:,Kmm) - : ', mask1=tmask ) |
---|
290 | ! |
---|
291 | IF( ln_timing ) CALL timing_stop('ssh_atf') |
---|
292 | ! |
---|
293 | END SUBROUTINE ssh_atf |
---|
294 | |
---|
295 | SUBROUTINE wAimp( kt, Kmm ) |
---|
296 | !!---------------------------------------------------------------------- |
---|
297 | !! *** ROUTINE wAimp *** |
---|
298 | !! |
---|
299 | !! ** Purpose : compute the Courant number and partition vertical velocity |
---|
300 | !! if a proportion needs to be treated implicitly |
---|
301 | !! |
---|
302 | !! ** Method : - |
---|
303 | !! |
---|
304 | !! ** action : ww : now vertical velocity (to be handled explicitly) |
---|
305 | !! : wi : now vertical velocity (for implicit treatment) |
---|
306 | !! |
---|
307 | !! Reference : Shchepetkin, A. F. (2015): An adaptive, Courant-number-dependent |
---|
308 | !! implicit scheme for vertical advection in oceanic modeling. |
---|
309 | !! Ocean Modelling, 91, 38-69. |
---|
310 | !!---------------------------------------------------------------------- |
---|
311 | INTEGER, INTENT(in) :: kt ! time step |
---|
312 | INTEGER, INTENT(in) :: Kmm ! time level index |
---|
313 | ! |
---|
314 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
315 | REAL(wp) :: zCu, zcff, z1_e3t ! local scalars |
---|
316 | REAL(wp) , PARAMETER :: Cu_min = 0.15_wp ! local parameters |
---|
317 | REAL(wp) , PARAMETER :: Cu_max = 0.30_wp ! local parameters |
---|
318 | REAL(wp) , PARAMETER :: Cu_cut = 2._wp*Cu_max - Cu_min ! local parameters |
---|
319 | REAL(wp) , PARAMETER :: Fcu = 4._wp*Cu_max*(Cu_max-Cu_min) ! local parameters |
---|
320 | !!---------------------------------------------------------------------- |
---|
321 | ! |
---|
322 | IF( ln_timing ) CALL timing_start('wAimp') |
---|
323 | ! |
---|
324 | IF( kt == nit000 ) THEN |
---|
325 | IF(lwp) WRITE(numout,*) |
---|
326 | IF(lwp) WRITE(numout,*) 'wAimp : Courant number-based partitioning of now vertical velocity ' |
---|
327 | IF(lwp) WRITE(numout,*) '~~~~~ ' |
---|
328 | wi(:,:,:) = 0._wp |
---|
329 | ENDIF |
---|
330 | ! |
---|
331 | ! Calculate Courant numbers |
---|
332 | IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN |
---|
333 | DO_3D_00_00( 1, jpkm1 ) |
---|
334 | z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm) |
---|
335 | ! 2*rn_Dt and not rDt (for restartability) |
---|
336 | Cu_adv(ji,jj,jk) = 2._wp * rn_Dt * ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) ) & |
---|
337 | & + ( MAX( e2u(ji ,jj)*e3u(ji ,jj,jk,Kmm)*uu(ji ,jj,jk,Kmm) + un_td(ji ,jj,jk), 0._wp ) - & |
---|
338 | & MIN( e2u(ji-1,jj)*e3u(ji-1,jj,jk,Kmm)*uu(ji-1,jj,jk,Kmm) + un_td(ji-1,jj,jk), 0._wp ) ) & |
---|
339 | & * r1_e1e2t(ji,jj) & |
---|
340 | & + ( MAX( e1v(ji,jj )*e3v(ji,jj ,jk,Kmm)*vv(ji,jj ,jk,Kmm) + vn_td(ji,jj ,jk), 0._wp ) - & |
---|
341 | & MIN( e1v(ji,jj-1)*e3v(ji,jj-1,jk,Kmm)*vv(ji,jj-1,jk,Kmm) + vn_td(ji,jj-1,jk), 0._wp ) ) & |
---|
342 | & * r1_e1e2t(ji,jj) & |
---|
343 | & ) * z1_e3t |
---|
344 | END_3D |
---|
345 | ELSE |
---|
346 | DO_3D_00_00( 1, jpkm1 ) |
---|
347 | z1_e3t = 1._wp / e3t(ji,jj,jk,Kmm) |
---|
348 | ! 2*rn_Dt and not rDt (for restartability) |
---|
349 | Cu_adv(ji,jj,jk) = 2._wp * rn_Dt * ( ( MAX( ww(ji,jj,jk) , 0._wp ) - MIN( ww(ji,jj,jk+1) , 0._wp ) ) & |
---|
350 | & + ( MAX( e2u(ji ,jj)*e3u(ji ,jj,jk,Kmm)*uu(ji ,jj,jk,Kmm), 0._wp ) - & |
---|
351 | & MIN( e2u(ji-1,jj)*e3u(ji-1,jj,jk,Kmm)*uu(ji-1,jj,jk,Kmm), 0._wp ) ) & |
---|
352 | & * r1_e1e2t(ji,jj) & |
---|
353 | & + ( MAX( e1v(ji,jj )*e3v(ji,jj ,jk,Kmm)*vv(ji,jj ,jk,Kmm), 0._wp ) - & |
---|
354 | & MIN( e1v(ji,jj-1)*e3v(ji,jj-1,jk,Kmm)*vv(ji,jj-1,jk,Kmm), 0._wp ) ) & |
---|
355 | & * r1_e1e2t(ji,jj) & |
---|
356 | & ) * z1_e3t |
---|
357 | END_3D |
---|
358 | ENDIF |
---|
359 | CALL lbc_lnk( 'sshwzv', Cu_adv, 'T', 1. ) |
---|
360 | ! |
---|
361 | CALL iom_put("Courant",Cu_adv) |
---|
362 | ! |
---|
363 | IF( MAXVAL( Cu_adv(:,:,:) ) > Cu_min ) THEN ! Quick check if any breaches anywhere |
---|
364 | DO_3DS_11_11( jpkm1, 2, -1 ) |
---|
365 | ! |
---|
366 | zCu = MAX( Cu_adv(ji,jj,jk) , Cu_adv(ji,jj,jk-1) ) |
---|
367 | ! alt: |
---|
368 | ! IF ( ww(ji,jj,jk) > 0._wp ) THEN |
---|
369 | ! zCu = Cu_adv(ji,jj,jk) |
---|
370 | ! ELSE |
---|
371 | ! zCu = Cu_adv(ji,jj,jk-1) |
---|
372 | ! ENDIF |
---|
373 | ! |
---|
374 | IF( zCu <= Cu_min ) THEN !<-- Fully explicit |
---|
375 | zcff = 0._wp |
---|
376 | ELSEIF( zCu < Cu_cut ) THEN !<-- Mixed explicit |
---|
377 | zcff = ( zCu - Cu_min )**2 |
---|
378 | zcff = zcff / ( Fcu + zcff ) |
---|
379 | ELSE !<-- Mostly implicit |
---|
380 | zcff = ( zCu - Cu_max )/ zCu |
---|
381 | ENDIF |
---|
382 | zcff = MIN(1._wp, zcff) |
---|
383 | ! |
---|
384 | wi(ji,jj,jk) = zcff * ww(ji,jj,jk) |
---|
385 | ww(ji,jj,jk) = ( 1._wp - zcff ) * ww(ji,jj,jk) |
---|
386 | ! |
---|
387 | Cu_adv(ji,jj,jk) = zcff ! Reuse array to output coefficient below and in stp_ctl |
---|
388 | END_3D |
---|
389 | Cu_adv(:,:,1) = 0._wp |
---|
390 | ELSE |
---|
391 | ! Fully explicit everywhere |
---|
392 | Cu_adv(:,:,:) = 0._wp ! Reuse array to output coefficient below and in stp_ctl |
---|
393 | wi (:,:,:) = 0._wp |
---|
394 | ENDIF |
---|
395 | CALL iom_put("wimp",wi) |
---|
396 | CALL iom_put("wi_cff",Cu_adv) |
---|
397 | CALL iom_put("wexp",ww) |
---|
398 | ! |
---|
399 | IF( ln_timing ) CALL timing_stop('wAimp') |
---|
400 | ! |
---|
401 | END SUBROUTINE wAimp |
---|
402 | !!====================================================================== |
---|
403 | END MODULE sshwzv |
---|