[13769] | 1 | MODULE stprk3 |
---|
[12983] | 2 | !!====================================================================== |
---|
[13769] | 3 | !! *** MODULE stprk3 *** |
---|
[12983] | 4 | !! Time-stepping : manager of the shallow water equation time stepping |
---|
[13604] | 5 | !! 3rd order Runge-Kutta scheme |
---|
[12983] | 6 | !!====================================================================== |
---|
| 7 | !! History : NEMO ! 2020-03 (A. Nasser, G. Madec) Original code from 4.0.2 |
---|
[13604] | 8 | !! - ! 2020-10 (S. Techene, G. Madec) cleanning |
---|
[12983] | 9 | !!---------------------------------------------------------------------- |
---|
| 10 | |
---|
| 11 | !!---------------------------------------------------------------------- |
---|
[13604] | 12 | !! stp_RK3 : RK3 Shallow Water Eq. time-stepping |
---|
[12983] | 13 | !!---------------------------------------------------------------------- |
---|
[13604] | 14 | USE stp_oce ! modules used in nemo_init and stp_RK3 |
---|
[12983] | 15 | ! |
---|
[13604] | 16 | USE domqco ! quasi-eulerian coordinate |
---|
| 17 | USE phycst ! physical constants |
---|
| 18 | USE usrdef_nam ! user defined namelist parameters |
---|
[12983] | 19 | |
---|
| 20 | IMPLICIT NONE |
---|
| 21 | PRIVATE |
---|
| 22 | |
---|
| 23 | PUBLIC stp_RK3 ! called by nemogcm.F90 |
---|
[13604] | 24 | |
---|
| 25 | ! !** time level indices **! |
---|
| 26 | INTEGER, PUBLIC :: Nbb, Nnn, Naa, Nrhs !: used by nemo_init |
---|
[12983] | 27 | |
---|
| 28 | !! * Substitutions |
---|
| 29 | # include "do_loop_substitute.h90" |
---|
| 30 | # include "domzgr_substitute.h90" |
---|
| 31 | !!---------------------------------------------------------------------- |
---|
| 32 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
| 33 | !! $Id: step.F90 12614 2020-03-26 14:59:52Z gm $ |
---|
| 34 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
| 35 | !!---------------------------------------------------------------------- |
---|
| 36 | CONTAINS |
---|
| 37 | |
---|
| 38 | SUBROUTINE stp_RK3( kstp ) |
---|
| 39 | !!---------------------------------------------------------------------- |
---|
| 40 | !! *** ROUTINE stp_RK3 *** |
---|
| 41 | !! |
---|
[13604] | 42 | !! ** Purpose : - RK3 Time stepping scheme for shallow water Eq. |
---|
[12983] | 43 | !! |
---|
[13604] | 44 | !! ** Method : 3rd order time stepping scheme which has 3 stages |
---|
| 45 | !! * Update calendar and forcings |
---|
| 46 | !! stage 1 : n ==> n+1/3 using variables at n |
---|
| 47 | !! - Compute the rhs of momentum |
---|
| 48 | !! - Time step ssh at Naa (n+1/3) |
---|
| 49 | !! - Time step u,v at Naa (n+1/3) |
---|
| 50 | !! - Swap time indices |
---|
| 51 | !! stage 2 : n ==> n+1/2 using variables at n and n+1/3 |
---|
| 52 | !! - Compute the rhs of momentum |
---|
| 53 | !! - Time step ssh at Naa (n+1/2) |
---|
| 54 | !! - Time step u,v at Naa (n+1/2) |
---|
| 55 | !! - Swap time indices |
---|
| 56 | !! stage 3 : n ==> n+1 using variables at n and n+1/2 |
---|
| 57 | !! - Compute the rhs of momentum |
---|
| 58 | !! - Time step ssh at Naa (n+1) |
---|
| 59 | !! - Time step u,v at Naa (n+1) |
---|
| 60 | !! - Swap time indices |
---|
| 61 | !! * Outputs and diagnostics |
---|
| 62 | !! |
---|
| 63 | !! NB: in stages 1 and 2 lateral mixing and forcing are not taken |
---|
| 64 | !! into account in the momentum RHS execpt if key_RK3all is used |
---|
[12983] | 65 | !!---------------------------------------------------------------------- |
---|
[13604] | 66 | INTEGER, INTENT(in ) :: kstp ! ocean time-step index |
---|
| 67 | ! |
---|
[12983] | 68 | INTEGER :: ji, jj, jk ! dummy loop indice |
---|
| 69 | INTEGER :: indic ! error indicator if < 0 |
---|
| 70 | REAL(wp):: z1_2rho0, z5_6, z3_4 ! local scalars |
---|
[13769] | 71 | REAL(wp):: zue3a, zue3b, zua, zrhs_u ! local scalars |
---|
| 72 | REAL(wp):: zve3a, zve3b, zva, zrhs_v ! - - |
---|
[12983] | 73 | !! --------------------------------------------------------------------- |
---|
| 74 | ! |
---|
| 75 | IF( ln_timing ) CALL timing_start('stp_RK3') |
---|
| 76 | ! |
---|
| 77 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
| 78 | ! model timestep |
---|
| 79 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 80 | ! |
---|
| 81 | IF ( kstp == nit000 ) ww(:,:,:) = 0._wp ! initialize vertical velocity one for all to zero |
---|
| 82 | |
---|
| 83 | ! |
---|
| 84 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
| 85 | ! update I/O and calendar |
---|
| 86 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 87 | indic = 0 ! reset to no error condition |
---|
| 88 | |
---|
[13604] | 89 | IF( kstp == nit000 ) THEN ! initialize IOM context |
---|
| 90 | CALL iom_init( cxios_context, ld_closedef=.FALSE. ) ! for model grid (including possible AGRIF zoom) |
---|
[12983] | 91 | CALL iom_init_closedef |
---|
| 92 | ENDIF |
---|
| 93 | IF( kstp /= nit000 ) CALL day( kstp ) ! Calendar (day was already called at nit000 in day_init) |
---|
| 94 | CALL iom_setkt( kstp - nit000 + 1, cxios_context ) ! tell IOM we are at time step kstp |
---|
| 95 | |
---|
| 96 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
[13604] | 97 | ! Update external forcing (SWE: surface boundary condition only) |
---|
[12983] | 98 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 99 | |
---|
[13604] | 100 | CALL sbc ( kstp, Nbb, Nnn ) ! Sea Boundary Condition |
---|
| 101 | |
---|
[12983] | 102 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
[13604] | 103 | ! Ocean physics update (SWE: eddy viscosity only) |
---|
[12983] | 104 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[13604] | 105 | |
---|
[12983] | 106 | IF( l_ldfdyn_time ) CALL ldf_dyn( kstp, Nbb ) ! eddy viscosity coeff. |
---|
| 107 | |
---|
| 108 | !====================================================================== |
---|
| 109 | !====================================================================== |
---|
| 110 | ! ===== RK3 ===== |
---|
| 111 | !====================================================================== |
---|
| 112 | !====================================================================== |
---|
| 113 | |
---|
| 114 | |
---|
| 115 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
[13604] | 116 | ! RK3 1st stage Ocean dynamics : u, v, ssh |
---|
[12983] | 117 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 118 | rDt = rn_Dt / 3._wp |
---|
| 119 | r1_Dt = 1._wp / rDt |
---|
[13604] | 120 | ! |
---|
| 121 | ! !== RHS of the momentum Eq. ==! |
---|
| 122 | ! |
---|
| 123 | uu(:,:,:,Nrhs) = 0._wp ! set dynamics trends to zero |
---|
| 124 | vv(:,:,:,Nrhs) = 0._wp |
---|
[12983] | 125 | |
---|
[13604] | 126 | CALL dyn_adv( kstp, Nbb, Nbb, uu, vv, Nrhs ) ! advection (VF or FF) ==> RHS |
---|
| 127 | CALL dyn_vor( kstp, Nbb, uu, vv, Nrhs ) ! vorticity ==> RHS |
---|
[12983] | 128 | #if defined key_RK3all |
---|
[13604] | 129 | CALL dyn_ldf( kstp, Nbb, Nbb, uu, vv, Nrhs ) ! lateral mixing |
---|
[12983] | 130 | #endif |
---|
[13769] | 131 | z5_6 = 5._wp/6._wp |
---|
| 132 | DO_3D( 0, 0, 0, 0, 1, jpkm1 ) |
---|
[13604] | 133 | ! ! horizontal pressure gradient |
---|
[13769] | 134 | zrhs_u = - grav * ( ssh(ji+1,jj,Nbb) - ssh(ji,jj,Nbb) ) * r1_e1u(ji,jj) |
---|
| 135 | zrhs_v = - grav * ( ssh(ji,jj+1,Nbb) - ssh(ji,jj,Nbb) ) * r1_e2v(ji,jj) |
---|
| 136 | #if defined key_RK3all |
---|
| 137 | ! ! wind stress and layer friction |
---|
| 138 | zrhs_u = zrhs_u + r1_rho0 * ( z5_6*utau_b(ji,jj) + (1._wp - z5_6)*utau(ji,jj) ) / e3u(ji,jj,jk,Nbb) & |
---|
| 139 | & - rn_rfr * uu(ji,jj,jk,Nbb) |
---|
| 140 | zrhs_v = zrhs_v + r1_rho0 * ( z5_6*vtau_b(ji,jj) + (1._wp - z5_6)*vtau(ji,jj) ) / e3v(ji,jj,jk,Nbb) & |
---|
| 141 | & - rn_rfr * vv(ji,jj,jk,Nbb) |
---|
| 142 | #endif |
---|
| 143 | ! ! ==> RHS |
---|
| 144 | uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + zrhs_u |
---|
| 145 | vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + zrhs_v |
---|
[12983] | 146 | END_3D |
---|
| 147 | ! |
---|
[13604] | 148 | ! !== Time stepping of ssh Eq. ==! (and update r3_Naa) |
---|
| 149 | ! |
---|
| 150 | CALL ssh_nxt( kstp, Nbb, Nbb, ssh, Naa ) ! after ssh |
---|
[13769] | 151 | ! ! after ssh/h_0 ratio |
---|
[13604] | 152 | CALL dom_qco_r3c( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) ) |
---|
| 153 | ! |
---|
| 154 | ! !== Time stepping of momentum Eq. ==! |
---|
| 155 | ! |
---|
| 156 | IF( ln_dynadv_vec ) THEN ! vector invariant form : applied on velocity |
---|
[13769] | 157 | DO_3D( 0, 0, 0, 0, 1,jpkm1) |
---|
[12983] | 158 | uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
| 159 | vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
| 160 | END_3D |
---|
| 161 | ELSE |
---|
[13769] | 162 | DO_3D( 0, 0, 0, 0, 1,jpkm1) ! flux form : applied on thickness weighted velocity |
---|
[13604] | 163 | zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb) |
---|
| 164 | zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb) |
---|
| 165 | zue3a = zue3b + rDt * e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
| 166 | zve3a = zve3b + rDt * e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
| 167 | ! |
---|
| 168 | uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa) |
---|
| 169 | vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa) |
---|
[12983] | 170 | END_3D |
---|
| 171 | ENDIF |
---|
[13604] | 172 | ! |
---|
[14433] | 173 | CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. ) |
---|
[14834] | 174 | IF (nn_hls==2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.) |
---|
[13769] | 175 | ! |
---|
[13604] | 176 | ! !== Swap time levels ==! |
---|
[12983] | 177 | Nrhs= Nnn |
---|
| 178 | Nnn = Naa |
---|
| 179 | Naa = Nrhs |
---|
| 180 | |
---|
| 181 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
| 182 | ! RK3 2nd stage Ocean dynamics : hdiv, ssh, e3, u, v, w |
---|
| 183 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 184 | rDt = rn_Dt / 2._wp |
---|
| 185 | r1_Dt = 1._wp / rDt |
---|
[13604] | 186 | ! |
---|
| 187 | ! !== RHS of the momentum Eq. ==! |
---|
| 188 | ! |
---|
| 189 | uu(:,:,:,Nrhs) = 0._wp ! set dynamics trends to zero |
---|
| 190 | vv(:,:,:,Nrhs) = 0._wp |
---|
| 191 | CALL dyn_adv( kstp, Nbb, Nnn, uu, vv, Nrhs ) ! advection (VF or FF) ==> RHS |
---|
| 192 | CALL dyn_vor( kstp, Nnn, uu, vv, Nrhs ) ! vorticity ==> RHS |
---|
[12983] | 193 | #if defined key_RK3all |
---|
[13604] | 194 | CALL dyn_ldf( kstp, Nbb, Nbb, uu, vv, Nrhs ) ! lateral mixing |
---|
[12983] | 195 | #endif |
---|
| 196 | ! |
---|
[13769] | 197 | z3_4 = 3._wp/4._wp |
---|
[13295] | 198 | DO_3D( 0, 0, 0, 0, 1, jpkm1 ) |
---|
[13604] | 199 | ! ! horizontal pressure gradient |
---|
[13769] | 200 | zrhs_u = - grav * ( ssh(ji+1,jj,Nnn) - ssh(ji,jj,Nnn) ) * r1_e1u(ji,jj) |
---|
| 201 | zrhs_v = - grav * ( ssh(ji,jj+1,Nnn) - ssh(ji,jj,Nnn) ) * r1_e2v(ji,jj) |
---|
[12983] | 202 | #if defined key_RK3all |
---|
[13769] | 203 | ! ! wind stress and layer friction |
---|
| 204 | zrhs_u = zrhs_u + r1_rho0 * ( z3_4*utau_b(ji,jj) + (1._wp - z3_4)*utau(ji,jj) ) / e3u(ji,jj,jk,Nnn) & |
---|
| 205 | & - rn_rfr * uu(ji,jj,jk,Nbb) |
---|
| 206 | zrhs_v = zrhs_v + r1_rho0 * ( z3_4*vtau_b(ji,jj) + (1._wp - z3_4)*vtau(ji,jj) ) / e3v(ji,jj,jk,Nnn) & |
---|
| 207 | & - rn_rfr * vv(ji,jj,jk,Nbb) |
---|
| 208 | #endif |
---|
| 209 | ! ! ==> RHS |
---|
| 210 | uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + zrhs_u |
---|
| 211 | vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + zrhs_v |
---|
[12983] | 212 | END_3D |
---|
[13604] | 213 | ! |
---|
| 214 | ! !== Time stepping of ssh Eq. ==! (and update r3_Naa) |
---|
| 215 | ! |
---|
| 216 | CALL ssh_nxt( kstp, Nbb, Nnn, ssh, Naa ) ! after ssh |
---|
[13769] | 217 | ! ! after ssh/h_0 ratio |
---|
[13604] | 218 | CALL dom_qco_r3c( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) ) |
---|
| 219 | ! |
---|
| 220 | ! !== Time stepping of momentum Eq. ==! |
---|
| 221 | ! |
---|
| 222 | IF( ln_dynadv_vec ) THEN ! vector invariant form : applied on velocity |
---|
[13769] | 223 | DO_3D( 0, 0, 0, 0, 1,jpkm1) |
---|
[12983] | 224 | uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
| 225 | vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
| 226 | END_3D |
---|
| 227 | ELSE |
---|
[13769] | 228 | DO_3D( 0, 0, 0, 0, 1,jpkm1) ! flux form : applied on thickness weighted velocity |
---|
[13604] | 229 | zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb) |
---|
| 230 | zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb) |
---|
| 231 | zue3a = zue3b + rDt * e3u(ji,jj,jk,Nnn) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
| 232 | zve3a = zve3b + rDt * e3v(ji,jj,jk,Nnn) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
| 233 | ! |
---|
| 234 | uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa) |
---|
| 235 | vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa) |
---|
[12983] | 236 | END_3D |
---|
| 237 | ENDIF |
---|
[13604] | 238 | ! |
---|
[14433] | 239 | CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. ) |
---|
[14834] | 240 | IF (nn_hls==2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.) |
---|
[13769] | 241 | ! |
---|
[13604] | 242 | ! !== Swap time levels ==! |
---|
[12983] | 243 | Nrhs= Nnn |
---|
| 244 | Nnn = Naa |
---|
| 245 | Naa = Nrhs |
---|
| 246 | |
---|
| 247 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
| 248 | ! RK3 3rd stage Ocean dynamics : hdiv, ssh, e3, u, v, w |
---|
| 249 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
| 250 | rDt = rn_Dt |
---|
| 251 | r1_Dt = 1._wp / rDt |
---|
[13604] | 252 | ! |
---|
| 253 | ! !== RHS of the momentum Eq. ==! |
---|
| 254 | ! |
---|
| 255 | uu(:,:,:,Nrhs) = 0._wp ! set dynamics trends to zero |
---|
| 256 | vv(:,:,:,Nrhs) = 0._wp |
---|
| 257 | ! |
---|
| 258 | CALL dyn_adv( kstp, Nbb, Nnn, uu, vv, Nrhs ) ! advection (VF or FF) ==> RHS |
---|
| 259 | CALL dyn_vor( kstp, Nnn, uu, vv, Nrhs ) ! vorticity ==> RHS |
---|
| 260 | CALL dyn_ldf( kstp, Nbb, Nnn, uu, vv, Nrhs ) ! lateral mixing |
---|
[12983] | 261 | |
---|
[13769] | 262 | z1_2rho0 = 0.5_wp * r1_rho0 |
---|
| 263 | DO_3D( 0, 0, 0, 0, 1,jpkm1 ) |
---|
[13604] | 264 | ! ! horizontal pressure gradient |
---|
[13769] | 265 | zrhs_u = - grav * ( ssh(ji+1,jj,Nnn) - ssh(ji,jj,Nnn) ) * r1_e1u(ji,jj) |
---|
| 266 | zrhs_v = - grav * ( ssh(ji,jj+1,Nnn) - ssh(ji,jj,Nnn) ) * r1_e2v(ji,jj) |
---|
| 267 | ! ! wind stress and layer friction |
---|
| 268 | zrhs_u = zrhs_u + z1_2rho0 * ( utau_b(ji,jj) + utau(ji,jj) ) / e3u(ji,jj,jk,Nnn) & |
---|
| 269 | & - rn_rfr * uu(ji,jj,jk,Nbb) |
---|
| 270 | zrhs_v = zrhs_v + z1_2rho0 * ( vtau_b(ji,jj) + vtau(ji,jj) ) / e3v(ji,jj,jk,Nnn) & |
---|
| 271 | & - rn_rfr * vv(ji,jj,jk,Nbb) |
---|
| 272 | ! ! ==> RHS |
---|
| 273 | uu(ji,jj,jk,Nrhs) = uu(ji,jj,jk,Nrhs) + zrhs_u |
---|
| 274 | vv(ji,jj,jk,Nrhs) = vv(ji,jj,jk,Nrhs) + zrhs_v |
---|
[12983] | 275 | END_3D |
---|
[13604] | 276 | ! |
---|
| 277 | ! !== Time stepping of ssh Eq. ==! (and update r3_Naa) |
---|
| 278 | ! |
---|
| 279 | CALL ssh_nxt( kstp, Nbb, Nnn, ssh, Naa ) ! after ssh |
---|
[13769] | 280 | ! ! after ssh/h_0 ratio |
---|
[13604] | 281 | CALL dom_qco_r3c( ssh(:,:,Naa), r3t(:,:,Naa), r3u(:,:,Naa), r3v(:,:,Naa), r3f(:,:) ) |
---|
| 282 | ! |
---|
| 283 | ! !== Time stepping of momentum Eq. ==! |
---|
| 284 | ! |
---|
| 285 | IF( ln_dynadv_vec ) THEN ! vector invariant form : applied on velocity |
---|
[13769] | 286 | DO_3D( 0, 0, 0, 0, 1,jpkm1) |
---|
| 287 | uu(ji,jj,jk,Naa) = uu(ji,jj,jk,Nbb) + rDt * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
| 288 | vv(ji,jj,jk,Naa) = vv(ji,jj,jk,Nbb) + rDt * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
[12983] | 289 | END_3D |
---|
| 290 | ! |
---|
[13604] | 291 | ELSE ! flux form : applied on thickness weighted velocity |
---|
[13769] | 292 | DO_3D( 0, 0, 0, 0, 1,jpkm1) |
---|
[12983] | 293 | zue3b = e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nbb) |
---|
| 294 | zve3b = e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nbb) |
---|
[13604] | 295 | zue3a = zue3b + rDt * e3u(ji,jj,jk,Nbb) * uu(ji,jj,jk,Nrhs) * umask(ji,jj,jk) |
---|
| 296 | zve3a = zve3b + rDt * e3v(ji,jj,jk,Nbb) * vv(ji,jj,jk,Nrhs) * vmask(ji,jj,jk) |
---|
[12983] | 297 | ! |
---|
[13604] | 298 | uu(ji,jj,jk,Naa) = zue3a / e3u(ji,jj,jk,Naa) |
---|
| 299 | vv(ji,jj,jk,Naa) = zve3a / e3v(ji,jj,jk,Naa) |
---|
[12983] | 300 | END_3D |
---|
| 301 | ENDIF |
---|
[13604] | 302 | ! |
---|
[14433] | 303 | CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. ) |
---|
[14834] | 304 | IF (nn_hls==2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.) |
---|
[13769] | 305 | ! |
---|
[13604] | 306 | ! !== Swap time levels ==! |
---|
| 307 | ! |
---|
[12983] | 308 | Nrhs = Nbb |
---|
| 309 | Nbb = Naa |
---|
| 310 | Naa = Nrhs |
---|
[13604] | 311 | |
---|
[12983] | 312 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
[14179] | 313 | ! diagnostics and outputs at Nbb (i.e. the just computed time step) |
---|
[12983] | 314 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[13604] | 315 | |
---|
[14179] | 316 | IF( ln_diacfl ) CALL dia_cfl ( kstp, Nbb ) ! Courant number diagnostics |
---|
| 317 | CALL dia_wri ( kstp, Nbb ) ! ocean model: outputs |
---|
[12983] | 318 | ! |
---|
[14179] | 319 | IF( lrst_oce ) CALL rst_write ( kstp, Nbb, Nbb ) ! write output ocean restart file |
---|
[12983] | 320 | |
---|
| 321 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
| 322 | ! Control |
---|
| 323 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
[14318] | 324 | CALL stp_ctl ( kstp , Nbb ) |
---|
[13604] | 325 | |
---|
[12983] | 326 | IF( kstp == nit000 ) THEN ! 1st time step only |
---|
| 327 | CALL iom_close( numror ) ! close input ocean restart file |
---|
| 328 | IF(lwm) CALL FLUSH ( numond ) ! flush output namelist oce |
---|
| 329 | IF(lwm .AND. numoni /= -1 ) CALL FLUSH ( numoni ) ! flush output namelist ice (if exist) |
---|
| 330 | ENDIF |
---|
| 331 | |
---|
| 332 | ! |
---|
[14239] | 333 | #if defined key_xios |
---|
[12983] | 334 | IF( kstp == nitend .OR. indic < 0 ) THEN |
---|
[13604] | 335 | CALL iom_context_finalize( cxios_context ) |
---|
[12983] | 336 | ENDIF |
---|
| 337 | #endif |
---|
| 338 | ! |
---|
| 339 | IF( ln_timing ) CALL timing_stop('stp_RK3') |
---|
| 340 | ! |
---|
| 341 | END SUBROUTINE stp_RK3 |
---|
[13604] | 342 | |
---|
[12983] | 343 | !!====================================================================== |
---|
[13769] | 344 | END MODULE stprk3 |
---|