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

source: branches/2014/dev_r4650_UKMO2_ice_shelves/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 4666

Last change on this file since 4666 was 4666, checked in by mathiot, 10 years ago

#1331 : add ISOMIP config files + ice shelf code

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