New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
stepLF.F90 in NEMO/branches/2020/dev_r12527_Gurvan_ShallowWater/cfgs/penzAM98/MY_SRC – NEMO

source: NEMO/branches/2020/dev_r12527_Gurvan_ShallowWater/cfgs/penzAM98/MY_SRC/stepLF.F90 @ 13562

Last change on this file since 13562 was 13562, checked in by gm, 4 years ago

zgr_pse created only for NS coast

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