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

source: branches/UKMO/r6232_tracer_advection/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 9295

Last change on this file since 9295 was 9295, checked in by jcastill, 6 years ago

Remove svn keywords

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