New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
step.F90 in branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC – NEMO

source: branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 7795

Last change on this file since 7795 was 7256, checked in by cbricaud, 7 years ago

phaze NEMO routines in CRS branch with nemo_v3_6_STABLE branch at rev 7213 (09-09-2016) (merge -r 5519:7213 )

  • Property svn:keywords set to Id
File size: 25.0 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   use wrk_nemo
36
37   IMPLICIT NONE
38   PRIVATE
39
40   PUBLIC   stp   ! called by opa.F90
41
42   !! * Substitutions
43#  include "domzgr_substitute.h90"
44!!gm   #  include "zdfddm_substitute.h90"
45   !!----------------------------------------------------------------------
46   !! NEMO/OPA 3.7 , NEMO Consortium (2014)
47   !! $Id$
48   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
49   !!----------------------------------------------------------------------
50CONTAINS
51
52#if defined key_agrif
53   RECURSIVE SUBROUTINE stp( )
54      INTEGER             ::   kstp   ! ocean time-step index
55#else
56   SUBROUTINE stp( kstp )
57      INTEGER, INTENT(in) ::   kstp   ! ocean time-step index
58#endif
59      !!----------------------------------------------------------------------
60      !!                     ***  ROUTINE stp  ***
61      !!
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.)
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
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      !!----------------------------------------------------------------------
75      INTEGER ::   jk       ! dummy loop indice
76      INTEGER ::   indic    ! error indicator if < 0
77      INTEGER ::   kcall    ! optional integer argument (dom_vvl_sf_nxt)
78      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zrhop,zrhd
79
80      !! ---------------------------------------------------------------------
81
82#if defined key_agrif
83      kstp = nit000 + Agrif_Nb_Step()
84      IF ( lk_agrif_debug ) THEN
85         IF ( Agrif_Root() .and. lwp) Write(*,*) '---'
86         IF (lwp) Write(*,*) 'Grid Number',Agrif_Fixed(),' time step ',kstp, 'int tstep',Agrif_NbStepint()
87      ENDIF
88
89      IF ( kstp == (nit000 + 1) ) lk_agrif_fstep = .FALSE.
90
91# if defined key_iomput
92      IF( Agrif_Nbstepint() == 0 )   CALL iom_swap( cxios_context )
93# endif
94#endif
95                             indic = 0           ! reset to no error condition
96      IF( kstp == nit000 ) THEN
97         ! must be done after nemo_init for AGRIF+XIOS+OASIS
98                      CALL iom_init(      cxios_context          )  ! iom_put initialization
99         IF( lk_crs ) CALL iom_init( TRIM(cxios_context)//"_crs" )  ! initialize context for coarse grid
100      ENDIF
101
102      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
103                             CALL iom_setkt( kstp - nit000 + 1,      cxios_context          )   ! tell iom we are at time step kstp
104      IF( lk_crs     )       CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" )   ! tell iom we are at time step kstp
105
106      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
107      ! Update data, open boundaries, surface boundary condition (including sea-ice)
108      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
109      IF( lk_tide    )   CALL sbc_tide( kstp )
110      IF( lk_bdy     )  THEN
111         IF( ln_apr_dyn) CALL sbc_apr( kstp )   ! bdy_dta needs ssh_ib
112                         CALL bdy_dta ( kstp, time_offset=+1 )   ! update dynamic & tracer data at open boundaries
113      ENDIF
114                         CALL sbc    ( kstp )         ! Sea Boundary Condition (including sea-ice)
115                                                      ! clem: moved here for bdy ice purpose
116      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
117      ! Update stochastic parameters and random T/S fluctuations
118      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
119       IF( ln_sto_eos ) CALL sto_par( kstp )          ! Stochastic parameters
120       IF( ln_sto_eos ) CALL sto_pts( tsn  )          ! Random T/S fluctuations
121
122      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
123      ! Ocean physics update                (ua, va, tsa used as workspace)
124      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
125      !  THERMODYNAMICS
126                         CALL eos_rab( tsb, rab_b )       ! before local thermal/haline expension ratio at T-points
127                         CALL eos_rab( tsn, rab_n )       ! now    local thermal/haline expension ratio at T-points
128                         CALL bn2    ( tsb, rab_b, rn2b ) ! before Brunt-Vaisala frequency
129                         CALL bn2    ( tsn, rab_n, rn2  ) ! now    Brunt-Vaisala frequency
130      !
131      !  VERTICAL PHYSICS
132                         CALL zdf_bfr( kstp )         ! bottom friction (if quadratic)
133      !                                               ! Vertical eddy viscosity and diffusivity coefficients
134      IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz
135      IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz
136      IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz
137      IF( lk_zdfkpp  )   CALL zdf_kpp( kstp )            ! KPP closure scheme for Kz
138      IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value)
139         avt (:,:,:) = rn_avt0 * wmask (:,:,:)
140         avmu(:,:,:) = rn_avm0 * wumask(:,:,:)
141         avmv(:,:,:) = rn_avm0 * wvmask(:,:,:)
142      ENDIF
143      IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths
144         DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:) * tmask(:,:,jk)   ;   END DO
145      ENDIF
146      IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity
147
148      IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing
149
150      IF( lk_zdfddm .AND. .NOT. lk_zdfkpp )   &
151         &               CALL zdf_ddm( kstp )         ! double diffusive mixing
152
153                         CALL zdf_mxl( kstp )         ! mixed layer depth
154
155                                                      ! write TKE or GLS information in the restart file
156      IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' )
157      IF( lrst_oce .AND. lk_zdfgls )   CALL gls_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(  lk_asminc .AND. ln_asmiau .AND. &
209             & ln_dyninc       )  CALL dyn_asm_inc  ( kstp )   ! apply dynamics assimilation increment
210          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! subtract Neptune velocities (simplified)
211          IF( lk_bdy           )  CALL bdy_dyn3d_dmp( kstp )   ! bdy damping trends
212                                  CALL dyn_adv      ( kstp )   ! advection (vector or flux form)
213                                  CALL dyn_vor      ( kstp )   ! vorticity term including Coriolis
214                                  CALL dyn_ldf      ( kstp )   ! lateral mixing
215          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! add Neptune velocities (simplified)
216#if defined key_agrif
217          IF(.NOT. Agrif_Root())  CALL Agrif_Sponge_dyn        ! momentum sponge
218#endif
219                                  CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
220                                  CALL dyn_spg( kstp, indic )  ! surface pressure gradient
221
222                                  ua_sv(:,:,:) = ua(:,:,:)     ! Save trends (barotropic trend has been fully updated at this stage)
223                                  va_sv(:,:,:) = va(:,:,:)
224
225                                  CALL div_cur( kstp )         ! Horizontal divergence & Relative vorticity (2nd call in time-split case)
226          IF( lk_vvl     )        CALL dom_vvl_sf_nxt( kstp, kcall=2 )  ! after vertical scale factors (update depth average component)
227                                  CALL wzv           ( kstp )  ! now cross-level velocity
228      ENDIF
229
230      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
231      ! diagnostics and outputs             (ua, va, tsa used as workspace)
232      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
233      IF( lk_floats  )      CALL flo_stp( kstp )         ! drifting Floats
234      IF( lk_diahth  )      CALL dia_hth( kstp )         ! Thermocline depth (20 degres isotherm depth)
235      IF( .NOT. ln_cpl )    CALL dia_fwb( kstp )         ! Fresh water budget diagnostics
236      IF( lk_diadct  )      CALL dia_dct( kstp )         ! Transports
237      IF( lk_diaar5  )      CALL dia_ar5( kstp )         ! ar5 diag
238      IF( lk_diaharm )      CALL dia_harm( kstp )        ! Tidal harmonic analysis
239                            CALL dia_wri( kstp )         ! ocean model: outputs
240      !
241      IF( lk_crs     )      CALL crs_fld( kstp )         ! ocean model: online field coarsening & output
242
243      CALL wrk_alloc( jpi, jpj, jpk, zrhop, zrhd )
244      CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) )                 ! now in situ and potential density
245      zrhop(:,:,jpk) = 0._wp
246      CALL iom_put( 'rhop', zrhop )
247
248      CALL wrk_dealloc( jpi, jpj, jpk, zrhop, zrhd )
249
250#if defined key_top
251      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
252      ! Passive Tracer Model
253      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
254      IF( lk_crs )   THEN
255                         CALL dom_grid_crs
256
257                         CALL zdf_mxl_crs(kstp)
258                         CALL eos_crs ( tsn_crs, rhd_crs, rhop_crs, fsdept_crs(:,:,:) ) ! now in situ density for hpg computation
259
260         IF( ln_zps )    CALL zps_hde_crs( kstp, jpts, tsb_crs, gtsu_crs, gtsv_crs, rhd_crs, gru_crs, grv_crs )
261
262         IF( lk_ldfslp .AND.  .NOT. ln_traldf_grif )  &
263                         CALL ldf_slp_crs( kstp, rhd_crs, rb2_crs )
264                         CALL dom_grid_glo
265
266      ENDIF
267
268      IF( ln_crs_top )   CALL dom_grid_crs
269
270                         CALL trc_stp( kstp )         ! time-stepping
271
272      IF( ln_crs_top )   CALL dom_grid_glo
273
274#endif
275
276
277      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
278      ! Active tracers                              (ua, va used as workspace)
279      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
280                             tsa(:,:,:,:) = 0.e0            ! set tracer trends to zero
281
282      IF(  lk_asminc .AND. ln_asmiau .AND. &
283         & ln_trainc     )   CALL tra_asm_inc( kstp )       ! apply tracer assimilation increment
284                             CALL tra_sbc    ( kstp )       ! surface boundary condition
285      IF( ln_traqsr      )   CALL tra_qsr    ( kstp )       ! penetrative solar radiation qsr
286      IF( ln_trabbc      )   CALL tra_bbc    ( kstp )       ! bottom heat flux
287      IF( lk_trabbl      )   CALL tra_bbl    ( kstp )       ! advective (and/or diffusive) bottom boundary layer scheme
288      IF( ln_tradmp      )   CALL tra_dmp    ( kstp )       ! internal damping trends
289      IF( lk_bdy         )   CALL bdy_tra_dmp( kstp )       ! bdy damping trends
290                             CALL tra_adv    ( kstp )       ! horizontal & vertical advection
291      IF( lk_zdfkpp      )   CALL tra_kpp    ( kstp )       ! KPP non-local tracer fluxes
292                             CALL tra_ldf    ( kstp )       ! lateral mixing
293
294      IF( ln_diaptr      )   CALL dia_ptr                   ! Poleward adv/ldf TRansports diagnostics
295
296#if defined key_agrif
297      IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_tra          ! tracers sponge
298#endif
299                             CALL tra_zdf    ( kstp )       ! vertical mixing and after tracer fields
300
301      IF( ln_dynhpg_imp  ) THEN                             ! semi-implicit hpg (time stepping then eos)
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                             CALL eos    ( tsa, rhd, rhop, fsdept_n(:,:,:) )  ! Time-filtered in situ density for hpg computation
305            IF( ln_zps .AND. .NOT. ln_isfcav)                                &
306               &             CALL zps_hde    ( kstp, jpts, tsa, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
307               &                                           rhd, gru , grv    )  ! of t, s, rd at the last ocean level
308            IF( ln_zps .AND.       ln_isfcav)                                &
309               &             CALL zps_hde_isf( kstp, jpts, tsa, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
310               &                                           rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
311               &                                    gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
312      ELSE                                                  ! centered hpg  (eos then time stepping)
313         IF ( .NOT. lk_dynspg_ts ) THEN                     ! eos already called in time-split case
314                             CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) )  ! now in situ density for hpg computation
315         IF( ln_zps .AND. .NOT. ln_isfcav)                                   &
316               &             CALL zps_hde    ( kstp, jpts, tsn, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
317               &                                           rhd, gru , grv    )  ! of t, s, rd at the last ocean level
318         IF( ln_zps .AND.       ln_isfcav)                                   & 
319               &             CALL zps_hde_isf( kstp, jpts, tsn, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
320               &                                           rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
321               &                                    gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
322         ENDIF
323         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
324                             CALL tra_nxt( kstp )                ! tracer fields at next time step
325      ENDIF
326
327      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
328      ! Dynamics                                    (tsa used as workspace)
329      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
330      IF( lk_dynspg_ts   )  THEN
331                                                             ! revert to previously computed momentum tendencies
332                                                             ! (not using ua, va as temporary arrays during tracers' update could avoid that)
333                               ua(:,:,:) = ua_sv(:,:,:)
334                               va(:,:,:) = va_sv(:,:,:)
335                                                             ! Revert now divergence and rotational to previously computed ones
336                                                             !(needed because of the time swap in div_cur, at the beginning of each time step)
337                               hdivn(:,:,:) = hdivb(:,:,:)
338                               rotn(:,:,:)  = rotb(:,:,:) 
339
340                               CALL dyn_bfr( kstp )         ! bottom friction
341                               CALL dyn_zdf( kstp )         ! vertical diffusion
342      ELSE
343                               ua(:,:,:) = 0.e0             ! set dynamics trends to zero
344                               va(:,:,:) = 0.e0
345
346        IF(  lk_asminc .AND. ln_asmiau .AND. &
347           & ln_dyninc      )  CALL dyn_asm_inc( kstp )     ! apply dynamics assimilation increment
348        IF( ln_bkgwri )        CALL asm_bkg_wri( kstp )     ! output background fields
349        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! subtract Neptune velocities (simplified)
350        IF( lk_bdy          )  CALL bdy_dyn3d_dmp(kstp )    ! bdy damping trends
351                               CALL dyn_adv( kstp )         ! advection (vector or flux form)
352                               CALL dyn_vor( kstp )         ! vorticity term including Coriolis
353                               CALL dyn_ldf( kstp )         ! lateral mixing
354        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! add Neptune velocities (simplified)
355#if defined key_agrif
356        IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_dyn        ! momemtum sponge
357#endif
358                               CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
359                               CALL dyn_bfr( kstp )         ! bottom friction
360                               CALL dyn_zdf( kstp )         ! vertical diffusion
361                               CALL dyn_spg( kstp, indic )  ! surface pressure gradient
362      ENDIF
363                               CALL dyn_nxt( kstp )         ! lateral velocity at next time step
364
365                               CALL ssh_swp( kstp )         ! swap of sea surface height
366      IF( lk_vvl           )   CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors
367      !
368      IF( ln_diahsb        )   CALL dia_hsb( kstp )         ! - ML - global conservation diagnostics
369
370      IF( lrst_oce         )   CALL rst_write( kstp )       ! write output ocean restart file
371      IF( ln_sto_eos       )   CALL sto_rst_write( kstp )   ! write restart file for stochastic parameters
372
373#if defined key_agrif
374      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
375      ! AGRIF
376      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<     
377                               CALL Agrif_Integrate_ChildGrids( stp ) 
378
379      IF ( Agrif_NbStepint().EQ.0 ) THEN
380                               CALL Agrif_Update_Tra()      ! Update active tracers
381                               CALL Agrif_Update_Dyn()      ! Update momentum
382      ENDIF
383#endif
384      IF( lk_diaobs        )   CALL dia_obs( kstp )         ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
385
386      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
387      ! Control
388      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
389                               CALL stp_ctl( kstp, indic )
390      IF( indic < 0        )   THEN
391                               CALL ctl_stop( 'step: indic < 0' )
392                               CALL dia_wri_state( 'output.abort', kstp )
393      ENDIF
394      IF( kstp == nit000   )   THEN
395                 CALL iom_close( numror )     ! close input  ocean restart file
396         IF(lwm) CALL FLUSH    ( numond )     ! flush output namelist oce
397         IF( lwm.AND.numoni /= -1 ) CALL FLUSH    ( numoni )     ! flush output namelist ice
398      ENDIF
399
400      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
401      ! Coupled mode
402      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
403      IF( lk_oasis         )   CALL sbc_cpl_snd( kstp )     ! coupled mode : field exchanges
404      !
405#if defined key_iomput
406      IF( kstp == nitend .OR. indic < 0 ) THEN
407                      CALL iom_context_finalize(      cxios_context          ) ! needed for XIOS+AGRIF
408         IF( lk_crs ) CALL iom_context_finalize( trim(cxios_context)//"_crs" ) !
409      ENDIF
410#endif
411      !
412      IF( nn_timing == 1 .AND.  kstp == nit000  )   CALL timing_reset
413      !     
414      !
415   END SUBROUTINE stp
416
417   !!======================================================================
418END MODULE step
Note: See TracBrowser for help on using the repository browser.