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/2021/dev_r14273_HPC-02_Daley_Tiling/src/OCE – NEMO

source: NEMO/branches/2021/dev_r14273_HPC-02_Daley_Tiling/src/OCE/stpmlf.F90 @ 14636

Last change on this file since 14636 was 14636, checked in by hadcv, 3 years ago

#2600: Merge in dev_r14393_HPC-03_Mele_Comm_Cleanup [14538:14609]

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