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 @ 6772

Last change on this file since 6772 was 6772, checked in by cbricaud, 8 years ago

clean in coarsening branch

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