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

source: branches/UKMO/AMM15_v3_6_STABLE_package_collate/NEMOGCM/NEMO/OPA_SRC/step.F90 @ 11277

Last change on this file since 11277 was 11277, checked in by kingr, 5 years ago

Merged Juan's changes for running AMM15 woth wave coupling.
Corrected minor logic error to allow AMM7-uncoupled to reproduce earlier results.
Few line spacing changes to allow merging with OBS br and trunk rvn 5518.

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