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

source: branches/UKMO/dev_r5518_medusa_fix_restart/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 7865

Last change on this file since 7865 was 7865, checked in by marc, 7 years ago

Adding extra arguments to ZPS_HDE

File size: 24.6 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 dom_oce, ONLY: e3w_n, gdep3w_n, gdept_n
34   USE step_oce         ! time stepping definition modules
35   USE iom
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   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      !! ---------------------------------------------------------------------
79
80#if defined key_agrif
81      kstp = nit000 + Agrif_Nb_Step()
82!      IF ( Agrif_Root() .and. lwp) Write(*,*) '---'
83!      IF (lwp) Write(*,*) 'Grid Number',Agrif_Fixed(),' time step ',kstp
84      IF ( kstp == (nit000 + 1) ) lk_agrif_fstep = .FALSE.
85# if defined key_iomput
86      IF( Agrif_Nbstepint() == 0 )   CALL iom_swap( cxios_context )
87# endif
88#endif
89                             indic = 0           ! reset to no error condition
90      IF( kstp == nit000 ) THEN
91         ! must be done after nemo_init for AGRIF+XIOS+OASIS
92                      CALL iom_init(      cxios_context          )  ! iom_put initialization
93         IF( ln_crs ) CALL iom_init( TRIM(cxios_context)//"_crs" )  ! initialize context for coarse grid
94      ENDIF
95
96      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
97                             CALL iom_setkt( kstp - nit000 + 1,      cxios_context          )   ! tell iom we are at time step kstp
98      IF( ln_crs     )       CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" )   ! tell iom we are at time step kstp
99
100      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
101      ! Update data, open boundaries, surface boundary condition (including sea-ice)
102      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
103      IF( lk_tide    )   CALL sbc_tide( kstp )
104      IF( lk_bdy     )  THEN
105         IF( ln_apr_dyn) CALL sbc_apr( kstp )   ! bdy_dta needs ssh_ib
106                         CALL bdy_dta ( kstp, time_offset=+1 )   ! update dynamic & tracer data at open boundaries
107      ENDIF
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                        CALL sto_par( kstp )          ! Stochastic parameters
114
115      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
116      ! Ocean physics update                (ua, va, tsa used as workspace)
117      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
118      !  THERMODYNAMICS
119                         CALL eos_rab( tsb, rab_b )       ! before local thermal/haline expension ratio at T-points
120                         CALL eos_rab( tsn, rab_n )       ! now    local thermal/haline expension ratio at T-points
121                         CALL bn2    ( tsb, rab_b, rn2b ) ! before Brunt-Vaisala frequency
122                         CALL bn2    ( tsn, rab_n, rn2  ) ! now    Brunt-Vaisala frequency
123      !
124      !  VERTICAL PHYSICS
125                         CALL zdf_bfr( kstp )         ! bottom friction (if quadratic)
126      !                                               ! Vertical eddy viscosity and diffusivity coefficients
127      IF( lk_zdfric  )   CALL zdf_ric( kstp )            ! Richardson number dependent Kz
128      IF( lk_zdftke  )   CALL zdf_tke( kstp )            ! TKE closure scheme for Kz
129      IF( lk_zdfgls  )   CALL zdf_gls( kstp )            ! GLS closure scheme for Kz
130      IF( lk_zdfkpp  )   CALL zdf_kpp( kstp )            ! KPP closure scheme for Kz
131      IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value)
132         avt (:,:,:) = rn_avt0 * wmask (:,:,:)
133         avmu(:,:,:) = rn_avm0 * wumask(:,:,:)
134         avmv(:,:,:) = rn_avm0 * wvmask(:,:,:)
135      ENDIF
136      IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths
137         DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:) * tmask(:,:,jk)   ;   END DO
138      ENDIF
139      IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity
140
141      IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing
142
143      IF( lk_zdfddm .AND. .NOT. lk_zdfkpp )   &
144         &               CALL zdf_ddm( kstp )         ! double diffusive mixing
145
146                         CALL zdf_mxl( kstp )         ! mixed layer depth
147
148                                                      ! write TKE or GLS information in the restart file
149      IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' )
150      IF( lrst_oce .AND. lk_zdfgls )   CALL gls_rst( kstp, 'WRITE' )
151      !
152      !  LATERAL  PHYSICS
153      !
154      IF( lk_ldfslp ) THEN                            ! slope of lateral mixing
155         IF(ln_sto_eos ) CALL sto_pts( tsn )          ! Random T/S fluctuations
156                         CALL eos( tsb, rhd, gdept_0(:,:,:) )               ! before in situ density
157         IF( ln_zps .AND. .NOT. ln_isfcav)                               &
158            &            CALL zps_hde    ( kstp, jpts, tsb,              &
159            &                              e3w_n, gdept_n, gtsu, gtsv,   &  ! Partial steps: before horizontal gradient
160            &                                          rhd, gru , grv    )  ! of t, s, rd at the last ocean level
161         IF( ln_zps .AND.       ln_isfcav)                               &
162            &            CALL zps_hde_isf( kstp, jpts, tsb,              &
163            &                              e3w_n, gdept_n, gdep3w_n,     &
164            &                              gtsu, gtsv,                   &    ! Partial steps for top cell (ISF)
165            &                                          rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
166            &                                   gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the first ocean level
167         IF( ln_traldf_grif ) THEN                           ! before slope for Griffies operator
168                         CALL ldf_slp_grif( kstp )
169         ELSE
170                         CALL ldf_slp( kstp, rhd, rn2b )     ! before slope for Madec operator
171         ENDIF
172      ENDIF
173#if defined key_traldf_c2d
174      IF( lk_traldf_eiv )   CALL ldf_eiv( kstp )      ! eddy induced velocity coefficient
175#endif
176#if defined key_traldf_c3d && defined key_traldf_smag
177                          CALL ldf_tra_smag( kstp )      ! eddy induced velocity coefficient
178#  endif
179#if defined key_dynldf_c3d && defined key_dynldf_smag
180                          CALL ldf_dyn_smag( kstp )      ! eddy induced velocity coefficient
181#  endif
182
183      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
184      !  Ocean dynamics : hdiv, rot, ssh, e3, wn
185      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
186                         CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_cur)
187      IF( lk_vvl     )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors
188                         CALL wzv           ( kstp )  ! now cross-level velocity
189
190      IF( lk_dynspg_ts ) THEN 
191          ! In case the time splitting case, update almost all momentum trends here:
192          ! Note that the computation of vertical velocity above, hence "after" sea level
193          ! is necessary to compute momentum advection for the rhs of barotropic loop:
194            IF(ln_sto_eos ) CALL sto_pts( tsn )                             ! Random T/S fluctuations
195                            CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) ) ! now in situ density for hpg computation
196            IF( ln_zps .AND. .NOT. ln_isfcav)                               &
197               &            CALL zps_hde    ( kstp, jpts, tsn,              &
198               &                              e3w_n, gdept_n, gtsu, gtsv,   &    ! Partial steps: before horizontal gradient
199               &                                          rhd, gru , grv    )  ! of t, s, rd at the last ocean level
200            IF( ln_zps .AND.       ln_isfcav)                               &
201               &            CALL zps_hde_isf( kstp, jpts, tsn,              &
202               &                              e3w_n, gdept_n, gdep3w_n,     &
203               &                              gtsu, gtsv,                   &    ! Partial steps for top cell (ISF)
204               &                                          rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
205               &                                   gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
206
207                                  ua(:,:,:) = 0.e0             ! set dynamics trends to zero
208                                  va(:,:,:) = 0.e0
209          IF(  ln_asmiau .AND. &
210             & ln_dyninc       )  CALL dyn_asm_inc  ( kstp )   ! apply dynamics assimilation increment
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_wri( kstp )         ! ocean model: outputs
241      !
242      IF( ln_crs     )      CALL crs_fld( kstp )         ! ocean model: online field coarsening & output
243
244#if defined key_top
245      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
246      ! Passive Tracer Model
247      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
248                         CALL trc_stp( kstp )         ! time-stepping
249#endif
250
251
252      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
253      ! Active tracers                              (ua, va used as workspace)
254      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
255                             tsa(:,:,:,:) = 0.e0            ! set tracer trends to zero
256
257      IF(  ln_asmiau .AND. &
258         & ln_trainc     )   CALL tra_asm_inc( kstp )       ! apply tracer assimilation increment
259                             CALL tra_sbc    ( kstp )       ! surface boundary condition
260      IF( ln_traqsr      )   CALL tra_qsr    ( kstp )       ! penetrative solar radiation qsr
261      IF( ln_trabbc      )   CALL tra_bbc    ( kstp )       ! bottom heat flux
262      IF( lk_trabbl      )   CALL tra_bbl    ( kstp )       ! advective (and/or diffusive) bottom boundary layer scheme
263      IF( ln_tradmp      )   CALL tra_dmp    ( kstp )       ! internal damping trends
264      IF( lk_bdy         )   CALL bdy_tra_dmp( kstp )       ! bdy damping trends
265                             CALL tra_adv    ( kstp )       ! horizontal & vertical advection
266      IF( lk_zdfkpp      )   CALL tra_kpp    ( kstp )       ! KPP non-local tracer fluxes
267                             CALL tra_ldf    ( kstp )       ! lateral mixing
268
269      IF( ln_diaptr      )   CALL dia_ptr                   ! Poleward adv/ldf TRansports diagnostics
270
271#if defined key_agrif
272      IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_tra          ! tracers sponge
273#endif
274                             CALL tra_zdf    ( kstp )       ! vertical mixing and after tracer fields
275
276      IF( ln_dynhpg_imp  ) THEN                             ! semi-implicit hpg (time stepping then eos)
277         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
278                             CALL tra_nxt( kstp )                ! tracer fields at next time step
279            IF( ln_sto_eos ) CALL sto_pts( tsn )                 ! Random T/S fluctuations
280                             CALL eos    ( tsa, rhd, rhop, fsdept_n(:,:,:) )  ! Time-filtered in situ density for hpg computation
281            IF( ln_zps .AND. .NOT. ln_isfcav)                                &
282               &             CALL zps_hde    ( kstp, jpts, tsa,              &
283               &                               e3w_n, gdept_n, gtsu, gtsv,   &    ! Partial steps: before horizontal gradient
284               &                                           rhd, gru , grv    )  ! of t, s, rd at the last ocean level
285            IF( ln_zps .AND.       ln_isfcav)                                &
286               &             CALL zps_hde_isf( kstp, jpts, tsa,              &
287               &                               e3w_n, gdept_n, gdep3w_n,     &
288               &                               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            IF( ln_sto_eos ) CALL sto_pts( tsn )    ! Random T/S fluctuations
294                             CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) )  ! now in situ density for hpg computation
295         IF( ln_zps .AND. .NOT. ln_isfcav)                                   &
296               &             CALL zps_hde    ( kstp, jpts, tsn,              &
297               &                               e3w_n, gdept_n, gtsu, gtsv,   &    ! Partial steps: before horizontal gradient
298               &                                           rhd, gru , grv    )  ! of t, s, rd at the last ocean level
299         IF( ln_zps .AND.       ln_isfcav)                                   & 
300               &             CALL zps_hde_isf( kstp, jpts, tsn,              &
301               &                               e3w_n, gdept_n, gdep3w_n,     &
302               &                               gtsu, gtsv,                   &    ! Partial steps for top cell (ISF)
303               &                                           rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
304               &                                    gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
305         ENDIF
306         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
307                             CALL tra_nxt( kstp )                ! tracer fields at next time step
308      ENDIF
309
310      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
311      ! Dynamics                                    (tsa used as workspace)
312      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
313      IF( lk_dynspg_ts   )  THEN
314                                                             ! revert to previously computed momentum tendencies
315                                                             ! (not using ua, va as temporary arrays during tracers' update could avoid that)
316                               ua(:,:,:) = ua_sv(:,:,:)
317                               va(:,:,:) = va_sv(:,:,:)
318                                                             ! Revert now divergence and rotational to previously computed ones
319                                                             !(needed because of the time swap in div_cur, at the beginning of each time step)
320                               hdivn(:,:,:) = hdivb(:,:,:)
321                               rotn(:,:,:)  = rotb(:,:,:) 
322
323                               CALL dyn_bfr( kstp )         ! bottom friction
324                               CALL dyn_zdf( kstp )         ! vertical diffusion
325      ELSE
326                               ua(:,:,:) = 0.e0             ! set dynamics trends to zero
327                               va(:,:,:) = 0.e0
328
329        IF(  ln_asmiau .AND. &
330           & ln_dyninc      )  CALL dyn_asm_inc( kstp )     ! apply dynamics assimilation increment
331        IF( ln_bkgwri )        CALL asm_bkg_wri( kstp )     ! output background fields
332        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! subtract Neptune velocities (simplified)
333        IF( lk_bdy          )  CALL bdy_dyn3d_dmp(kstp )    ! bdy damping trends
334                               CALL dyn_adv( kstp )         ! advection (vector or flux form)
335                               CALL dyn_vor( kstp )         ! vorticity term including Coriolis
336                               CALL dyn_ldf( kstp )         ! lateral mixing
337        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! add Neptune velocities (simplified)
338#if defined key_agrif
339        IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_dyn        ! momemtum sponge
340#endif
341                               CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
342                               CALL dyn_bfr( kstp )         ! bottom friction
343                               CALL dyn_zdf( kstp )         ! vertical diffusion
344                               CALL dyn_spg( kstp, indic )  ! surface pressure gradient
345      ENDIF
346                               CALL dyn_nxt( kstp )         ! lateral velocity at next time step
347
348                               CALL ssh_swp( kstp )         ! swap of sea surface height
349      IF( lk_vvl           )   CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors
350
351      IF( ln_diahsb        )   CALL dia_hsb( kstp )         ! - ML - global conservation diagnostics
352      IF( lk_diaobs  )         CALL dia_obs( kstp )         ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
353
354      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
355      ! Control and restarts
356      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
357                               CALL stp_ctl( kstp, indic )
358      IF( indic < 0        )   THEN
359                               CALL ctl_stop( 'step: indic < 0' )
360                               CALL dia_wri_state( 'output.abort', kstp )
361      ENDIF
362      IF( kstp == nit000   )   THEN
363                 CALL iom_close( numror )     ! close input  ocean restart file
364         IF(lwm) CALL FLUSH    ( numond )     ! flush output namelist oce
365         IF( lwm.AND.numoni /= -1 ) CALL FLUSH    ( numoni )     ! flush output namelist ice
366      ENDIF
367      IF( lrst_oce         )   CALL rst_write    ( kstp )   ! write output ocean restart file
368
369      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
370      ! Coupled mode
371      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
372      IF( lk_oasis         )   CALL sbc_cpl_snd( kstp )     ! coupled mode : field exchanges
373      !
374#if defined key_iomput
375      IF( kstp == nitend .OR. indic < 0 ) THEN
376                      CALL iom_context_finalize(      cxios_context          ) ! needed for XIOS+AGRIF
377         IF( ln_crs ) CALL iom_context_finalize( trim(cxios_context)//"_crs" ) !
378      ENDIF
379#endif
380      !
381      IF( nn_timing == 1 .AND.  kstp == nit000  )   CALL timing_reset
382      !
383   END SUBROUTINE stp
384
385   !!======================================================================
386END MODULE step
Note: See TracBrowser for help on using the repository browser.