MODULE steplf !!====================================================================== !! *** MODULE step *** !! Time-stepping : manager of the ocean, tracer and ice time stepping !!====================================================================== !! History : OPA ! 1991-03 (G. Madec) Original code !! - ! 1991-11 (G. Madec) !! - ! 1992-06 (M. Imbard) add a first output record !! - ! 1996-04 (G. Madec) introduction of dynspg !! - ! 1996-04 (M.A. Foujols) introduction of passive tracer !! 8.0 ! 1997-06 (G. Madec) new architecture of call !! 8.2 ! 1997-06 (G. Madec, M. Imbard, G. Roullet) free surface !! - ! 1999-02 (G. Madec, N. Grima) hpg implicit !! - ! 2000-07 (J-M Molines, M. Imbard) Open Bondary Conditions !! NEMO 1.0 ! 2002-06 (G. Madec) free form, suppress macro-tasking !! - ! 2004-08 (C. Talandier) New trends organization !! - ! 2005-01 (C. Ethe) Add the KPP closure scheme !! - ! 2005-11 (G. Madec) Reorganisation of tra and dyn calls !! - ! 2006-01 (L. Debreu, C. Mazauric) Agrif implementation !! - ! 2006-07 (S. Masson) restart using iom !! 3.2 ! 2009-02 (G. Madec, R. Benshila) reintroduicing z*-coordinate !! - ! 2009-06 (S. Masson, G. Madec) TKE restart compatible with key_cpl !! 3.3 ! 2010-05 (K. Mogensen, A. Weaver, M. Martin, D. Lea) Assimilation interface !! - ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA !! 3.4 ! 2011-04 (G. Madec, C. Ethe) Merge of dtatem and dtasal !! 3.6 ! 2012-07 (J. Simeon, G. Madec. C. Ethe) Online coarsening of outputs !! 3.6 ! 2014-04 (F. Roquet, G. Madec) New equations of state !! 3.6 ! 2014-10 (E. Clementi, P. Oddo) Add Qiao vertical mixing in case of waves !! 3.7 ! 2014-10 (G. Madec) LDF simplication !! - ! 2014-12 (G. Madec) remove KPP scheme !! - ! 2015-11 (J. Chanut) free surface simplification (remove filtered free surface) !! 4.0 ! 2017-05 (G. Madec) introduction of the vertical physics manager (zdfphy) !! 4.1 ! 2019-08 (A. Coward, D. Storkey) rewrite in preparation for new timestepping scheme !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! stplf : OPA system time-stepping !!---------------------------------------------------------------------- USE step_oce ! time stepping definition modules ! USE iom ! xIOs server USE domqe USE traatfqco ! time filtering (tra_atf_qco routine) USE dynatfqco ! time filtering (dyn_atf_qco routine) USE dynspg_ts ! surface pressure gradient: split-explicit scheme (define un_adv) USE bdydyn ! ocean open boundary conditions (define bdy_dyn) IMPLICIT NONE PRIVATE PUBLIC stplf ! called by nemogcm.F90 !!---------------------------------------------------------------------- !! time level indices !!---------------------------------------------------------------------- INTEGER, PUBLIC :: Nbb, Nnn, Naa, Nrhs !! used by nemo_init !!---------------------------------------------------------------------- !! NEMO/OCE 4.0 , NEMO Consortium (2018) !! $Id: step.F90 12377 2020-02-12 14:39:06Z acc $ !! Software governed by the CeCILL license (see ./LICENSE) !!---------------------------------------------------------------------- CONTAINS #if defined key_agrif RECURSIVE SUBROUTINE stplf( ) INTEGER :: kstp ! ocean time-step index #else SUBROUTINE stplf( kstp ) INTEGER, INTENT(in) :: kstp ! ocean time-step index #endif !!---------------------------------------------------------------------- !! *** ROUTINE stplf *** !! !! ** Purpose : - Time stepping of OPA (momentum and active tracer eqs.) !! - Time stepping of SI3 (dynamic and thermodynamic eqs.) !! - Time stepping of TRC (passive tracer eqs.) !! !! ** Method : -1- Update forcings and data !! -2- Update ocean physics !! -3- Compute the t and s trends !! -4- Update t and s !! -5- Compute the momentum trends !! -6- Update the horizontal velocity !! -7- Compute the diagnostics variables (rd,N2, hdiv,w) !! -8- Outputs and diagnostics !!---------------------------------------------------------------------- INTEGER :: ji, jj, jk ! dummy loop indice INTEGER :: indic ! error indicator if < 0 !!gm kcall can be removed, I guess INTEGER :: kcall ! optional integer argument (dom_vvl_sf_nxt) !! --------------------------------------------------------------------- #if defined key_agrif kstp = nit000 + Agrif_Nb_Step() Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs ! agrif_oce module copies of time level indices IF( lk_agrif_debug ) THEN IF( Agrif_Root() .and. lwp) WRITE(*,*) '---' IF(lwp) WRITE(*,*) 'Grid Number', Agrif_Fixed(),' time step ', kstp, 'int tstep', Agrif_NbStepint() ENDIF IF( kstp == nit000 + 1 ) lk_agrif_fstep = .FALSE. # if defined key_iomput IF( Agrif_Nbstepint() == 0 ) CALL iom_swap( cxios_context ) # endif #endif ! IF( ln_timing ) CALL timing_start('stplf') ! !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! update I/O and calendar !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< indic = 0 ! reset to no error condition IF( kstp == nit000 ) THEN ! initialize IOM context (must be done after nemo_init for AGRIF+XIOS+OASIS) CALL iom_init( cxios_context, ld_closedef=.FALSE. ) ! for model grid (including passible AGRIF zoom) IF( lk_diamlr ) CALL dia_mlr_iom_init ! with additional setup for multiple-linear-regression analysis CALL iom_init_closedef IF( ln_crs ) CALL iom_init( TRIM(cxios_context)//"_crs" ) ! for coarse grid ENDIF IF( kstp /= nit000 ) CALL day( kstp ) ! Calendar (day was already called at nit000 in day_init) CALL iom_setkt( kstp - nit000 + 1, cxios_context ) ! tell IOM we are at time step kstp IF( ln_crs ) CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" ) ! tell IOM we are at time step kstp !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Update external forcing (tides, open boundaries, ice shelf interaction and surface boundary condition (including sea-ice) !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< IF( ln_tide ) CALL tide_update( kstp ) ! update tide potential IF( ln_apr_dyn ) CALL sbc_apr ( kstp ) ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib) IF( ln_bdy ) CALL bdy_dta ( kstp, Nnn ) ! update dynamic & tracer data at open boundaries IF( ln_isf ) CALL isf_stp ( kstp, Nnn ) CALL sbc ( kstp, Nbb, Nnn ) ! Sea Boundary Condition (including sea-ice) ! !!st !!!!!!!!!!!!!!!!!!!!!!! ! ! emp = 0._wp ! emp_b = 0._wp ! qns = 0._wp ! qsr = 0._wp ! qns_b = 0._wp ! ! !!st !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Update stochastic parameters and random T/S fluctuations !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> IF( ln_sto_eos ) CALL sto_par( kstp ) ! Stochastic parameters IF( ln_sto_eos ) CALL sto_pts( ts(:,:,:,:,Nnn) ) ! Random T/S fluctuations !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Ocean physics update !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< ! THERMODYNAMICS CALL eos_rab( ts(:,:,:,:,Nbb), rab_b, Nnn ) ! before local thermal/haline expension ratio at T-points CALL eos_rab( ts(:,:,:,:,Nnn), rab_n, Nnn ) ! now local thermal/haline expension ratio at T-points CALL bn2 ( ts(:,:,:,:,Nbb), rab_b, rn2b, Nnn ) ! before Brunt-Vaisala frequency CALL bn2 ( ts(:,:,:,:,Nnn), rab_n, rn2, Nnn ) ! now Brunt-Vaisala frequency ! VERTICAL PHYSICS CALL zdf_phy( kstp, Nbb, Nnn, Nrhs ) ! vertical physics update (top/bot drag, avt, avs, avm + MLD) ! LATERAL PHYSICS ! IF( l_ldfslp ) THEN ! slope of lateral mixing CALL eos( ts(:,:,:,:,Nbb), rhd, gdept_0(:,:,:) ) ! before in situ density IF( ln_zps .AND. .NOT. ln_isfcav) & & CALL zps_hde ( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, & ! Partial steps: before horizontal gradient & rhd, gru , grv ) ! of t, s, rd at the last ocean level IF( ln_zps .AND. ln_isfcav) & & CALL zps_hde_isf( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, gtui, gtvi, & ! Partial steps for top cell (ISF) & rhd, gru , grv , grui, grvi ) ! of t, s, rd at the first ocean level IF( ln_traldf_triad ) THEN CALL ldf_slp_triad( kstp, Nbb, Nnn ) ! before slope for triad operator ELSE CALL ldf_slp ( kstp, rhd, rn2b, Nbb, Nnn ) ! before slope for standard operator ENDIF ENDIF ! ! eddy diffusivity coeff. IF( l_ldftra_time .OR. l_ldfeiv_time ) CALL ldf_tra( kstp, Nbb, Nnn ) ! and/or eiv coeff. IF( l_ldfdyn_time ) CALL ldf_dyn( kstp, Nbb ) ! eddy viscosity coeff. !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Ocean dynamics : hdiv, ssh, e3, u, v, w !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< CALL ssh_nxt ( kstp, Nbb, Nnn, ssh, Naa ) ! after ssh (includes call to div_hor) IF( .NOT.ln_linssh ) CALL dom_qe_r3c ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa) ) IF( .NOT.ln_linssh ) CALL dom_h_nxt ( kstp, Nbb, Nnn, Naa ) ! after vertical scale factors !IF( .NOT.ln_linssh ) CALL dom_qe_sf_nxt ( kstp, Nbb, Nnn, Naa ) ! after vertical scale factors CALL wzv ( kstp, Nbb, Nnn, ww, Naa ) ! now cross-level velocity IF( ln_zad_Aimp ) CALL wAimp ( kstp, Nnn ) ! Adaptive-implicit vertical advection partitioning CALL eos ( ts(:,:,:,:,Nnn), rhd, rhop, gdept(:,:,:,Nnn) ) ! now in situ density for hpg computation uu(:,:,:,Nrhs) = 0._wp ! set dynamics trends to zero vv(:,:,:,Nrhs) = 0._wp IF( lk_asminc .AND. ln_asmiau .AND. ln_dyninc ) & & CALL dyn_asm_inc ( kstp, Nbb, Nnn, uu, vv, Nrhs ) ! apply dynamics assimilation increment IF( ln_bdy ) CALL bdy_dyn3d_dmp ( kstp, Nbb, uu, vv, Nrhs ) ! bdy damping trends #if defined key_agrif IF(.NOT. Agrif_Root()) & & CALL Agrif_Sponge_dyn ! momentum sponge #endif CALL dyn_adv( kstp, Nbb, Nnn , uu, vv, Nrhs ) ! advection (VF or FF) ==> RHS CALL dyn_vor( kstp, Nnn , uu, vv, Nrhs ) ! vorticity ==> RHS CALL dyn_ldf( kstp, Nbb, Nnn , uu, vv, Nrhs ) ! lateral mixing IF( ln_zdfosm ) CALL dyn_osm( kstp, Nnn , uu, vv, Nrhs ) ! OSMOSIS non-local velocity fluxes ==> RHS CALL dyn_hpg( kstp, Nnn , uu, vv, Nrhs ) ! horizontal gradient of Hydrostatic pressure CALL dyn_spg( kstp, Nbb, Nnn, Nrhs, uu, vv, ssh, uu_b, vv_b, Naa ) ! surface pressure gradient ! With split-explicit free surface, since now transports have been updated and ssh(:,:,Nrhs) as well IF( ln_dynspg_ts ) THEN ! vertical scale factors and vertical velocity need to be updated CALL div_hor ( kstp, Nbb, Nnn ) ! Horizontal divergence (2nd call in time-split case) IF(.NOT.ln_linssh) CALL dom_qe_r3c ( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) ) IF(.NOT.ln_linssh) CALL dom_qe_sf_nxt ( kstp, Nbb, Nnn, Naa, kcall=2 ) ! after vertical scale factors (update depth average component) !IF(.NOT.ln_linssh) CALL dom_h_nxt ( kstp, Nbb, Nnn, Naa, kcall=2 ) ! after vertical scale factors (update depth average component) ENDIF CALL dyn_zdf ( kstp, Nbb, Nnn, Nrhs, uu, vv, Naa ) ! vertical diffusion IF( ln_dynspg_ts ) THEN ! vertical scale factors and vertical velocity need to be updated CALL wzv ( kstp, Nbb, Nnn, ww, Naa ) ! now cross-level velocity IF( ln_zad_Aimp ) CALL wAimp ( kstp, Nnn ) ! Adaptive-implicit vertical advection partitioning ENDIF !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! cool skin !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< IF ( ln_diurnal ) CALL diurnal_layers( kstp ) !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! diagnostics and outputs !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< IF( ln_floats ) CALL flo_stp ( kstp, Nbb, Nnn ) ! drifting Floats IF( ln_diacfl ) CALL dia_cfl ( kstp, Nnn ) ! Courant number diagnostics CALL dia_hth ( kstp, Nnn ) ! Thermocline depth (20 degres isotherm depth) IF( ln_diadct ) CALL dia_dct ( kstp, Nnn ) ! Transports CALL dia_ar5 ( kstp, Nnn ) ! ar5 diag CALL dia_ptr ( kstp, Nnn ) ! Poleward adv/ldf TRansports diagnostics CALL dia_wri ( kstp, Nnn ) ! ocean model: outputs IF( ln_crs ) CALL crs_fld ( kstp, Nnn ) ! ocean model: online field coarsening & output IF( lk_diadetide ) CALL dia_detide( kstp ) ! Weights computation for daily detiding of model diagnostics IF( lk_diamlr ) CALL dia_mlr ! Update time used in multiple-linear-regression analysis #if defined key_top !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Passive Tracer Model !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< CALL trc_stp ( kstp, Nbb, Nnn, Nrhs, Naa ) ! time-stepping #endif !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Active tracers !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< ts(:,:,:,:,Nrhs) = 0._wp ! set tracer trends to zero IF( lk_asminc .AND. ln_asmiau .AND. & & ln_trainc ) CALL tra_asm_inc( kstp, Nbb, Nnn, ts, Nrhs ) ! apply tracer assimilation increment CALL tra_sbc ( kstp, Nnn, ts, Nrhs ) ! surface boundary condition IF( ln_traqsr ) CALL tra_qsr ( kstp, Nnn, ts, Nrhs ) ! penetrative solar radiation qsr IF( ln_isf ) CALL tra_isf ( kstp, Nnn, ts, Nrhs ) ! ice shelf heat flux IF( ln_trabbc ) CALL tra_bbc ( kstp, Nnn, ts, Nrhs ) ! bottom heat flux IF( ln_trabbl ) CALL tra_bbl ( kstp, Nbb, Nnn, ts, Nrhs ) ! advective (and/or diffusive) bottom boundary layer scheme IF( ln_tradmp ) CALL tra_dmp ( kstp, Nbb, Nnn, ts, Nrhs ) ! internal damping trends IF( ln_bdy ) CALL bdy_tra_dmp( kstp, Nbb, ts, Nrhs ) ! bdy damping trends #if defined key_agrif IF(.NOT. Agrif_Root()) & & CALL Agrif_Sponge_tra ! tracers sponge #endif CALL tra_adv ( kstp, Nbb, Nnn, ts, Nrhs ) ! hor. + vert. advection ==> RHS IF( ln_zdfosm ) CALL tra_osm ( kstp, Nnn, ts, Nrhs ) ! OSMOSIS non-local tracer fluxes ==> RHS IF( lrst_oce .AND. ln_zdfosm ) & & CALL osm_rst ( kstp, Nnn, 'WRITE' ) ! write OSMOSIS outputs + ww (so must do here) to restarts CALL tra_ldf ( kstp, Nbb, Nnn, ts, Nrhs ) ! lateral mixing CALL tra_zdf ( kstp, Nbb, Nnn, Nrhs, ts, Naa ) ! vertical mixing and after tracer fields IF( ln_zdfnpc ) CALL tra_npc ( kstp, Nnn, Nrhs, ts, Naa ) ! update after fields by non-penetrative convection !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Set boundary conditions, time filter and swap time levels !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< !!jc1: For agrif, it would be much better to finalize tracers/momentum here (e.g. bdy conditions) and move the swap !! (and time filtering) after Agrif update. Then restart would be done after and would contain updated fields. !! If so: !! (i) no need to call agrif update at initialization time !! (ii) no need to update "before" fields !! !! Apart from creating new tra_swp/dyn_swp routines, this however: !! (i) makes boundary conditions at initialization time computed from updated fields which is not the case between !! two restarts => restartability issue. One can circumvent this, maybe, by assuming "interface separation", !! e.g. a shift of the feedback interface inside child domain. !! (ii) requires that all restart outputs of updated variables by agrif (e.g. passive tracers/tke/barotropic arrays) are done at the same !! place. !! !!jc2: dynnxt must be the latest call. e3t(:,:,:,Nbb) are indeed updated in that routine CALL zdyn_ts ( Nnn, Naa, e3u, e3v, uu, vv ) ! barotrope ajustment CALL finalize_sbc ( kstp, Nbb, Naa, uu, vv, ts ) ! boundary condifions CALL ssh_atf ( kstp, Nbb, Nnn, Naa, ssh ) ! time filtering of "now" sea surface height CALL dom_qe_r3c ( ssh(:,:,Nnn), r3t_f, r3u_f, r3v_f ) ! "now" ssh/h_0 ratio from filtrered ssh CALL tra_atf_qco ( kstp, Nbb, Nnn, Naa, ts ) ! time filtering of "now" tracer arrays CALL dyn_atf_qco ( kstp, Nbb, Nnn, Naa, uu, vv, e3t, e3u, e3v ) ! time filtering of "now" velocities and scale factors r3t(:,:,Nnn) = r3t_f(:,:) r3u(:,:,Nnn) = r3u_f(:,:) r3v(:,:,Nnn) = r3v_f(:,:) ! ! Swap time levels Nrhs = Nbb Nbb = Nnn Nnn = Naa Naa = Nrhs ! IF(.NOT.ln_linssh) CALL dom_qe_sf_update( kstp, Nbb, Nnn, Naa ) ! recompute vertical scale factors ! IF( ln_diahsb ) CALL dia_hsb ( kstp, Nbb, Nnn ) ! - ML - global conservation diagnostics !!gm : This does not only concern the dynamics ==>>> add a new title !!gm2: why ouput restart before AGRIF update? !! !!jc: That would be better, but see comment above !! IF( lrst_oce ) CALL rst_write ( kstp, Nbb, Nnn ) ! write output ocean restart file IF( ln_sto_eos ) CALL sto_rst_write( kstp ) ! write restart file for stochastic parameters #if defined key_agrif !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! AGRIF !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs ! agrif_oce module copies of time level indices CALL Agrif_Integrate_ChildGrids( stplf ) ! allows to finish all the Child Grids before updating IF( Agrif_NbStepint() == 0 ) THEN CALL Agrif_update_all( ) ! Update all components ENDIF #endif IF( ln_diaobs ) CALL dia_obs ( kstp, Nnn ) ! obs-minus-model (assimilation) diagnostics (call after dynamics update) !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Control !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< CALL stp_ctl ( kstp, Nbb, Nnn, indic ) IF( kstp == nit000 ) THEN ! 1st time step only CALL iom_close( numror ) ! close input ocean restart file IF(lwm) CALL FLUSH ( numond ) ! flush output namelist oce IF(lwm .AND. numoni /= -1 ) CALL FLUSH ( numoni ) ! flush output namelist ice (if exist) ENDIF !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> ! Coupled mode !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< !!gm why lk_oasis and not lk_cpl ???? IF( lk_oasis ) CALL sbc_cpl_snd( kstp, Nbb, Nnn ) ! coupled mode : field exchanges ! #if defined key_iomput IF( kstp == nitend .OR. indic < 0 ) THEN CALL iom_context_finalize( cxios_context ) ! needed for XIOS+AGRIF IF(lrxios) CALL iom_context_finalize( crxios_context ) IF( ln_crs ) CALL iom_context_finalize( trim(cxios_context)//"_crs" ) ! ENDIF #endif ! IF( ln_timing ) CALL timing_stop('stplf') ! END SUBROUTINE stplf SUBROUTINE zdyn_ts (Kmm, Kaa, pe3u, pe3v, puu, pvv) !!---------------------------------------------------------------------- !! *** ROUTINE zdyn_ts *** !! !! ** Purpose : Finalize after horizontal velocity. !! !! ** Method : * Ensure after velocities transport matches time splitting !! estimate (ln_dynspg_ts=T) !! !! ** Action : puu(Kmm),pvv(Kmm),puu(Kaa),pvv(Kaa) now and after horizontal velocity !!---------------------------------------------------------------------- INTEGER , INTENT(in ) :: Kmm, Kaa ! before and after time level indices REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv ! velocities REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(in ) :: pe3u, pe3v ! scale factors ! REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zue, zve ! INTEGER :: jk ! dummy loop indices !!---------------------------------------------------------------------- IF ( ln_dynspg_ts ) THEN ALLOCATE( zue(jpi,jpj) , zve(jpi,jpj) ) ! Ensure below that barotropic velocities match time splitting estimate ! Compute actual transport and replace it with ts estimate at "after" time step zue(:,:) = pe3u(:,:,1,Kaa) * puu(:,:,1,Kaa) * umask(:,:,1) zve(:,:) = pe3v(:,:,1,Kaa) * pvv(:,:,1,Kaa) * vmask(:,:,1) DO jk = 2, jpkm1 zue(:,:) = zue(:,:) + pe3u(:,:,jk,Kaa) * puu(:,:,jk,Kaa) * umask(:,:,jk) zve(:,:) = zve(:,:) + pe3v(:,:,jk,Kaa) * pvv(:,:,jk,Kaa) * vmask(:,:,jk) END DO DO jk = 1, jpkm1 puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kaa) - zue(:,:) * r1_hu(:,:,Kaa) + uu_b(:,:,Kaa) ) * umask(:,:,jk) pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kaa) - zve(:,:) * r1_hv(:,:,Kaa) + vv_b(:,:,Kaa) ) * vmask(:,:,jk) END DO ! IF( .NOT.ln_bt_fw ) THEN ! Remove advective velocity from "now velocities" ! prior to asselin filtering ! In the forward case, this is done below after asselin filtering ! so that asselin contribution is removed at the same time DO jk = 1, jpkm1 puu(:,:,jk,Kmm) = ( puu(:,:,jk,Kmm) - un_adv(:,:)*r1_hu(:,:,Kmm) + uu_b(:,:,Kmm) )*umask(:,:,jk) pvv(:,:,jk,Kmm) = ( pvv(:,:,jk,Kmm) - vn_adv(:,:)*r1_hv(:,:,Kmm) + vv_b(:,:,Kmm) )*vmask(:,:,jk) END DO ENDIF ! DEALLOCATE( zue, zve ) ! ENDIF ! END SUBROUTINE zdyn_ts SUBROUTINE finalize_sbc (kt, Kbb, Kaa, puu, pvv, pts) !!---------------------------------------------------------------------- !! *** ROUTINE finalize_sbc *** !! !! ** Purpose : Apply the boundary condition on the after velocity !! !! ** Method : * Apply lateral boundary conditions on after velocity !! at the local domain boundaries through lbc_lnk call, !! at the one-way open boundaries (ln_bdy=T), !! at the AGRIF zoom boundaries (lk_agrif=T) !! !! ** Action : puu(Kaa),pvv(Kaa) after horizontal velocity and tracers !!---------------------------------------------------------------------- INTEGER , INTENT(in ) :: kt ! ocean time-step index INTEGER , INTENT(in ) :: Kbb, Kaa ! before and after time level indices REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv ! velocities to be time filtered REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts ! active tracers ! ! Update after tracer and velocity on domain lateral boundaries ! #if defined key_agrif CALL Agrif_tra ! AGRIF zoom boundaries CALL Agrif_dyn( kt ) !* AGRIF zoom boundaries #endif ! ! local domain boundaries (T-point, unchanged sign) CALL lbc_lnk_multi( 'finalize_sbc', puu(:,:,:, Kaa), 'U', -1., pvv(:,:,: ,Kaa), 'V', -1. & & , pts(:,:,:,jp_tem,Kaa), 'T', 1., pts(:,:,:,jp_sal,Kaa), 'T', 1. ) !* local domain boundaries ! ! !* BDY open boundaries IF( ln_bdy ) THEN CALL bdy_tra( kt, Kbb, pts, Kaa ) IF( ln_dynspg_exp ) CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa ) IF( ln_dynspg_ts ) CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa, dyn3d_only=.true. ) ENDIF ! END SUBROUTINE finalize_sbc ! !!====================================================================== END MODULE steplf