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

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 4369

Last change on this file since 4369 was 4369, checked in by jchanut, 10 years ago

Change time step order, ticket #1232

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