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

source: branches/UKMO/dev_1d_bugfixes_tocommit/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 13191

Last change on this file since 13191 was 13191, checked in by jwhile, 4 years ago

Updates for 1d runnig

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