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

source: branches/UKMO/r6232_HZG_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/step.F90

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

Changes to stop reading and using the pressure forcing file when coupling to the atmospheric mean sea level pressure

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