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

source: trunk/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 4673

Last change on this file since 4673 was 4624, checked in by acc, 10 years ago

#1305. Fix slow start-up problems on some systems by introducing and using lwm logical to restrict output of merged namelists to the first (or only) processor. lwm is true only on the first processor regardless of ln_ctl. Small changes to all flavours of nemogcm.F90 are also required to write namctl and namcfg after the call to mynode which now opens output.namelist.dyn and writes nammpp.

  • Property svn:keywords set to Id
File size: 20.8 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( lk_diahth  )   CALL dia_hth( kstp )         ! Thermocline depth (20 degres isotherm depth)
205      IF( lk_diafwb  )   CALL dia_fwb( kstp )         ! Fresh water budget diagnostics
206      IF( ln_diaptr  )   CALL dia_ptr( kstp )         ! Poleward TRansports diagnostics
207      IF( lk_diadct  )   CALL dia_dct( kstp )         ! Transports
208      IF( lk_diaar5  )   CALL dia_ar5( kstp )         ! ar5 diag
209      IF( lk_diaharm )   CALL dia_harm( kstp )        ! Tidal harmonic analysis
210                         CALL dia_wri( kstp )         ! ocean model: outputs
211      !
212      IF( ln_crs     )   CALL crs_fld( kstp )         ! ocean model: online field coarsening & output
213
214
215#if defined key_top
216      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
217      ! Passive Tracer Model
218      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
219                         CALL trc_stp( kstp )         ! time-stepping
220#endif
221
222      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
223      ! Active tracers                              (ua, va used as workspace)
224      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
225                             tsa(:,:,:,:) = 0.e0            ! set tracer trends to zero
226
227      IF(  ln_asmiau .AND. &
228         & ln_trainc     )   CALL tra_asm_inc( kstp )       ! apply tracer assimilation increment
229                             CALL tra_sbc    ( kstp )       ! surface boundary condition
230      IF( ln_traqsr      )   CALL tra_qsr    ( kstp )       ! penetrative solar radiation qsr
231      IF( ln_trabbc      )   CALL tra_bbc    ( kstp )       ! bottom heat flux
232      IF( lk_trabbl      )   CALL tra_bbl    ( kstp )       ! advective (and/or diffusive) bottom boundary layer scheme
233      IF( ln_tradmp      )   CALL tra_dmp    ( kstp )       ! internal damping trends
234      IF( lk_bdy         )   CALL bdy_tra_dmp( kstp )       ! bdy damping trends
235                             CALL tra_adv    ( kstp )       ! horizontal & vertical advection
236      IF( lk_zdfkpp      )   CALL tra_kpp    ( kstp )       ! KPP non-local tracer fluxes
237                             CALL tra_ldf    ( kstp )       ! lateral mixing
238#if defined key_agrif
239      IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_tra          ! tracers sponge
240#endif
241                             CALL tra_zdf    ( kstp )       ! vertical mixing and after tracer fields
242
243      IF( ln_dynhpg_imp  ) THEN                             ! semi-implicit hpg (time stepping then eos)
244         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
245                             CALL tra_nxt( kstp )                ! tracer fields at next time step
246                             CALL eos    ( tsa, rhd, rhop, fsdept_n(:,:,:) )  ! Time-filtered in situ density for hpg computation
247         IF( ln_zps      )   CALL zps_hde( kstp, jpts, tsa, gtsu, gtsv,  &    ! zps: time filtered hor. derivative
248            &                                          rhd, gru , grv  )      ! of t, s, rd at the last ocean level
249
250      ELSE                                                  ! centered hpg  (eos then time stepping)
251         IF ( .NOT. lk_dynspg_ts ) THEN                     ! eos already called in time-split case
252                                CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) )  ! now in situ density for hpg computation
253            IF( ln_zps      )   CALL zps_hde( kstp, jpts, tsn, gtsu, gtsv,  &    ! zps: now hor. derivative
254            &                                          rhd, gru , grv  )      ! of t, s, rd at the last ocean level
255         ENDIF
256         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
257                             CALL tra_nxt( kstp )                ! tracer fields at next time step
258      ENDIF
259
260      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
261      ! Dynamics                                    (tsa used as workspace)
262      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
263      IF( lk_dynspg_ts   )  THEN
264                                                             ! revert to previously computed momentum tendencies
265                                                             ! (not using ua, va as temporary arrays during tracers' update could avoid that)
266                               ua(:,:,:) = ua_sv(:,:,:)
267                               va(:,:,:) = va_sv(:,:,:)
268                                                             ! Revert now divergence and rotational to previously computed ones
269                                                             !(needed because of the time swap in div_cur, at the beginning of each time step)
270                               hdivn(:,:,:) = hdivb(:,:,:)
271                               rotn(:,:,:)  = rotb(:,:,:) 
272
273                               CALL dyn_bfr( kstp )         ! bottom friction
274                               CALL dyn_zdf( kstp )         ! vertical diffusion
275      ELSE
276                               ua(:,:,:) = 0.e0             ! set dynamics trends to zero
277                               va(:,:,:) = 0.e0
278
279        IF(  ln_asmiau .AND. &
280           & ln_dyninc      )  CALL dyn_asm_inc( kstp )     ! apply dynamics assimilation increment
281        IF( ln_bkgwri )        CALL asm_bkg_wri( kstp )     ! output background fields
282        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! subtract Neptune velocities (simplified)
283        IF( lk_bdy          )  CALL bdy_dyn3d_dmp(kstp )    ! bdy damping trends
284                               CALL dyn_adv( kstp )         ! advection (vector or flux form)
285                               CALL dyn_vor( kstp )         ! vorticity term including Coriolis
286                               CALL dyn_ldf( kstp )         ! lateral mixing
287        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! add Neptune velocities (simplified)
288#if defined key_agrif
289        IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_dyn        ! momemtum sponge
290#endif
291                               CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
292                               CALL dyn_bfr( kstp )         ! bottom friction
293                               CALL dyn_zdf( kstp )         ! vertical diffusion
294                               CALL dyn_spg( kstp, indic )  ! surface pressure gradient
295      ENDIF
296                               CALL dyn_nxt( kstp )         ! lateral velocity at next time step
297
298                               CALL ssh_swp( kstp )         ! swap of sea surface height
299      IF( lk_vvl           )   CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors
300
301      IF( ln_diahsb        )   CALL dia_hsb( kstp )         ! - ML - global conservation diagnostics
302      IF( lk_diaobs  )         CALL dia_obs( kstp )         ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
303
304      IF( lrst_oce .AND. ln_diahsb )   CALL dia_hsb_rst( kstp, 'WRITE' )
305      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
306      ! Control and restarts
307      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
308                               CALL stp_ctl( kstp, indic )
309      IF( indic < 0        )   THEN
310                               CALL ctl_stop( 'step: indic < 0' )
311                               CALL dia_wri_state( 'output.abort', kstp )
312      ENDIF
313      IF( kstp == nit000   )   THEN
314                 CALL iom_close( numror )     ! close input  ocean restart file
315         IF(lwm) CALL FLUSH    ( numond )     ! flush output namelist oce
316         IF(lwm) CALL FLUSH    ( numoni )     ! flush output namelist ice   
317      ENDIF
318      IF( lrst_oce         )   CALL rst_write    ( kstp )   ! write output ocean restart file
319
320      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
321      ! Trends                              (ua, va, tsa used as workspace)
322      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
323      IF( nstop == 0 ) THEN
324         IF( lk_trddyn     )   CALL trd_dwr( kstp )         ! trends: dynamics
325         IF( lk_trdtra     )   CALL trd_twr( kstp )         ! trends: active tracers
326         IF( lk_trdmld     )   CALL trd_mld( kstp )         ! trends: Mixed-layer
327         IF( lk_trdvor     )   CALL trd_vor( kstp )         ! trends: vorticity budget
328      ENDIF
329
330      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
331      ! Coupled mode
332      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
333      IF( lk_cpl           )   CALL sbc_cpl_snd( kstp )     ! coupled mode : field exchanges
334      !
335#if defined key_iomput
336      IF( kstp == nitend .OR. indic < 0 ) THEN
337                      CALL iom_context_finalize( "nemo"     ) ! needed for XIOS+AGRIF
338         IF( ln_crs ) CALL iom_context_finalize( "nemo_crs" ) !
339      ENDIF
340#endif
341      !
342      IF( nn_timing == 1 .AND.  kstp == nit000  )   CALL timing_reset
343      !
344   END SUBROUTINE stp
345
346   !!======================================================================
347END MODULE step
Note: See TracBrowser for help on using the repository browser.