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.
step.F90 in branches/2014/dev_r4650_UKMO12_CFL_diags_take2/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/2014/dev_r4650_UKMO12_CFL_diags_take2/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 4714

Last change on this file since 4714 was 4714, checked in by timgraham, 10 years ago

Added diacfl module and added a nn_diacfl to namctl to control it

  • Property svn:keywords set to Id
File size: 20.9 KB
Line 
1MODULE step
2   !!======================================================================
3   !!                       ***  MODULE step  ***
4   !! Time-stepping    : manager of the ocean, tracer and ice time stepping
5   !!======================================================================
6   !! History :  OPA  !  1991-03  (G. Madec)  Original code
7   !!             -   !  1991-11  (G. Madec)
8   !!             -   !  1992-06  (M. Imbard)  add a first output record
9   !!             -   !  1996-04  (G. Madec)  introduction of dynspg
10   !!             -   !  1996-04  (M.A. Foujols)  introduction of passive tracer
11   !!            8.0  !  1997-06  (G. Madec)  new architecture of call
12   !!            8.2  !  1997-06  (G. Madec, M. Imbard, G. Roullet)  free surface
13   !!             -   !  1999-02  (G. Madec, N. Grima)  hpg implicit
14   !!             -   !  2000-07  (J-M Molines, M. Imbard)  Open Bondary Conditions
15   !!   NEMO     1.0  !  2002-06  (G. Madec)  free form, suppress macro-tasking
16   !!             -   !  2004-08  (C. Talandier) New trends organization
17   !!             -   !  2005-01  (C. Ethe) Add the KPP closure scheme
18   !!             -   !  2005-11  (G. Madec)  Reorganisation of tra and dyn calls
19   !!             -   !  2006-01  (L. Debreu, C. Mazauric)  Agrif implementation
20   !!             -   !  2006-07  (S. Masson)  restart using iom
21   !!            3.2  !  2009-02  (G. Madec, R. Benshila)  reintroduicing z*-coordinate
22   !!             -   !  2009-06  (S. Masson, G. Madec)  TKE restart compatible with key_cpl
23   !!            3.3  !  2010-05  (K. Mogensen, A. Weaver, M. Martin, D. Lea) Assimilation interface
24   !!             -   !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA
25   !!            3.4  !  2011-04  (G. Madec, C. Ethe) Merge of dtatem and dtasal
26   !!                 !  2012-07  (J. Simeon, G. Madec. C. Ethe) Online coarsening of outputs
27   !!----------------------------------------------------------------------
28
29   !!----------------------------------------------------------------------
30   !!   stp             : OPA system time-stepping
31   !!----------------------------------------------------------------------
32   USE step_oce         ! time stepping definition modules
33
34   IMPLICIT NONE
35   PRIVATE
36
37   PUBLIC   stp   ! called by opa.F90
38
39   !! * Substitutions
40#  include "domzgr_substitute.h90"
41#  include "zdfddm_substitute.h90"
42   !!----------------------------------------------------------------------
43   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
44   !! $Id$
45   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
46   !!----------------------------------------------------------------------
47CONTAINS
48
49#if defined key_agrif
50   SUBROUTINE stp( )
51      INTEGER             ::   kstp   ! ocean time-step index
52#else
53   SUBROUTINE stp( kstp )
54      INTEGER, INTENT(in) ::   kstp   ! ocean time-step index
55#endif
56      !!----------------------------------------------------------------------
57      !!                     ***  ROUTINE stp  ***
58      !!
59      !! ** Purpose : - Time stepping of OPA (momentum and active tracer eqs.)
60      !!              - Time stepping of LIM (dynamic and thermodynamic eqs.)
61      !!              - Tme stepping  of TRC (passive tracer eqs.)
62      !!
63      !! ** Method  : -1- Update forcings and data
64      !!              -2- Update ocean physics
65      !!              -3- Compute the t and s trends
66      !!              -4- Update t and s
67      !!              -5- Compute the momentum trends
68      !!              -6- Update the horizontal velocity
69      !!              -7- Compute the diagnostics variables (rd,N2, div,cur,w)
70      !!              -8- Outputs and diagnostics
71      !!----------------------------------------------------------------------
72      INTEGER ::   jk       ! dummy loop indice
73      INTEGER ::   indic    ! error indicator if < 0
74      INTEGER ::   kcall    ! optional integer argument (dom_vvl_sf_nxt)
75      !! ---------------------------------------------------------------------
76
77#if defined key_agrif
78      kstp = nit000 + Agrif_Nb_Step()
79!      IF ( Agrif_Root() .and. lwp) Write(*,*) '---'
80!      IF (lwp) Write(*,*) 'Grid Number',Agrif_Fixed(),' time step ',kstp
81      IF ( kstp == (nit000 + 1) ) lk_agrif_fstep = .FALSE.
82# if defined key_iomput
83      IF( Agrif_Nbstepint() == 0 )   CALL iom_swap( "nemo" )
84# endif
85#endif
86                             indic = 0           ! reset to no error condition
87      IF( kstp == nit000 ) THEN
88                      CALL iom_init( "nemo" )      ! iom_put initialization (must be done after nemo_init for AGRIF+XIOS+OASIS)
89         IF( ln_crs ) CALL iom_init( "nemo_crs" )  ! initialize context for coarse grid
90      ENDIF
91
92      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
93                             CALL iom_setkt( kstp - nit000 + 1, "nemo"     )   ! say to iom that we are at time step kstp
94      IF( ln_crs     )       CALL iom_setkt( kstp - nit000 + 1, "nemo_crs" )   ! say to iom that we are at time step kstp
95
96      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
97      ! Update data, open boundaries, surface boundary condition (including sea-ice)
98      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
99      IF( lk_tide    )   CALL sbc_tide( kstp )
100      IF( lk_bdy     )   CALL bdy_dta ( kstp, time_offset=+1 )   ! update dynamic & tracer data at open boundaries
101
102                         CALL sbc    ( kstp )         ! Sea Boundary Condition (including sea-ice)
103                                                      ! clem: moved here for bdy ice purpose
104
105      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
106      ! Ocean physics update                (ua, va, tsa used as workspace)
107      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
108                         CALL bn2( tsb, rn2b )        ! before Brunt-Vaisala frequency
109                         CALL bn2( tsn, rn2  )        ! now    Brunt-Vaisala frequency
110      !
111      !  VERTICAL PHYSICS
112                         CALL zdf_bfr( kstp )         ! bottom friction (if quadratic)
113      !                                               ! Vertical eddy viscosity and diffusivity coefficients
114      IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz
115      IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz
116      IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz
117      IF( lk_zdfkpp  )   CALL zdf_kpp( kstp )            ! KPP closure scheme for Kz
118      IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value)
119         avt (:,:,:) = rn_avt0 * tmask(:,:,:)
120         avmu(:,:,:) = rn_avm0 * umask(:,:,:)
121         avmv(:,:,:) = rn_avm0 * vmask(:,:,:)
122      ENDIF
123      IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths
124         DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:) * tmask(:,:,jk)   ;   END DO
125      ENDIF
126      IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity
127
128      IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing
129
130      IF( lk_zdfddm .AND. .NOT. lk_zdfkpp )   &
131         &               CALL zdf_ddm( kstp )         ! double diffusive mixing
132
133                         CALL zdf_mxl( kstp )         ! mixed layer depth
134
135                                                      ! write TKE or GLS information in the restart file
136      IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' )
137      IF( lrst_oce .AND. lk_zdfgls )   CALL gls_rst( kstp, 'WRITE' )
138      !
139      !  LATERAL  PHYSICS
140      !
141      IF( lk_ldfslp ) THEN                            ! slope of lateral mixing
142                         CALL eos( tsb, rhd, gdept_0(:,:,:) )             ! before in situ density
143         IF( ln_zps )    CALL zps_hde( kstp, jpts, tsb, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
144            &                                      rhd, gru , grv  )      ! of t, s, rd at the last ocean level
145         IF( ln_traldf_grif ) THEN                           ! before slope for Griffies operator
146                         CALL ldf_slp_grif( kstp )
147         ELSE
148                         CALL ldf_slp( kstp, rhd, rn2b )     ! before slope for Madec operator
149         ENDIF
150      ENDIF
151#if defined key_traldf_c2d
152      IF( lk_traldf_eiv )   CALL ldf_eiv( kstp )      ! eddy induced velocity coefficient
153#endif
154#if defined key_traldf_c3d && key_traldf_smag
155                          CALL ldf_tra_smag( kstp )      ! eddy induced velocity coefficient
156#  endif
157#if defined key_dynldf_c3d && key_dynldf_smag
158                          CALL ldf_dyn_smag( kstp )      ! eddy induced velocity coefficient
159#  endif
160
161      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
162      !  Ocean dynamics : hdiv, rot, ssh, e3, wn
163      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
164                         CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_cur)
165      IF( lk_vvl     )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors
166                         CALL wzv           ( kstp )  ! now cross-level velocity
167
168      IF( lk_dynspg_ts ) THEN 
169          ! In case the time splitting case, update almost all momentum trends here:
170          ! Note that the computation of vertical velocity above, hence "after" sea level
171          ! is necessary to compute momentum advection for the rhs of barotropic loop:
172                                  CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) ) ! now in situ density for hpg computation
173          IF( ln_zps      )       CALL zps_hde( kstp, jpts, tsn, gtsu, gtsv,  &   ! zps: now hor. derivative
174                &                                          rhd, gru , grv  )      ! of t, s, rd at the last ocean level
175
176                                  ua(:,:,:) = 0.e0             ! set dynamics trends to zero
177                                  va(:,:,:) = 0.e0
178          IF(  ln_asmiau .AND. &
179             & ln_dyninc       )  CALL dyn_asm_inc  ( kstp )   ! apply dynamics assimilation increment
180          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! subtract Neptune velocities (simplified)
181          IF( lk_bdy           )  CALL bdy_dyn3d_dmp( kstp )   ! bdy damping trends
182                                  CALL dyn_adv      ( kstp )   ! advection (vector or flux form)
183                                  CALL dyn_vor      ( kstp )   ! vorticity term including Coriolis
184                                  CALL dyn_ldf      ( kstp )   ! lateral mixing
185          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! add Neptune velocities (simplified)
186#if defined key_agrif
187          IF(.NOT. Agrif_Root())  CALL Agrif_Sponge_dyn        ! momentum sponge
188#endif
189                                  CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
190                                  CALL dyn_spg( kstp, indic )  ! surface pressure gradient
191
192                                  ua_sv(:,:,:) = ua(:,:,:)     ! Save trends (barotropic trend has been fully updated at this stage)
193                                  va_sv(:,:,:) = va(:,:,:)
194
195                                  CALL div_cur( kstp )         ! Horizontal divergence & Relative vorticity (2nd call in time-split case)
196          IF( lk_vvl     )        CALL dom_vvl_sf_nxt( kstp, kcall=2 )  ! after vertical scale factors (update depth average component)
197                                  CALL wzv           ( kstp )  ! now cross-level velocity
198      ENDIF
199
200      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
201      ! diagnostics and outputs             (ua, va, tsa used as workspace)
202      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
203      IF( lk_floats  )   CALL flo_stp( kstp )         ! drifting Floats
204      IF( nn_diacfl == 1 )   CALL dia_cfl( kstp )         ! Courant number diagnostics
205      IF( lk_diahth  )   CALL dia_hth( kstp )         ! Thermocline depth (20 degres isotherm depth)
206      IF( lk_diafwb  )   CALL dia_fwb( kstp )         ! Fresh water budget diagnostics
207      IF( ln_diaptr  )   CALL dia_ptr( kstp )         ! Poleward TRansports diagnostics
208      IF( lk_diadct  )   CALL dia_dct( kstp )         ! Transports
209      IF( lk_diaar5  )   CALL dia_ar5( kstp )         ! ar5 diag
210      IF( lk_diaharm )   CALL dia_harm( kstp )        ! Tidal harmonic analysis
211                         CALL dia_wri( kstp )         ! ocean model: outputs
212      !
213      IF( ln_crs     )   CALL crs_fld( kstp )         ! ocean model: online field coarsening & output
214
215
216#if defined key_top
217      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
218      ! Passive Tracer Model
219      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
220                         CALL trc_stp( kstp )         ! time-stepping
221#endif
222
223      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
224      ! Active tracers                              (ua, va used as workspace)
225      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
226                             tsa(:,:,:,:) = 0.e0            ! set tracer trends to zero
227
228      IF(  ln_asmiau .AND. &
229         & ln_trainc     )   CALL tra_asm_inc( kstp )       ! apply tracer assimilation increment
230                             CALL tra_sbc    ( kstp )       ! surface boundary condition
231      IF( ln_traqsr      )   CALL tra_qsr    ( kstp )       ! penetrative solar radiation qsr
232      IF( ln_trabbc      )   CALL tra_bbc    ( kstp )       ! bottom heat flux
233      IF( lk_trabbl      )   CALL tra_bbl    ( kstp )       ! advective (and/or diffusive) bottom boundary layer scheme
234      IF( ln_tradmp      )   CALL tra_dmp    ( kstp )       ! internal damping trends
235      IF( lk_bdy         )   CALL bdy_tra_dmp( kstp )       ! bdy damping trends
236                             CALL tra_adv    ( kstp )       ! horizontal & vertical advection
237      IF( lk_zdfkpp      )   CALL tra_kpp    ( kstp )       ! KPP non-local tracer fluxes
238                             CALL tra_ldf    ( kstp )       ! lateral mixing
239#if defined key_agrif
240      IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_tra          ! tracers sponge
241#endif
242                             CALL tra_zdf    ( kstp )       ! vertical mixing and after tracer fields
243
244      IF( ln_dynhpg_imp  ) THEN                             ! semi-implicit hpg (time stepping then eos)
245         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
246                             CALL tra_nxt( kstp )                ! tracer fields at next time step
247                             CALL eos    ( tsa, rhd, rhop, fsdept_n(:,:,:) )  ! Time-filtered in situ density for hpg computation
248         IF( ln_zps      )   CALL zps_hde( kstp, jpts, tsa, gtsu, gtsv,  &    ! zps: time filtered hor. derivative
249            &                                          rhd, gru , grv  )      ! of t, s, rd at the last ocean level
250
251      ELSE                                                  ! centered hpg  (eos then time stepping)
252         IF ( .NOT. lk_dynspg_ts ) THEN                     ! eos already called in time-split case
253                                CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) )  ! now in situ density for hpg computation
254            IF( ln_zps      )   CALL zps_hde( kstp, jpts, tsn, gtsu, gtsv,  &    ! zps: now hor. derivative
255            &                                          rhd, gru , grv  )      ! of t, s, rd at the last ocean level
256         ENDIF
257         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
258                             CALL tra_nxt( kstp )                ! tracer fields at next time step
259      ENDIF
260
261      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
262      ! Dynamics                                    (tsa used as workspace)
263      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
264      IF( lk_dynspg_ts   )  THEN
265                                                             ! revert to previously computed momentum tendencies
266                                                             ! (not using ua, va as temporary arrays during tracers' update could avoid that)
267                               ua(:,:,:) = ua_sv(:,:,:)
268                               va(:,:,:) = va_sv(:,:,:)
269                                                             ! Revert now divergence and rotational to previously computed ones
270                                                             !(needed because of the time swap in div_cur, at the beginning of each time step)
271                               hdivn(:,:,:) = hdivb(:,:,:)
272                               rotn(:,:,:)  = rotb(:,:,:) 
273
274                               CALL dyn_bfr( kstp )         ! bottom friction
275                               CALL dyn_zdf( kstp )         ! vertical diffusion
276      ELSE
277                               ua(:,:,:) = 0.e0             ! set dynamics trends to zero
278                               va(:,:,:) = 0.e0
279
280        IF(  ln_asmiau .AND. &
281           & ln_dyninc      )  CALL dyn_asm_inc( kstp )     ! apply dynamics assimilation increment
282        IF( ln_bkgwri )        CALL asm_bkg_wri( kstp )     ! output background fields
283        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! subtract Neptune velocities (simplified)
284        IF( lk_bdy          )  CALL bdy_dyn3d_dmp(kstp )    ! bdy damping trends
285                               CALL dyn_adv( kstp )         ! advection (vector or flux form)
286                               CALL dyn_vor( kstp )         ! vorticity term including Coriolis
287                               CALL dyn_ldf( kstp )         ! lateral mixing
288        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! add Neptune velocities (simplified)
289#if defined key_agrif
290        IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_dyn        ! momemtum sponge
291#endif
292                               CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
293                               CALL dyn_bfr( kstp )         ! bottom friction
294                               CALL dyn_zdf( kstp )         ! vertical diffusion
295                               CALL dyn_spg( kstp, indic )  ! surface pressure gradient
296      ENDIF
297                               CALL dyn_nxt( kstp )         ! lateral velocity at next time step
298
299                               CALL ssh_swp( kstp )         ! swap of sea surface height
300      IF( lk_vvl           )   CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors
301
302      IF( ln_diahsb        )   CALL dia_hsb( kstp )         ! - ML - global conservation diagnostics
303      IF( lk_diaobs  )         CALL dia_obs( kstp )         ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
304
305      IF( lrst_oce .AND. ln_diahsb )   CALL dia_hsb_rst( kstp, 'WRITE' )
306      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
307      ! Control and restarts
308      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
309                               CALL stp_ctl( kstp, indic )
310      IF( indic < 0        )   THEN
311                               CALL ctl_stop( 'step: indic < 0' )
312                               CALL dia_wri_state( 'output.abort', kstp )
313      ENDIF
314      IF( kstp == nit000   )   THEN
315                 CALL iom_close( numror )     ! close input  ocean restart file
316         IF(lwm) CALL FLUSH    ( numond )     ! flush output namelist oce
317         IF(lwm) CALL FLUSH    ( numoni )     ! flush output namelist ice   
318      ENDIF
319      IF( lrst_oce         )   CALL rst_write    ( kstp )   ! write output ocean restart file
320
321      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
322      ! Trends                              (ua, va, tsa used as workspace)
323      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
324      IF( nstop == 0 ) THEN
325         IF( lk_trddyn     )   CALL trd_dwr( kstp )         ! trends: dynamics
326         IF( lk_trdtra     )   CALL trd_twr( kstp )         ! trends: active tracers
327         IF( lk_trdmld     )   CALL trd_mld( kstp )         ! trends: Mixed-layer
328         IF( lk_trdvor     )   CALL trd_vor( kstp )         ! trends: vorticity budget
329      ENDIF
330
331      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
332      ! Coupled mode
333      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
334      IF( lk_cpl           )   CALL sbc_cpl_snd( kstp )     ! coupled mode : field exchanges
335      !
336#if defined key_iomput
337      IF( kstp == nitend .OR. indic < 0 ) THEN
338                      CALL iom_context_finalize( "nemo"     ) ! needed for XIOS+AGRIF
339         IF( ln_crs ) CALL iom_context_finalize( "nemo_crs" ) !
340      ENDIF
341#endif
342      !
343      IF( nn_timing == 1 .AND.  kstp == nit000  )   CALL timing_reset
344      !
345   END SUBROUTINE stp
346
347   !!======================================================================
348END MODULE step
Note: See TracBrowser for help on using the repository browser.