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

source: branches/2015/dev_r5094_UKMO_ISFCLEAN/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 5098

Last change on this file since 5098 was 5098, checked in by mathiot, 9 years ago

add wmask, wumask, wvmask and restore loop order and add flag to ignore specific isf code if no isf

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