1 | MODULE step |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE step *** |
---|
4 | !! Time-stepping : manager of the shallow water equation time stepping |
---|
5 | !!====================================================================== |
---|
6 | !! History : NEMO ! 2020-03 (A. Nasser, G. Madec) Original code from 4.0.2 |
---|
7 | !!---------------------------------------------------------------------- |
---|
8 | #if defined key_qco |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | !! 'key_qco' EMPTY MODULE Quasi-Eulerian vertical coordonate |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | #else |
---|
13 | !!---------------------------------------------------------------------- |
---|
14 | !! stp : Shallow Water time-stepping |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | USE step_oce ! time stepping definition modules |
---|
17 | USE phycst ! physical constants |
---|
18 | USE usrdef_nam |
---|
19 | USE lib_mpp ! MPP library |
---|
20 | USE iom ! xIOs server |
---|
21 | |
---|
22 | IMPLICIT NONE |
---|
23 | PRIVATE |
---|
24 | |
---|
25 | PUBLIC stp ! called by nemogcm.F90 |
---|
26 | |
---|
27 | !!---------------------------------------------------------------------- |
---|
28 | !! time level indices |
---|
29 | !!---------------------------------------------------------------------- |
---|
30 | INTEGER, PUBLIC :: Nbb, Nnn, Naa, Nrhs !! used by nemo_init |
---|
31 | |
---|
32 | !! * Substitutions |
---|
33 | # include "do_loop_substitute.h90" |
---|
34 | !!---------------------------------------------------------------------- |
---|
35 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
36 | !! $Id$ |
---|
37 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
38 | !!---------------------------------------------------------------------- |
---|
39 | CONTAINS |
---|
40 | |
---|
41 | #if defined key_agrif |
---|
42 | RECURSIVE SUBROUTINE stp( ) |
---|
43 | INTEGER :: kstp ! ocean time-step index |
---|
44 | #else |
---|
45 | SUBROUTINE stp( kstp ) |
---|
46 | INTEGER, INTENT(in) :: kstp ! ocean time-step index |
---|
47 | #endif |
---|
48 | !!---------------------------------------------------------------------- |
---|
49 | !! *** ROUTINE stp *** |
---|
50 | !! |
---|
51 | !! ** Purpose : - Time stepping of shallow water (SHW) (momentum and ssh eqs.) |
---|
52 | !! |
---|
53 | !! ** Method : -1- Update forcings |
---|
54 | !! -2- Update the ssh at Naa |
---|
55 | !! -3- Compute the momentum trends (Nrhs) |
---|
56 | !! -4- Update the horizontal velocity |
---|
57 | !! -5- Apply Asselin time filter to uu,vv,ssh |
---|
58 | !! -6- Outputs and diagnostics |
---|
59 | !!---------------------------------------------------------------------- |
---|
60 | INTEGER :: ji, jj, jk ! dummy loop indice |
---|
61 | INTEGER :: indic ! error indicator if < 0 |
---|
62 | !!gm kcall can be removed, I guess |
---|
63 | INTEGER :: kcall ! optional integer argument (dom_vvl_sf_nxt) |
---|
64 | REAL(wp):: z1_2rho0 ! local scalars |
---|
65 | |
---|
66 | REAL(wp) :: zue3a, zue3n, zue3b ! local scalars |
---|
67 | REAL(wp) :: zve3a, zve3n, zve3b ! - - |
---|
68 | REAL(wp) :: ze3t_tf, ze3u_tf, ze3v_tf, zua, zva |
---|
69 | REAL(wp), DIMENSION(jpi,jpj) :: zssh |
---|
70 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3u, ze3v |
---|
71 | !! --------------------------------------------------------------------- |
---|
72 | #if defined key_agrif |
---|
73 | kstp = nit000 + Agrif_Nb_Step() |
---|
74 | Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs ! agrif_oce module copies of time level indices |
---|
75 | IF( lk_agrif_debug ) THEN |
---|
76 | IF( Agrif_Root() .and. lwp) WRITE(*,*) '---' |
---|
77 | IF(lwp) WRITE(*,*) 'Grid Number', Agrif_Fixed(),' time step ', kstp, 'int tstep', Agrif_NbStepint() |
---|
78 | ENDIF |
---|
79 | IF( kstp == nit000 + 1 ) lk_agrif_fstep = .FALSE. |
---|
80 | # if defined key_iomput |
---|
81 | IF( Agrif_Nbstepint() == 0 ) CALL iom_swap( cxios_context ) |
---|
82 | # endif |
---|
83 | #endif |
---|
84 | ! |
---|
85 | IF( ln_timing ) CALL timing_start('stp') |
---|
86 | ! |
---|
87 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
88 | ! model timestep |
---|
89 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
90 | ! |
---|
91 | IF( l_1st_euler ) THEN ! start or restart with Euler 1st time-step |
---|
92 | rDt = rn_Dt |
---|
93 | r1_Dt = 1._wp / rDt |
---|
94 | ENDIF |
---|
95 | |
---|
96 | IF ( kstp == nit000 ) ww(:,:,:) = 0._wp ! initialize vertical velocity one for all to zero |
---|
97 | |
---|
98 | ! |
---|
99 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
100 | ! update I/O and calendar |
---|
101 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
102 | indic = 0 ! reset to no error condition |
---|
103 | |
---|
104 | IF( kstp == nit000 ) THEN ! initialize IOM context (must be done after nemo_init for AGRIF+XIOS+OASIS) |
---|
105 | CALL iom_init( cxios_context, ld_closedef=.FALSE. ) ! for model grid (including passible AGRIF zoom) |
---|
106 | IF( lk_diamlr ) CALL dia_mlr_iom_init ! with additional setup for multiple-linear-regression analysis |
---|
107 | CALL iom_init_closedef |
---|
108 | IF( ln_crs ) CALL iom_init( TRIM(cxios_context)//"_crs" ) ! for coarse grid |
---|
109 | ENDIF |
---|
110 | IF( kstp /= nit000 ) CALL day( kstp ) ! Calendar (day was already called at nit000 in day_init) |
---|
111 | CALL iom_setkt( kstp - nit000 + 1, cxios_context ) ! tell IOM we are at time step kstp |
---|
112 | IF( ln_crs ) CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" ) ! tell IOM we are at time step kstp |
---|
113 | |
---|
114 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
115 | ! Update external forcing (tides, open boundaries, ice shelf interaction and surface boundary condition (including sea-ice) |
---|
116 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
117 | IF( ln_tide ) CALL tide_update( kstp ) ! update tide potential |
---|
118 | IF( ln_apr_dyn ) CALL sbc_apr ( kstp ) ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib) |
---|
119 | IF( ln_bdy ) CALL bdy_dta ( kstp, Nnn ) ! update dynamic & tracer data at open boundaries |
---|
120 | CALL sbc ( kstp, Nbb, Nnn ) ! Sea Boundary Condition |
---|
121 | |
---|
122 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
123 | ! Ocean physics update |
---|
124 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
125 | ! LATERAL PHYSICS |
---|
126 | ! ! eddy diffusivity coeff. |
---|
127 | IF( l_ldfdyn_time ) CALL ldf_dyn( kstp, Nbb ) ! eddy viscosity coeff. |
---|
128 | |
---|
129 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
130 | ! Ocean dynamics : hdiv, ssh, e3, u, v, w |
---|
131 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
132 | |
---|
133 | CALL ssh_nxt ( kstp, Nbb, Nnn, ssh, Naa ) ! after ssh (includes call to div_hor) |
---|
134 | uu(:,:,:,Nrhs) = 0._wp ! set dynamics trends to zero |
---|
135 | vv(:,:,:,Nrhs) = 0._wp |
---|
136 | |
---|
137 | IF( .NOT.ln_linssh ) CALL dom_vvl_sf_nxt( kstp, Nbb, Nnn, Naa ) ! after vertical scale factors |
---|
138 | |
---|
139 | IF( ln_bdy ) CALL bdy_dyn3d_dmp ( kstp, Nbb, uu, vv, Nrhs ) ! bdy damping trends |
---|
140 | |
---|
141 | #if defined key_agrif |
---|
142 | IF(.NOT. Agrif_Root()) & |
---|
143 | & CALL Agrif_Sponge_dyn ! momentum sponge |
---|
144 | #endif |
---|
145 | |
---|
146 | !!an - calcul du gradient de pression horizontal (explicit) |
---|
147 | zssh(:,:) = ssh(:,:,Nnn) |
---|
148 | # if defined key_bvp |
---|
149 | ! gradient and divergence not penalised |
---|
150 | zssh(:,:) = zssh(:,:) * r1_rpo(:,:,1) |
---|
151 | #endif |
---|
152 | DO_3D_00_00( 1, jpkm1 ) |
---|
153 | uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) - grav * ( zssh(ji+1,jj) - zssh(ji,jj) ) * r1_e1u(ji,jj) |
---|
154 | vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) - grav * ( zssh(ji,jj+1) - zssh(ji,jj) ) * r1_e2v(ji,jj) |
---|
155 | END_3D |
---|
156 | ! |
---|
157 | |
---|
158 | ! IF( kstp == nit000 .AND. lwp ) THEN |
---|
159 | ! WRITE(numout,*) |
---|
160 | ! WRITE(numout,*) 'step.F90 : classic script used' |
---|
161 | ! WRITE(numout,*) '~~~~~~~' |
---|
162 | ! IF( ln_dynvor_ens_adVO .OR. ln_dynvor_ens_adKE .OR. ln_dynvor_ens_adKEVO & |
---|
163 | ! & .OR. ln_dynvor_ene_adVO .OR. ln_dynvor_ene_adKE .OR. ln_dynvor_ene_adKEVO ) THEN |
---|
164 | ! CALL ctl_stop('STOP','step : alternative direction asked but classis step' ) |
---|
165 | ! ENDIF |
---|
166 | ! ENDIF |
---|
167 | !!an |
---|
168 | ! CALL dyn_adv( kstp, Nbb, Nnn , uu, vv, Nrhs ) ! advection (VF or FF) ==> RHS |
---|
169 | ! |
---|
170 | ! CALL dyn_vor( kstp, Nbb, Nnn , uu, vv, Nrhs ) ! vorticity ==> RHS |
---|
171 | ! |
---|
172 | !!an In dynvor, dynkegAD is called even if not AD, so we keep the same step.F90 |
---|
173 | |
---|
174 | CALL dyn_vor( kstp, Nbb, Nnn , uu, vv, Nrhs) ! vorticity ==> RHS |
---|
175 | |
---|
176 | CALL dyn_ldf( kstp, Nbb, Nnn , uu, vv, Nrhs ) ! lateral mixing |
---|
177 | |
---|
178 | ! add wind stress forcing and layer linear friction to the RHS |
---|
179 | ze3u(:,:,:) = e3u(:,:,:,Nnn) |
---|
180 | ze3v(:,:,:) = e3v(:,:,:,Nnn) |
---|
181 | # if defined key_bvp |
---|
182 | !ze3u(:,:,:) = ze3u(:,:,:) * r1_rpo(:,:,:) |
---|
183 | !ze3v(:,:,:) = ze3v(:,:,:) * r1_rpo(:,:,:) |
---|
184 | #endif |
---|
185 | z1_2rho0 = 0.5_wp * r1_rho0 |
---|
186 | DO_3D_00_00(1,jpkm1) |
---|
187 | uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + z1_2rho0 * ( utau_b(ji,jj) + utau(ji,jj) ) / ze3u(ji,jj,jk) & |
---|
188 | & - rn_rfr * uu(ji,jj,jk,Nbb) |
---|
189 | vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + z1_2rho0 * ( vtau_b(ji,jj) + vtau(ji,jj) ) / ze3v(ji,jj,jk) & |
---|
190 | & - rn_rfr * vv(ji,jj,jk,Nbb) |
---|
191 | END_3D |
---|
192 | ! |
---|
193 | # if defined key_bvp |
---|
194 | ! Add frictionnal term - sigma * u |
---|
195 | ! can be done via sigmaT (mean_ij) |
---|
196 | DO_3D_00_00(1,jpkm1) |
---|
197 | uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) - bmpu(ji,jj,jk) * uu(ji,jj,jk,Nbb) |
---|
198 | vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) - bmpv(ji,jj,jk) * vv(ji,jj,jk,Nbb) |
---|
199 | END_3D |
---|
200 | #endif |
---|
201 | ! |
---|
202 | |
---|
203 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
204 | ! Leap-Frog time splitting + Robert-Asselin time filter on u,v,e3 |
---|
205 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
206 | |
---|
207 | !! what about IF( .NOT.ln_linssh ) ? |
---|
208 | !!an futur module dyn_nxt (a la place de dyn_atf) |
---|
209 | |
---|
210 | IF( ln_dynadv_vec ) THEN ! vector invariant form : applied on velocity |
---|
211 | IF( l_1st_euler ) THEN ! Euler time stepping (no Asselin filter) |
---|
212 | DO_3D_00_00(1,jpkm1) |
---|
213 | uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
214 | vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
215 | END_3D |
---|
216 | ELSE ! Leap Frog time stepping + Asselin filter |
---|
217 | DO_3D_11_11(1,jpkm1) |
---|
218 | zua = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
219 | zva = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
220 | ! ! Asselin time filter on u,v (Nnn) |
---|
221 | uu(ji,jj,jk,Nnn) = uu(ji,jj,jk,Nnn) + rn_atfp * (uu(ji,jj,jk,Nbb) - 2._wp * uu(ji,jj,jk,Nnn) + zua) |
---|
222 | vv(ji,jj,jk,Nnn) = vv(ji,jj,jk,Nnn) + rn_atfp * (vv(ji,jj,jk,Nbb) - 2._wp * vv(ji,jj,jk,Nnn) + zva) |
---|
223 | ! |
---|
224 | ze3u_tf = e3u(ji,jj,jk,Nnn) + rn_atfp * ( e3u(ji,jj,jk,Nbb) - 2._wp * e3u(ji,jj,jk,Nnn) + e3u(ji,jj,jk,Naa) ) |
---|
225 | ze3v_tf = e3v(ji,jj,jk,Nnn) + rn_atfp * ( e3v(ji,jj,jk,Nbb) - 2._wp * e3v(ji,jj,jk,Nnn) + e3v(ji,jj,jk,Naa) ) |
---|
226 | ze3t_tf = e3t(ji,jj,jk,Nnn) + rn_atfp * ( e3t(ji,jj,jk,Nbb) - 2._wp * e3t(ji,jj,jk,Nnn) + e3t(ji,jj,jk,Naa) ) |
---|
227 | ! |
---|
228 | e3u(ji,jj,jk,Nnn) = ze3u_tf |
---|
229 | e3v(ji,jj,jk,Nnn) = ze3v_tf |
---|
230 | e3t(ji,jj,jk,Nnn) = ze3t_tf |
---|
231 | ! |
---|
232 | uu(ji,jj,jk,Naa) = zua |
---|
233 | vv(ji,jj,jk,Naa) = zva |
---|
234 | END_3D |
---|
235 | ENDIF |
---|
236 | ! |
---|
237 | ELSE ! flux form : applied on thickness weighted velocity |
---|
238 | IF( l_1st_euler ) THEN ! Euler time stepping (no Asselin filter) |
---|
239 | DO_3D_00_00(1,jpkm1) |
---|
240 | zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb) |
---|
241 | zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb) |
---|
242 | ! ! LF time stepping |
---|
243 | zue3a = zue3b + rDt * e3u(ji,jj,jk,Nrhs) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
244 | zve3a = zve3b + rDt * e3v(ji,jj,jk,Nrhs) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
245 | ! |
---|
246 | uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa) |
---|
247 | vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa) |
---|
248 | END_3D |
---|
249 | ELSE ! Leap Frog time stepping + Asselin filter |
---|
250 | DO_3D_11_11(1,jpkm1) |
---|
251 | zue3n = e3u(ji,jj,jk,Nnn) * uu(ji,jj,jk,Nnn) |
---|
252 | zve3n = e3v(ji,jj,jk,Nnn) * vv(ji,jj,jk,Nnn) |
---|
253 | zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb) |
---|
254 | zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb) |
---|
255 | ! ! LF time stepping |
---|
256 | zue3a = zue3b + rDt * e3u(ji,jj,jk,Nrhs) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
257 | zve3a = zve3b + rDt * e3v(ji,jj,jk,Nrhs) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
258 | ! ! Asselin time filter on e3u/v/t |
---|
259 | ze3u_tf = e3u(ji,jj,jk,Nnn) + rn_atfp * ( e3u(ji,jj,jk,Nbb) - 2._wp * e3u(ji,jj,jk,Nnn) + e3u(ji,jj,jk,Naa) ) |
---|
260 | ze3v_tf = e3v(ji,jj,jk,Nnn) + rn_atfp * ( e3v(ji,jj,jk,Nbb) - 2._wp * e3v(ji,jj,jk,Nnn) + e3v(ji,jj,jk,Naa) ) |
---|
261 | ze3t_tf = e3t(ji,jj,jk,Nnn) + rn_atfp * ( e3t(ji,jj,jk,Nbb) - 2._wp * e3t(ji,jj,jk,Nnn) + e3t(ji,jj,jk,Naa) ) |
---|
262 | ! ! Asselin time filter on u,v (Nnn) |
---|
263 | uu(ji,jj,jk,Nnn) = ( zue3n + rn_atfp * ( zue3b - 2._wp * zue3n + zue3a ) ) / ze3u_tf |
---|
264 | vv(ji,jj,jk,Nnn) = ( zve3n + rn_atfp * ( zve3b - 2._wp * zve3n + zve3a ) ) / ze3v_tf |
---|
265 | ! |
---|
266 | e3u(ji,jj,jk,Nnn) = ze3u_tf |
---|
267 | e3v(ji,jj,jk,Nnn) = ze3v_tf |
---|
268 | e3t(ji,jj,jk,Nnn) = ze3t_tf |
---|
269 | ! |
---|
270 | uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa) |
---|
271 | vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa) |
---|
272 | END_3D |
---|
273 | ENDIF |
---|
274 | ENDIF |
---|
275 | |
---|
276 | |
---|
277 | CALL lbc_lnk_multi( 'stp', uu(:,:,:,Nnn), 'U', -1., vv(:,:,:,Nnn), 'V', -1., & !* local domain boundaries |
---|
278 | & uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. ) |
---|
279 | |
---|
280 | !!an |
---|
281 | |
---|
282 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
283 | ! Set boundary conditions, time filter and swap time levels |
---|
284 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
285 | |
---|
286 | !!an TO BE ADDED : dyn_nxt |
---|
287 | !! CALL dyn_atf ( kstp, Nbb, Nnn, Naa, uu, vv, e3t, e3u, e3v ) ! time filtering of "now" velocities and scale factors |
---|
288 | !!an TO BE ADDED : a simplifier |
---|
289 | !! CALL ssh_atf ( kstp, Nbb, Nnn, Naa, ssh ) ! time filtering of "now" sea surface height |
---|
290 | IF ( .NOT.( l_1st_euler ) ) THEN ! Only do time filtering for leapfrog timesteps |
---|
291 | ! ! filtering "now" field |
---|
292 | ssh(:,:,Nnn) = ssh(:,:,Nnn) + rn_atfp * ( ssh(:,:,Nbb) - 2 * ssh(:,:,Nnn) + ssh(:,:,Naa) ) |
---|
293 | ENDIF |
---|
294 | !!an |
---|
295 | |
---|
296 | |
---|
297 | ! Swap time levels |
---|
298 | Nrhs = Nbb |
---|
299 | Nbb = Nnn |
---|
300 | Nnn = Naa |
---|
301 | Naa = Nrhs |
---|
302 | ! |
---|
303 | CALL dom_vvl_sf_update( kstp, Nbb, Nnn, Naa ) ! recompute vertical scale factors |
---|
304 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
305 | ! diagnostics and outputs |
---|
306 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
307 | IF( ln_floats ) CALL flo_stp ( kstp, Nbb, Nnn ) ! drifting Floats |
---|
308 | IF( ln_diacfl ) CALL dia_cfl ( kstp, Nnn ) ! Courant number diagnostics |
---|
309 | |
---|
310 | CALL dia_wri ( kstp, Nnn ) ! ocean model: outputs |
---|
311 | ! |
---|
312 | IF( lrst_oce ) CALL rst_write ( kstp, Nbb, Nnn ) ! write output ocean restart file |
---|
313 | |
---|
314 | |
---|
315 | #if defined key_agrif |
---|
316 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
317 | ! AGRIF |
---|
318 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
319 | Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs ! agrif_oce module copies of time level indices |
---|
320 | CALL Agrif_Integrate_ChildGrids( stp ) ! allows to finish all the Child Grids before updating |
---|
321 | |
---|
322 | IF( Agrif_NbStepint() == 0 ) THEN |
---|
323 | CALL Agrif_update_all( ) ! Update all components |
---|
324 | ENDIF |
---|
325 | #endif |
---|
326 | IF( ln_diaobs ) CALL dia_obs ( kstp, Nnn ) ! obs-minus-model (assimilation) diagnostics (call after dynamics update) |
---|
327 | |
---|
328 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
329 | ! Control |
---|
330 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
331 | CALL stp_ctl ( kstp, Nbb, Nnn, indic ) |
---|
332 | |
---|
333 | |
---|
334 | IF( kstp == nit000 ) THEN ! 1st time step only |
---|
335 | CALL iom_close( numror ) ! close input ocean restart file |
---|
336 | IF(lwm) CALL FLUSH ( numond ) ! flush output namelist oce |
---|
337 | IF(lwm .AND. numoni /= -1 ) CALL FLUSH ( numoni ) ! flush output namelist ice (if exist) |
---|
338 | ENDIF |
---|
339 | |
---|
340 | ! |
---|
341 | #if defined key_iomput |
---|
342 | IF( kstp == nitend .OR. indic < 0 ) THEN |
---|
343 | CALL iom_context_finalize( cxios_context ) ! needed for XIOS+AGRIF |
---|
344 | IF(lrxios) CALL iom_context_finalize( crxios_context ) |
---|
345 | ENDIF |
---|
346 | #endif |
---|
347 | ! |
---|
348 | IF( l_1st_euler ) THEN ! recover Leap-frog timestep |
---|
349 | rDt = 2._wp * rn_Dt |
---|
350 | r1_Dt = 1._wp / rDt |
---|
351 | l_1st_euler = .FALSE. |
---|
352 | ENDIF |
---|
353 | ! |
---|
354 | IF( ln_timing ) CALL timing_stop('stp') |
---|
355 | ! |
---|
356 | END SUBROUTINE stp |
---|
357 | #endif |
---|
358 | ! |
---|
359 | !!====================================================================== |
---|
360 | END MODULE step |
---|