source: branches/UKMO/dev_r5518_STOPACK/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 11384

Last change on this file since 11384 was 11384, checked in by mattmartin, 2 years ago

Included Andrea Storto's STOPACK code into NEMO3.6 branch.

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