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_r13333_KERNEL-08_techene_gm_HPG_SPG/src/OCE – NEMO

source: NEMO/branches/2020/dev_r13333_KERNEL-08_techene_gm_HPG_SPG/src/OCE/stpMLF.F90 @ 13362

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

trunk: Mid-year merge, merge back KERNEL-06_techene_e3

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