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

source: branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 5105

Last change on this file since 5105 was 5105, checked in by cbricaud, 9 years ago

bug correction

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