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

source: branches/UKMO/r5518_INGV1_WAVE-coupling/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 7099

Last change on this file since 7099 was 7099, checked in by jcastill, 8 years ago

Changes as in r7078 of the original ING branch: branches/2015/dev_r5936_INGV1_WAVE

  • Property svn:keywords set to Id
File size: 24.3 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   !!            3.6  !  2014-10  (E. Clementi, P. Oddo) Add Qiao vertical mixing in case of waves
29   !!----------------------------------------------------------------------
30
31   !!----------------------------------------------------------------------
32   !!   stp             : OPA system time-stepping
33   !!----------------------------------------------------------------------
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 ::   ji,jj,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( ln_zdfqiao )   THEN
132        CALL zdf_qiao(kstp )                             ! Qiao vertical mixing 
133        DO jk = 1, jpkm1 
134          DO jj = 1, jpj 
135            DO ji = 1, jpi 
136              avmu(ji,jj,jk) = (avmu(ji,jj,jk) + QBvu(ji,jj,jk)) * umask(ji,jj,jk) 
137              avmv(ji,jj,jk) = (avmv(ji,jj,jk) + QBvv(ji,jj,jk)) * vmask(ji,jj,jk) 
138              avt( ji,jj,jk) = (avt( ji,jj,jk) + QBv(ji,jj,jk))  * tmask(ji,jj,jk) 
139            END DO
140          END DO
141        END DO
142      ENDIF 
143      !
144      IF( lk_zdfcst  ) THEN                              ! Constant Kz (reset avt, avm[uv] to the background value)
145         avt (:,:,:) = rn_avt0 * wmask (:,:,:)
146         avmu(:,:,:) = rn_avm0 * wumask(:,:,:)
147         avmv(:,:,:) = rn_avm0 * wvmask(:,:,:)
148      ENDIF
149      IF( ln_rnf_mouth ) THEN                         ! increase diffusivity at rivers mouths
150         DO jk = 2, nkrnf   ;   avt(:,:,jk) = avt(:,:,jk) + 2.e0 * rn_avt_rnf * rnfmsk(:,:) * tmask(:,:,jk)   ;   END DO
151      ENDIF
152      IF( ln_zdfevd  )   CALL zdf_evd( kstp )         ! enhanced vertical eddy diffusivity
153
154      IF( lk_zdftmx  )   CALL zdf_tmx( kstp )         ! tidal vertical mixing
155
156      IF( lk_zdfddm .AND. .NOT. lk_zdfkpp )   &
157         &               CALL zdf_ddm( kstp )         ! double diffusive mixing
158
159                         CALL zdf_mxl( kstp )         ! mixed layer depth
160
161                                                      ! write TKE or GLS information in the restart file
162      IF( lrst_oce .AND. lk_zdftke )   CALL tke_rst( kstp, 'WRITE' )
163      IF( lrst_oce .AND. lk_zdfgls )   CALL gls_rst( kstp, 'WRITE' )
164      !
165      !  LATERAL  PHYSICS
166      !
167      IF( lk_ldfslp ) THEN                            ! slope of lateral mixing
168         IF(ln_sto_eos ) CALL sto_pts( tsn )          ! Random T/S fluctuations
169                         CALL eos( tsb, rhd, gdept_0(:,:,:) )               ! before in situ density
170         IF( ln_zps .AND. .NOT. ln_isfcav)                               &
171            &            CALL zps_hde    ( kstp, jpts, tsb, gtsu, gtsv,  &  ! Partial steps: before horizontal gradient
172            &                                          rhd, gru , grv    )  ! of t, s, rd at the last ocean level
173         IF( ln_zps .AND.       ln_isfcav)                               &
174            &            CALL zps_hde_isf( kstp, jpts, tsb, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
175            &                                          rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
176            &                                   gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the first ocean level
177         IF( ln_traldf_grif ) THEN                           ! before slope for Griffies operator
178                         CALL ldf_slp_grif( kstp )
179         ELSE
180                         CALL ldf_slp( kstp, rhd, rn2b )     ! before slope for Madec operator
181         ENDIF
182      ENDIF
183#if defined key_traldf_c2d
184      IF( lk_traldf_eiv )   CALL ldf_eiv( kstp )      ! eddy induced velocity coefficient
185#endif
186#if defined key_traldf_c3d && defined key_traldf_smag
187                          CALL ldf_tra_smag( kstp )      ! eddy induced velocity coefficient
188#  endif
189#if defined key_dynldf_c3d && defined key_dynldf_smag
190                          CALL ldf_dyn_smag( kstp )      ! eddy induced velocity coefficient
191#  endif
192
193      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
194      !  Ocean dynamics : hdiv, rot, ssh, e3, wn
195      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
196                         CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_cur)
197      IF( lk_vvl     )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors
198                         CALL wzv           ( kstp )  ! now cross-level velocity
199
200      IF( lk_dynspg_ts ) THEN 
201          ! In case the time splitting case, update almost all momentum trends here:
202          ! Note that the computation of vertical velocity above, hence "after" sea level
203          ! is necessary to compute momentum advection for the rhs of barotropic loop:
204            IF(ln_sto_eos ) CALL sto_pts( tsn )                             ! Random T/S fluctuations
205                            CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) ) ! now in situ density for hpg computation
206            IF( ln_zps .AND. .NOT. ln_isfcav)                               &
207               &            CALL zps_hde    ( kstp, jpts, tsn, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
208               &                                          rhd, gru , grv    )  ! of t, s, rd at the last ocean level
209            IF( ln_zps .AND.       ln_isfcav)                               &
210               &            CALL zps_hde_isf( kstp, jpts, tsn, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
211               &                                          rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
212               &                                   gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
213
214                                  ua(:,:,:) = 0.e0             ! set dynamics trends to zero
215                                  va(:,:,:) = 0.e0
216          IF(  ln_asmiau .AND. &
217             & ln_dyninc       )  CALL dyn_asm_inc  ( kstp )   ! apply dynamics assimilation increment
218          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! subtract Neptune velocities (simplified)
219          IF( lk_bdy           )  CALL bdy_dyn3d_dmp( kstp )   ! bdy damping trends
220                                  CALL dyn_adv      ( kstp )   ! advection (vector or flux form)
221                                  CALL dyn_vor      ( kstp )   ! vorticity term including Coriolis
222                                  CALL dyn_ldf      ( kstp )   ! lateral mixing
223          IF( ln_stcor    )       CALL dyn_stcor    ( kstp )   ! Stokes-Coriolis forcing
224          IF( ln_neptsimp )       CALL dyn_nept_cor ( kstp )   ! add Neptune velocities (simplified)
225#if defined key_agrif
226          IF(.NOT. Agrif_Root())  CALL Agrif_Sponge_dyn        ! momentum sponge
227#endif
228                                  CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
229                                  CALL dyn_spg( kstp, indic )  ! surface pressure gradient
230
231                                  ua_sv(:,:,:) = ua(:,:,:)     ! Save trends (barotropic trend has been fully updated at this stage)
232                                  va_sv(:,:,:) = va(:,:,:)
233
234                                  CALL div_cur( kstp )         ! Horizontal divergence & Relative vorticity (2nd call in time-split case)
235          IF( lk_vvl     )        CALL dom_vvl_sf_nxt( kstp, kcall=2 )  ! after vertical scale factors (update depth average component)
236                                  CALL wzv           ( kstp )  ! now cross-level velocity
237      ENDIF
238
239      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
240      ! diagnostics and outputs             (ua, va, tsa used as workspace)
241      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
242      IF( lk_floats  )      CALL flo_stp( kstp )         ! drifting Floats
243      IF( lk_diahth  )      CALL dia_hth( kstp )         ! Thermocline depth (20 degres isotherm depth)
244      IF( .NOT. ln_cpl )    CALL dia_fwb( kstp )         ! Fresh water budget diagnostics
245      IF( lk_diadct  )      CALL dia_dct( kstp )         ! Transports
246      IF( lk_diaar5  )      CALL dia_ar5( kstp )         ! ar5 diag
247      IF( lk_diaharm )      CALL dia_harm( kstp )        ! Tidal harmonic analysis
248                            CALL dia_wri( kstp )         ! ocean model: outputs
249      !
250      IF( ln_crs     )      CALL crs_fld( kstp )         ! ocean model: online field coarsening & output
251
252#if defined key_top
253      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
254      ! Passive Tracer Model
255      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
256                         CALL trc_stp( kstp )         ! time-stepping
257#endif
258
259
260      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
261      ! Active tracers                              (ua, va used as workspace)
262      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
263                             tsa(:,:,:,:) = 0.e0            ! set tracer trends to zero
264
265      IF(  ln_asmiau .AND. &
266         & ln_trainc     )   CALL tra_asm_inc( kstp )       ! apply tracer assimilation increment
267                             CALL tra_sbc    ( kstp )       ! surface boundary condition
268      IF( ln_traqsr      )   CALL tra_qsr    ( kstp )       ! penetrative solar radiation qsr
269      IF( ln_trabbc      )   CALL tra_bbc    ( kstp )       ! bottom heat flux
270      IF( lk_trabbl      )   CALL tra_bbl    ( kstp )       ! advective (and/or diffusive) bottom boundary layer scheme
271      IF( ln_tradmp      )   CALL tra_dmp    ( kstp )       ! internal damping trends
272      IF( lk_bdy         )   CALL bdy_tra_dmp( kstp )       ! bdy damping trends
273                             CALL tra_adv    ( kstp )       ! horizontal & vertical advection
274      IF( lk_zdfkpp      )   CALL tra_kpp    ( kstp )       ! KPP non-local tracer fluxes
275                             CALL tra_ldf    ( kstp )       ! lateral mixing
276
277      IF( ln_diaptr      )   CALL dia_ptr                   ! Poleward adv/ldf TRansports diagnostics
278
279#if defined key_agrif
280      IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_tra          ! tracers sponge
281#endif
282                             CALL tra_zdf    ( kstp )       ! vertical mixing and after tracer fields
283
284      IF( ln_dynhpg_imp  ) THEN                             ! semi-implicit hpg (time stepping then eos)
285         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
286                             CALL tra_nxt( kstp )                ! tracer fields at next time step
287            IF( ln_sto_eos ) CALL sto_pts( tsn )                 ! Random T/S fluctuations
288                             CALL eos    ( tsa, rhd, rhop, fsdept_n(:,:,:) )  ! Time-filtered in situ density for hpg computation
289            IF( ln_zps .AND. .NOT. ln_isfcav)                                &
290               &             CALL zps_hde    ( kstp, jpts, tsa, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
291               &                                           rhd, gru , grv    )  ! of t, s, rd at the last ocean level
292            IF( ln_zps .AND.       ln_isfcav)                                &
293               &             CALL zps_hde_isf( kstp, jpts, tsa, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
294               &                                           rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
295               &                                    gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
296      ELSE                                                  ! centered hpg  (eos then time stepping)
297         IF ( .NOT. lk_dynspg_ts ) THEN                     ! eos already called in time-split case
298            IF( ln_sto_eos ) CALL sto_pts( tsn )    ! Random T/S fluctuations
299                             CALL eos    ( tsn, rhd, rhop, fsdept_n(:,:,:) )  ! now in situ density for hpg computation
300         IF( ln_zps .AND. .NOT. ln_isfcav)                                   &
301               &             CALL zps_hde    ( kstp, jpts, tsn, gtsu, gtsv,  &    ! Partial steps: before horizontal gradient
302               &                                           rhd, gru , grv    )  ! of t, s, rd at the last ocean level
303         IF( ln_zps .AND.       ln_isfcav)                                   & 
304               &             CALL zps_hde_isf( kstp, jpts, tsn, gtsu, gtsv,  &    ! Partial steps for top cell (ISF)
305               &                                           rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
306               &                                    gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
307         ENDIF
308         IF( ln_zdfnpc   )   CALL tra_npc( kstp )                ! update after fields by non-penetrative convection
309                             CALL tra_nxt( kstp )                ! tracer fields at next time step
310      ENDIF
311
312      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
313      ! Dynamics                                    (tsa used as workspace)
314      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
315      IF( lk_dynspg_ts   )  THEN
316                                                             ! revert to previously computed momentum tendencies
317                                                             ! (not using ua, va as temporary arrays during tracers' update could avoid that)
318                               ua(:,:,:) = ua_sv(:,:,:)
319                               va(:,:,:) = va_sv(:,:,:)
320                                                             ! Revert now divergence and rotational to previously computed ones
321                                                             !(needed because of the time swap in div_cur, at the beginning of each time step)
322                               hdivn(:,:,:) = hdivb(:,:,:)
323                               rotn(:,:,:)  = rotb(:,:,:) 
324
325                               CALL dyn_bfr( kstp )         ! bottom friction
326                               CALL dyn_zdf( kstp )         ! vertical diffusion
327      ELSE
328                               ua(:,:,:) = 0.e0             ! set dynamics trends to zero
329                               va(:,:,:) = 0.e0
330
331        IF(  ln_asmiau .AND. &
332           & ln_dyninc      )  CALL dyn_asm_inc( kstp )     ! apply dynamics assimilation increment
333        IF( ln_bkgwri )        CALL asm_bkg_wri( kstp )     ! output background fields
334        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! subtract Neptune velocities (simplified)
335        IF( lk_bdy          )  CALL bdy_dyn3d_dmp(kstp )    ! bdy damping trends
336                               CALL dyn_adv( kstp )         ! advection (vector or flux form)
337                               CALL dyn_vor( kstp )         ! vorticity term including Coriolis
338                               CALL dyn_ldf( kstp )         ! lateral mixing
339        IF( ln_neptsimp )      CALL dyn_nept_cor( kstp )    ! add Neptune velocities (simplified)
340#if defined key_agrif
341        IF(.NOT. Agrif_Root()) CALL Agrif_Sponge_dyn        ! momemtum sponge
342#endif
343                               CALL dyn_hpg( kstp )         ! horizontal gradient of Hydrostatic pressure
344                               CALL dyn_bfr( kstp )         ! bottom friction
345                               CALL dyn_zdf( kstp )         ! vertical diffusion
346                               CALL dyn_spg( kstp, indic )  ! surface pressure gradient
347      ENDIF
348                               CALL dyn_nxt( kstp )         ! lateral velocity at next time step
349
350                               CALL ssh_swp( kstp )         ! swap of sea surface height
351      IF( lk_vvl           )   CALL dom_vvl_sf_swp( kstp )  ! swap of vertical scale factors
352
353      IF( ln_diahsb        )   CALL dia_hsb( kstp )         ! - ML - global conservation diagnostics
354      IF( lk_diaobs  )         CALL dia_obs( kstp )         ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
355
356      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
357      ! Control and restarts
358      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
359                               CALL stp_ctl( kstp, indic )
360      IF( indic < 0        )   THEN
361                               CALL ctl_stop( 'step: indic < 0' )
362                               CALL dia_wri_state( 'output.abort', kstp )
363      ENDIF
364      IF( kstp == nit000   )   THEN
365                 CALL iom_close( numror )     ! close input  ocean restart file
366         IF(lwm) CALL FLUSH    ( numond )     ! flush output namelist oce
367         IF( lwm.AND.numoni /= -1 ) CALL FLUSH    ( numoni )     ! flush output namelist ice
368      ENDIF
369      IF( lrst_oce         )   CALL rst_write    ( kstp )   ! write output ocean restart file
370
371      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
372      ! Coupled mode
373      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
374      IF( lk_oasis         )   CALL sbc_cpl_snd( kstp )     ! coupled mode : field exchanges
375      !
376#if defined key_iomput
377      IF( kstp == nitend .OR. indic < 0 ) THEN
378                      CALL iom_context_finalize(      cxios_context          ) ! needed for XIOS+AGRIF
379         IF( ln_crs ) CALL iom_context_finalize( trim(cxios_context)//"_crs" ) !
380      ENDIF
381#endif
382      !
383      IF( nn_timing == 1 .AND.  kstp == nit000  )   CALL timing_reset
384      !
385   END SUBROUTINE stp
386
387   !!======================================================================
388END MODULE step
Note: See TracBrowser for help on using the repository browser.