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

source: branches/UKMO/AMM15_v3_6_STABLE_package_collate_coupling/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 10473

Last change on this file since 10473 was 10473, checked in by jcastill, 5 years ago

Merged branch UKMO/r6232_INGV1_WAVE-coupling@7620

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