source: NEMO/branches/2020/dev_r12377_KERNEL-06_techene_e3/src/OCE/stepLF.F90 @ 12731

Last change on this file since 12731 was 12731, checked in by techene, 11 months ago

replace h. and gde. in case key_qco is activated - quick and dirty

File size: 28.1 KB
Line 
1MODULE steplf
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   !!            3.6  !  2012-07  (J. Simeon, G. Madec. C. Ethe)  Online coarsening of outputs
27   !!            3.6  !  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   !!            3.7  !  2014-10  (G. Madec)  LDF simplication
30   !!             -   !  2014-12  (G. Madec) remove KPP scheme
31   !!             -   !  2015-11  (J. Chanut) free surface simplification (remove filtered free surface)
32   !!            4.0  !  2017-05  (G. Madec)  introduction of the vertical physics manager (zdfphy)
33   !!            4.1  !  2019-08  (A. Coward, D. Storkey) rewrite in preparation for new timestepping scheme
34   !!----------------------------------------------------------------------
35
36   !!----------------------------------------------------------------------
37   !!   stplf             : OPA system time-stepping
38   !!----------------------------------------------------------------------
39   USE step_oce         ! time stepping definition modules
40   !
41   USE iom              ! xIOs server
42   USE domqe
43   USE traatfqco         ! time filtering                   (tra_atf_qco routine)
44   USE dynatfqco         ! time filtering                   (dyn_atf_qco routine)
45   USE dynspg_ts         ! surface pressure gradient: split-explicit scheme (define un_adv)
46   USE bdydyn            ! ocean open boundary conditions (define bdy_dyn)
47
48   IMPLICIT NONE
49   PRIVATE
50
51   PUBLIC   stplf   ! called by nemogcm.F90
52
53   !!----------------------------------------------------------------------
54   !! time level indices
55   !!----------------------------------------------------------------------
56   INTEGER, PUBLIC :: Nbb, Nnn, Naa, Nrhs          !! used by nemo_init
57#  include "domzgr_substitute.h90"
58   !!----------------------------------------------------------------------
59   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
60   !! $Id: step.F90 12377 2020-02-12 14:39:06Z acc $
61   !! Software governed by the CeCILL license (see ./LICENSE)
62   !!----------------------------------------------------------------------
63CONTAINS
64
65#if defined key_agrif
66   RECURSIVE SUBROUTINE stplf( )
67      INTEGER             ::   kstp   ! ocean time-step index
68#else
69   SUBROUTINE stplf( kstp )
70      INTEGER, INTENT(in) ::   kstp   ! ocean time-step index
71#endif
72      !!----------------------------------------------------------------------
73      !!                     ***  ROUTINE stplf  ***
74      !!
75      !! ** Purpose : - Time stepping of OPA  (momentum and active tracer eqs.)
76      !!              - Time stepping of SI3 (dynamic and thermodynamic eqs.)
77      !!              - Time stepping of TRC  (passive tracer eqs.)
78      !!
79      !! ** Method  : -1- Update forcings and data
80      !!              -2- Update ocean physics
81      !!              -3- Compute the t and s trends
82      !!              -4- Update t and s
83      !!              -5- Compute the momentum trends
84      !!              -6- Update the horizontal velocity
85      !!              -7- Compute the diagnostics variables (rd,N2, hdiv,w)
86      !!              -8- Outputs and diagnostics
87      !!----------------------------------------------------------------------
88      INTEGER ::   ji, jj, jk   ! dummy loop indice
89      INTEGER ::   indic        ! error indicator if < 0
90!!gm kcall can be removed, I guess
91      INTEGER ::   kcall        ! optional integer argument (dom_vvl_sf_nxt)
92!!st patch
93      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zgdept  ! st patch
94      !! ---------------------------------------------------------------------
95#if defined key_agrif
96      kstp = nit000 + Agrif_Nb_Step()
97      Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs   ! agrif_oce module copies of time level indices
98      IF( lk_agrif_debug ) THEN
99         IF( Agrif_Root() .and. lwp)   WRITE(*,*) '---'
100         IF(lwp)   WRITE(*,*) 'Grid Number', Agrif_Fixed(),' time step ', kstp, 'int tstep', Agrif_NbStepint()
101      ENDIF
102      IF( kstp == nit000 + 1 )   lk_agrif_fstep = .FALSE.
103# if defined key_iomput
104      IF( Agrif_Nbstepint() == 0 )   CALL iom_swap( cxios_context )
105# endif
106#endif
107      !
108      IF( ln_timing )   CALL timing_start('stplf')
109      !
110      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
111      ! model timestep
112      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
113      !
114      IF( l_1st_euler ) THEN 
115         ! start or restart with Euler 1st time-step
116         rDt =  rn_Dt   
117         r1_Dt = 1._wp / rDt
118      ENDIF
119      !
120      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
121
122      ! update I/O and calendar
123      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
124                             indic = 0                ! reset to no error condition
125
126      IF( kstp == nit000 ) THEN                       ! initialize IOM context (must be done after nemo_init for AGRIF+XIOS+OASIS)
127                             CALL iom_init( cxios_context, ld_closedef=.FALSE. )   ! for model grid (including passible AGRIF zoom)
128         IF( lk_diamlr   )   CALL dia_mlr_iom_init    ! with additional setup for multiple-linear-regression analysis
129                             CALL iom_init_closedef
130         IF( ln_crs      )   CALL iom_init( TRIM(cxios_context)//"_crs" )  ! for coarse grid
131      ENDIF
132      IF( kstp /= nit000 )   CALL day( kstp )         ! Calendar (day was already called at nit000 in day_init)
133                             CALL iom_setkt( kstp - nit000 + 1,      cxios_context          )   ! tell IOM we are at time step kstp
134      IF( ln_crs         )   CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" )   ! tell IOM we are at time step kstp
135
136      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
137      ! Update external forcing (tides, open boundaries, ice shelf interaction and surface boundary condition (including sea-ice)
138      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
139      IF( ln_tide    )   CALL tide_update( kstp )                     ! update tide potential
140      IF( ln_apr_dyn )   CALL sbc_apr ( kstp )                        ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib)
141      IF( ln_bdy     )   CALL bdy_dta ( kstp, Nnn )                   ! update dynamic & tracer data at open boundaries
142      IF( ln_isf     )   CALL isf_stp ( kstp, Nnn )
143                         CALL sbc     ( kstp, Nbb, Nnn )              ! Sea Boundary Condition (including sea-ice)
144
145      !             !!st      !!!!!!!!!!!!!!!!!!!!!!!
146      !
147      ! emp = 0._wp
148      ! emp_b = 0._wp
149      ! qns = 0._wp
150      ! qsr = 0._wp
151      ! qns_b = 0._wp
152      !
153      !             !!st
154
155      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
156      ! Update stochastic parameters and random T/S fluctuations
157      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
158      IF( ln_sto_eos ) CALL sto_par( kstp )          ! Stochastic parameters
159      IF( ln_sto_eos ) CALL sto_pts( ts(:,:,:,:,Nnn)  )          ! Random T/S fluctuations
160
161      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
162      ! Ocean physics update
163      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
164      !  THERMODYNAMICS
165                         CALL eos_rab( ts(:,:,:,:,Nbb), rab_b, Nnn )       ! before local thermal/haline expension ratio at T-points
166                         CALL eos_rab( ts(:,:,:,:,Nnn), rab_n, Nnn )       ! now    local thermal/haline expension ratio at T-points
167                         CALL bn2    ( ts(:,:,:,:,Nbb), rab_b, rn2b, Nnn ) ! before Brunt-Vaisala frequency
168                         CALL bn2    ( ts(:,:,:,:,Nnn), rab_n, rn2, Nnn  ) ! now    Brunt-Vaisala frequency
169
170      !  VERTICAL PHYSICS
171                         CALL zdf_phy( kstp, Nbb, Nnn, Nrhs )   ! vertical physics update (top/bot drag, avt, avs, avm + MLD)
172
173      !  LATERAL  PHYSICS
174      !
175      IF( l_ldfslp ) THEN                             ! slope of lateral mixing
176                         CALL eos( ts(:,:,:,:,Nbb), rhd, gdept_0(:,:,:) )               ! before in situ density
177
178         IF( ln_zps .AND. .NOT. ln_isfcav)                                    &
179            &            CALL zps_hde    ( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv,  &  ! Partial steps: before horizontal gradient
180            &                                          rhd, gru , grv    )       ! of t, s, rd at the last ocean level
181
182         IF( ln_zps .AND.       ln_isfcav)                                                &
183            &            CALL zps_hde_isf( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, gtui, gtvi,  &  ! Partial steps for top cell (ISF)
184            &                                          rhd, gru , grv , grui, grvi   )       ! of t, s, rd at the first ocean level
185         IF( ln_traldf_triad ) THEN
186                         CALL ldf_slp_triad( kstp, Nbb, Nnn )             ! before slope for triad operator
187         ELSE
188                         CALL ldf_slp     ( kstp, rhd, rn2b, Nbb, Nnn )   ! before slope for standard operator
189         ENDIF
190      ENDIF
191      !                                                                        ! eddy diffusivity coeff.
192      IF( l_ldftra_time .OR. l_ldfeiv_time )   CALL ldf_tra( kstp, Nbb, Nnn )  !       and/or eiv coeff.
193      IF( l_ldfdyn_time                    )   CALL ldf_dyn( kstp, Nbb )       ! eddy viscosity coeff.
194
195      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
196      !  Ocean dynamics : hdiv, ssh, e3, u, v, w
197      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
198!!st patch
199!!st patch to be able to use substitution
200      DO jk = 1, jpk
201         zgdept(:,:,jk) = gdept(:,:,jk,Nnn)
202      END DO
203!!st end
204                            CALL ssh_nxt       ( kstp, Nbb, Nnn, ssh,  Naa )    ! after ssh (includes call to div_hor)
205      IF( .NOT.ln_linssh )  CALL dom_qe_r3c    ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa) )
206!!st      IF( .NOT.ln_linssh )  CALL dom_h_nxt     ( kstp, Nbb, Nnn,       Naa )    ! after vertical scale factors
207      !IF( .NOT.ln_linssh )  CALL dom_qe_sf_nxt ( kstp, Nbb, Nnn,      Naa )    ! after vertical scale factors
208                            CALL wzv           ( kstp, Nbb, Nnn, Naa, ww  )    ! Nnn cross-level velocity
209      IF( ln_zad_Aimp )     CALL wAimp         ( kstp,      Nnn           )  ! Adaptive-implicit vertical advection partitioning
210                            CALL eos    ( ts(:,:,:,:,Nnn), rhd, rhop, zgdept )  ! now in situ density for hpg computation
211
212
213                         uu(:,:,:,Nrhs) = 0._wp            ! set dynamics trends to zero
214                         vv(:,:,:,Nrhs) = 0._wp
215
216      IF(  lk_asminc .AND. ln_asmiau .AND. ln_dyninc )   &
217               &         CALL dyn_asm_inc   ( kstp, Nbb, Nnn, uu, vv, Nrhs )  ! apply dynamics assimilation increment
218      IF( ln_bdy     )   CALL bdy_dyn3d_dmp ( kstp, Nbb,      uu, vv, Nrhs )  ! bdy damping trends
219#if defined key_agrif
220      IF(.NOT. Agrif_Root())  &
221               &         CALL Agrif_Sponge_dyn        ! momentum sponge
222#endif
223                         CALL dyn_adv( kstp, Nbb, Nnn      , uu, vv, Nrhs )  ! advection (VF or FF)   ==> RHS
224                         CALL dyn_vor( kstp,      Nnn      , uu, vv, Nrhs )  ! vorticity              ==> RHS
225                         CALL dyn_ldf( kstp, Nbb, Nnn      , uu, vv, Nrhs )  ! lateral mixing
226      IF( ln_zdfosm  )   CALL dyn_osm( kstp,      Nnn      , uu, vv, Nrhs )  ! OSMOSIS non-local velocity fluxes ==> RHS
227                         CALL dyn_hpg( kstp,      Nnn      , uu, vv, Nrhs )  ! horizontal gradient of Hydrostatic pressure
228                         CALL dyn_spg( kstp, Nbb, Nnn, Nrhs, uu, vv, ssh, uu_b, vv_b, Naa )  ! surface pressure gradient
229
230                                                      ! With split-explicit free surface, since now transports have been updated and ssh(:,:,Nrhs) as well
231      IF( ln_dynspg_ts ) THEN                         ! vertical scale factors and vertical velocity need to be updated
232                            CALL div_hor    ( kstp, Nbb, Nnn )                ! Horizontal divergence  (2nd call in time-split case)
233         IF(.NOT.ln_linssh) CALL dom_qe_r3c ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) )
234         !IF(.NOT.ln_linssh) CALL dom_qe_sf_nxt ( kstp, Nbb, Nnn, Naa, kcall=2 )  ! after vertical scale factors (update depth average component)
235!!st         IF(.NOT.ln_linssh) CALL dom_h_nxt  ( kstp, Nbb, Nnn, Naa, kcall=2 )  ! after vertical scale factors (update depth average component)
236      ENDIF
237                            CALL dyn_zdf    ( kstp, Nbb, Nnn, Nrhs, uu, vv, Naa  )  ! vertical diffusion
238      IF( ln_dynspg_ts ) THEN                                                       ! vertical scale factors and vertical velocity need to be updated
239                            CALL wzv        ( kstp, Nbb, Nnn, Naa, ww )             ! Nnn cross-level velocity
240         IF( ln_zad_Aimp )  CALL wAimp      ( kstp,      Nnn )                      ! Adaptive-implicit vertical advection partitioning
241      ENDIF
242
243
244      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
245      ! cool skin
246      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
247      IF ( ln_diurnal )  CALL diurnal_layers( kstp )
248
249      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
250      ! diagnostics and outputs
251      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
252      IF( ln_floats  )   CALL flo_stp   ( kstp, Nbb, Nnn )      ! drifting Floats
253      IF( ln_diacfl  )   CALL dia_cfl   ( kstp,      Nnn )      ! Courant number diagnostics
254                         CALL dia_hth   ( kstp,      Nnn )      ! Thermocline depth (20 degres isotherm depth)
255      IF( ln_diadct  )   CALL dia_dct   ( kstp,      Nnn )      ! Transports
256                         CALL dia_ar5   ( kstp,      Nnn )      ! ar5 diag
257                         CALL dia_ptr   ( kstp,      Nnn )      ! Poleward adv/ldf TRansports diagnostics
258                         CALL dia_wri   ( kstp,      Nnn )      ! ocean model: outputs
259      IF( ln_crs     )   CALL crs_fld   ( kstp,      Nnn )      ! ocean model: online field coarsening & output
260      IF( lk_diadetide ) CALL dia_detide( kstp )                ! Weights computation for daily detiding of model diagnostics
261      IF( lk_diamlr  )   CALL dia_mlr                           ! Update time used in multiple-linear-regression analysis
262
263#if defined key_top
264      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
265      ! Passive Tracer Model
266      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
267                         CALL trc_stp       ( kstp, Nbb, Nnn, Nrhs, Naa )  ! time-stepping
268#endif
269
270      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
271      ! Active tracers
272      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
273                         ts(:,:,:,:,Nrhs) = 0._wp         ! set tracer trends to zero
274
275      IF(  lk_asminc .AND. ln_asmiau .AND. &
276         & ln_trainc )   CALL tra_asm_inc( kstp, Nbb, Nnn, ts, Nrhs )  ! apply tracer assimilation increment
277                         CALL tra_sbc    ( kstp,      Nnn, ts, Nrhs )  ! surface boundary condition
278      IF( ln_traqsr  )   CALL tra_qsr    ( kstp,      Nnn, ts, Nrhs )  ! penetrative solar radiation qsr
279      IF( ln_isf     )   CALL tra_isf    ( kstp,      Nnn, ts, Nrhs )  ! ice shelf heat flux
280      IF( ln_trabbc  )   CALL tra_bbc    ( kstp,      Nnn, ts, Nrhs )  ! bottom heat flux
281      IF( ln_trabbl  )   CALL tra_bbl    ( kstp, Nbb, Nnn, ts, Nrhs )  ! advective (and/or diffusive) bottom boundary layer scheme
282      IF( ln_tradmp  )   CALL tra_dmp    ( kstp, Nbb, Nnn, ts, Nrhs )  ! internal damping trends
283      IF( ln_bdy     )   CALL bdy_tra_dmp( kstp, Nbb,      ts, Nrhs )  ! bdy damping trends
284#if defined key_agrif
285      IF(.NOT. Agrif_Root())  &
286               &         CALL Agrif_Sponge_tra        ! tracers sponge
287#endif
288                         CALL tra_adv    ( kstp, Nbb, Nnn, ts, Nrhs )  ! hor. + vert. advection ==> RHS
289      IF( ln_zdfosm  )   CALL tra_osm    ( kstp,      Nnn, ts, Nrhs )  ! OSMOSIS non-local tracer fluxes ==> RHS
290      IF( lrst_oce .AND. ln_zdfosm ) &
291           &             CALL osm_rst    ( kstp,      Nnn, 'WRITE'  )  ! write OSMOSIS outputs + ww (so must do here) to restarts
292                         CALL tra_ldf    ( kstp, Nbb, Nnn, ts, Nrhs )  ! lateral mixing
293
294                         CALL tra_zdf    ( kstp, Nbb, Nnn, Nrhs, ts, Naa  )  ! vertical mixing and after tracer fields
295      IF( ln_zdfnpc  )   CALL tra_npc    ( kstp,      Nnn, Nrhs, ts, Naa  )  ! update after fields by non-penetrative convection
296
297      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
298      ! Set boundary conditions, time filter and swap time levels
299      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
300!!jc1: For agrif, it would be much better to finalize tracers/momentum here (e.g. bdy conditions) and move the swap
301!!    (and time filtering) after Agrif update. Then restart would be done after and would contain updated fields.
302!!    If so:
303!!    (i) no need to call agrif update at initialization time
304!!    (ii) no need to update "before" fields
305!!
306!!    Apart from creating new tra_swp/dyn_swp routines, this however:
307!!    (i) makes boundary conditions at initialization time computed from updated fields which is not the case between
308!!    two restarts => restartability issue. One can circumvent this, maybe, by assuming "interface separation",
309!!    e.g. a shift of the feedback interface inside child domain.
310!!    (ii) requires that all restart outputs of updated variables by agrif (e.g. passive tracers/tke/barotropic arrays) are done at the same
311!!    place.
312!!
313                         CALL zdyn_ts       ( Nnn, Naa, uu, vv )                   ! barotrope ajustment
314                         CALL finalize_sbc  ( kstp, Nbb, Naa, uu, vv, ts )                   ! boundary condifions
315                         CALL ssh_atf       ( kstp, Nbb, Nnn, Naa, ssh )                     ! time filtering of "now" sea surface height
316                         CALL dom_qe_r3c    ( ssh(:,:,Nnn), r3t_f, r3u_f, r3v_f )            ! "now" ssh/h_0 ratio from filtrered ssh
317                         CALL tra_atf_qco   ( kstp, Nbb, Nnn, Naa, ts )                      ! time filtering of "now" tracer arrays
318                         CALL dyn_atf_qco   ( kstp, Nbb, Nnn, Naa, uu, vv  )  ! time filtering of "now" velocities and scale factors
319                         r3t(:,:,Nnn) = r3t_f(:,:)
320                         r3u(:,:,Nnn) = r3u_f(:,:)
321                         r3v(:,:,Nnn) = r3v_f(:,:)
322
323      !
324      ! Swap time levels
325      Nrhs = Nbb
326      Nbb = Nnn
327      Nnn = Naa
328      Naa = Nrhs
329      !
330      !IF(.NOT.ln_linssh) CALL dom_qe_sf_update( kstp, Nbb, Nnn, Naa )  ! recompute vertical scale factors
331!!st      IF(.NOT.ln_linssh) CALL dom_h_update  ( kstp, Nbb, Nnn, Naa )  ! recompute vertical scale factors
332      !
333      IF( ln_diahsb  )   CALL dia_hsb       ( kstp, Nbb, Nnn )  ! - ML - global conservation diagnostics
334
335!!gm : This does not only concern the dynamics ==>>> add a new title
336!!gm2: why ouput restart before AGRIF update?
337!!
338!!jc: That would be better, but see comment above
339!!
340      IF( lrst_oce   )   CALL rst_write    ( kstp, Nbb, Nnn )   ! write output ocean restart file
341      IF( ln_sto_eos )   CALL sto_rst_write( kstp )   ! write restart file for stochastic parameters
342
343#if defined key_agrif
344      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
345      ! AGRIF
346      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
347                         Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs      ! agrif_oce module copies of time level indices
348                         CALL Agrif_Integrate_ChildGrids( stplf )       ! allows to finish all the Child Grids before updating
349
350                         IF( Agrif_NbStepint() == 0 ) THEN
351                            CALL Agrif_update_all( )                  ! Update all components
352                         ENDIF
353#endif
354      IF( ln_diaobs  )   CALL dia_obs      ( kstp, Nnn )      ! obs-minus-model (assimilation) diagnostics (call after dynamics update)
355
356      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
357      ! Control
358      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
359                         CALL stp_ctl      ( kstp, Nbb, Nnn, indic )
360
361      IF( kstp == nit000 ) THEN                          ! 1st time step only
362                                        CALL iom_close( numror )   ! close input  ocean restart file
363         IF(lwm)                        CALL FLUSH    ( numond )   ! flush output namelist oce
364         IF(lwm .AND. numoni /= -1 )    CALL FLUSH    ( numoni )   ! flush output namelist ice (if exist)
365      ENDIF
366
367      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>
368      ! Coupled mode
369      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
370!!gm why lk_oasis and not lk_cpl ????
371      IF( lk_oasis   )   CALL sbc_cpl_snd( kstp, Nbb, Nnn )        ! coupled mode : field exchanges
372      !
373#if defined key_iomput
374      IF( kstp == nitend .OR. indic < 0 ) THEN
375                      CALL iom_context_finalize(      cxios_context          ) ! needed for XIOS+AGRIF
376         IF(lrxios) CALL iom_context_finalize(      crxios_context          )
377         IF( ln_crs ) CALL iom_context_finalize( trim(cxios_context)//"_crs" ) !
378      ENDIF
379#endif
380      !
381      IF( l_1st_euler ) THEN         ! recover Leap-frog timestep
382         rDt = 2._wp * rn_Dt   
383         r1_Dt = 1._wp / rDt
384         l_1st_euler = .FALSE.     
385      ENDIF
386      !
387      IF( ln_timing )   CALL timing_stop('stplf')
388      !
389   END SUBROUTINE stplf
390
391
392   SUBROUTINE zdyn_ts (Kmm, Kaa, puu, pvv)
393      !!----------------------------------------------------------------------
394      !!                  ***  ROUTINE zdyn_ts  ***
395      !!
396      !! ** Purpose :   Finalize after horizontal velocity.
397      !!
398      !! ** Method  : * Ensure after velocities transport matches time splitting
399      !!             estimate (ln_dynspg_ts=T)
400      !!
401      !! ** Action :   puu(Kmm),pvv(Kmm),puu(Kaa),pvv(Kaa)   now and after horizontal velocity
402      !!----------------------------------------------------------------------
403      INTEGER                             , INTENT(in   ) :: Kmm, Kaa    ! before and after time level indices
404      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv         ! velocities
405      !
406      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   zue, zve
407      !
408      INTEGER  ::   jk   ! dummy loop indices
409      !!----------------------------------------------------------------------
410
411      IF ( ln_dynspg_ts ) THEN
412         ALLOCATE( zue(jpi,jpj)     , zve(jpi,jpj)     )
413         ! Ensure below that barotropic velocities match time splitting estimate
414         ! Compute actual transport and replace it with ts estimate at "after" time step
415         zue(:,:) = e3u(:,:,1,Kaa) * puu(:,:,1,Kaa) * umask(:,:,1)
416         zve(:,:) = e3v(:,:,1,Kaa) * pvv(:,:,1,Kaa) * vmask(:,:,1)
417         DO jk = 2, jpkm1
418            zue(:,:) = zue(:,:) + e3u(:,:,jk,Kaa) * puu(:,:,jk,Kaa) * umask(:,:,jk)
419            zve(:,:) = zve(:,:) + e3v(:,:,jk,Kaa) * pvv(:,:,jk,Kaa) * vmask(:,:,jk)
420         END DO
421         DO jk = 1, jpkm1
422            puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kaa) - zue(:,:) * r1_hu(:,:,Kaa) + uu_b(:,:,Kaa) ) * umask(:,:,jk)
423            pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kaa) - zve(:,:) * r1_hv(:,:,Kaa) + vv_b(:,:,Kaa) ) * vmask(:,:,jk)
424         END DO
425         !
426         IF( .NOT.ln_bt_fw ) THEN
427            ! Remove advective velocity from "now velocities"
428            ! prior to asselin filtering
429            ! In the forward case, this is done below after asselin filtering
430            ! so that asselin contribution is removed at the same time
431            DO jk = 1, jpkm1
432               puu(:,:,jk,Kmm) = ( puu(:,:,jk,Kmm) - un_adv(:,:)*r1_hu(:,:,Kmm) + uu_b(:,:,Kmm) )*umask(:,:,jk)
433               pvv(:,:,jk,Kmm) = ( pvv(:,:,jk,Kmm) - vn_adv(:,:)*r1_hv(:,:,Kmm) + vv_b(:,:,Kmm) )*vmask(:,:,jk)
434            END DO
435         ENDIF
436         !
437         DEALLOCATE( zue, zve )
438         !
439      ENDIF
440      !
441   END SUBROUTINE zdyn_ts
442
443
444   SUBROUTINE finalize_sbc (kt, Kbb, Kaa, puu, pvv, pts)
445      !!----------------------------------------------------------------------
446      !!                  ***  ROUTINE finalize_sbc  ***
447      !!
448      !! ** Purpose :   Apply the boundary condition on the after velocity
449      !!
450      !! ** Method  : * Apply lateral boundary conditions on after velocity
451      !!             at the local domain boundaries through lbc_lnk call,
452      !!             at the one-way open boundaries (ln_bdy=T),
453      !!             at the AGRIF zoom   boundaries (lk_agrif=T)
454      !!
455      !! ** Action :   puu(Kaa),pvv(Kaa)   after horizontal velocity and tracers
456      !!----------------------------------------------------------------------
457      INTEGER                             , INTENT(in   ) :: kt               ! ocean time-step index
458      INTEGER                             , INTENT(in   ) :: Kbb, Kaa    ! before and after time level indices
459      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv         ! velocities to be time filtered
460      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts            ! active tracers
461      !
462      ! Update after tracer and velocity on domain lateral boundaries
463      !
464#if defined key_agrif
465            CALL Agrif_tra                     ! AGRIF zoom boundaries
466            CALL Agrif_dyn( kt )               !* AGRIF zoom boundaries
467#endif
468      !                                        ! local domain boundaries  (T-point, unchanged sign)
469      CALL lbc_lnk_multi( 'finalize_sbc', puu(:,:,:,       Kaa), 'U', -1., pvv(:,:,:       ,Kaa), 'V', -1.   &
470                       &                , pts(:,:,:,jp_tem,Kaa), 'T',  1., pts(:,:,:,jp_sal,Kaa), 'T',  1. )     !* local domain boundaries
471      !
472      !                                        !* BDY open boundaries
473      IF( ln_bdy )   THEN
474                               CALL bdy_tra( kt, Kbb, pts,      Kaa )
475         IF( ln_dynspg_exp )   CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa )
476         IF( ln_dynspg_ts  )   CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa, dyn3d_only=.true. )
477      ENDIF
478      !
479   END SUBROUTINE finalize_sbc
480
481   !
482   !!======================================================================
483END MODULE steplf
Note: See TracBrowser for help on using the repository browser.