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