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_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/stpmlf.F90 @ 14286

Last change on this file since 14286 was 14286, checked in by mcastril, 3 years ago

Reformatting and allowing to use key_qco

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