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/UKMO/dev_r5107_hadgem3_cplseq/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/UKMO/dev_r5107_hadgem3_cplseq/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 5479

Last change on this file since 5479 was 5479, checked in by cguiavarch, 9 years ago

Changes for sequence of coupling calls in HadGEM3.

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