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 @ 4328

Last change on this file since 4328 was 4328, checked in by davestorkey, 10 years ago

Remove OBC module at NEMO 3.6. See ticket #1189.

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