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

source: branches/UKMO/dev_r5518_GO6_under_ice_relax_dr_hook/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 11738

Last change on this file since 11738 was 11738, checked in by marc, 5 years ago

The Dr Hook changes from my perl code.

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