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.
stpmlf.F90 in NEMO/trunk/src/OCE – NEMO

source: NEMO/trunk/src/OCE/stpmlf.F90 @ 14201

Last change on this file since 14201 was 14201, checked in by smasson, 4 years ago

trunk: phase step.F90 and stpmlf.F90...

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