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/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE – NEMO

source: NEMO/branches/2020/dev_r13327_KERNEL-06_2_techene_e3/src/OCE/stpMLF.F90 @ 13608

Last change on this file since 13608 was 13608, checked in by techene, 4 years ago

#2385 reordering and remove unnecessary USE - sette test not passed yet

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