[11050] | 1 | MODULE dynatf |
---|
[1502] | 2 | !!========================================================================= |
---|
[11050] | 3 | !! *** MODULE dynatf *** |
---|
| 4 | !! Ocean dynamics: time filtering |
---|
[1502] | 5 | !!========================================================================= |
---|
[1438] | 6 | !! History : OPA ! 1987-02 (P. Andrich, D. L Hostis) Original code |
---|
| 7 | !! ! 1990-10 (C. Levy, G. Madec) |
---|
| 8 | !! 7.0 ! 1993-03 (M. Guyon) symetrical conditions |
---|
| 9 | !! 8.0 ! 1997-02 (G. Madec & M. Imbard) opa, release 8.0 |
---|
| 10 | !! 8.2 ! 1997-04 (A. Weaver) Euler forward step |
---|
| 11 | !! - ! 1997-06 (G. Madec) lateral boudary cond., lbc routine |
---|
| 12 | !! NEMO 1.0 ! 2002-08 (G. Madec) F90: Free form and module |
---|
| 13 | !! - ! 2002-10 (C. Talandier, A-M. Treguier) Open boundary cond. |
---|
| 14 | !! 2.0 ! 2005-11 (V. Garnier) Surface pressure gradient organization |
---|
[14072] | 15 | !! 2.3 ! 2007-07 (D. Storkey) Calls to BDY routines. |
---|
[1502] | 16 | !! 3.2 ! 2009-06 (G. Madec, R.Benshila) re-introduce the vvl option |
---|
[2528] | 17 | !! 3.3 ! 2010-09 (D. Storkey, E.O'Dea) Bug fix for BDY module |
---|
[2723] | 18 | !! 3.3 ! 2011-03 (P. Oddo) Bug fix for time-splitting+(BDY-OBC) and not VVL |
---|
[4292] | 19 | !! 3.5 ! 2013-07 (J. Chanut) Compliant with time splitting changes |
---|
[6140] | 20 | !! 3.6 ! 2014-04 (G. Madec) add the diagnostic of the time filter trends |
---|
[5930] | 21 | !! 3.7 ! 2015-11 (J. Chanut) Free surface simplification |
---|
[11475] | 22 | !! 4.1 ! 2019-08 (A. Coward, D. Storkey) Rename dynnxt.F90 -> dynatf.F90. Now just does time filtering. |
---|
[1502] | 23 | !!------------------------------------------------------------------------- |
---|
[14072] | 24 | |
---|
[11050] | 25 | !!---------------------------------------------------------------------------------------------- |
---|
| 26 | !! dyn_atf : apply Asselin time filtering to "now" velocities and vertical scale factors |
---|
| 27 | !!---------------------------------------------------------------------------------------------- |
---|
[6140] | 28 | USE oce ! ocean dynamics and tracers |
---|
| 29 | USE dom_oce ! ocean space and time domain |
---|
| 30 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
[9023] | 31 | USE sbcrnf ! river runoffs |
---|
[6140] | 32 | USE phycst ! physical constants |
---|
| 33 | USE dynadv ! dynamics: vector invariant versus flux form |
---|
| 34 | USE dynspg_ts ! surface pressure gradient: split-explicit scheme |
---|
| 35 | USE domvvl ! variable volume |
---|
[13472] | 36 | USE bdy_oce , ONLY : ln_bdy |
---|
[6140] | 37 | USE bdydta ! ocean open boundary conditions |
---|
| 38 | USE bdydyn ! ocean open boundary conditions |
---|
| 39 | USE bdyvol ! ocean open boundary condition (bdy_vol routines) |
---|
| 40 | USE trd_oce ! trends: ocean variables |
---|
| 41 | USE trddyn ! trend manager: dynamics |
---|
| 42 | USE trdken ! trend manager: kinetic energy |
---|
[12150] | 43 | USE isf_oce , ONLY: ln_isf ! ice shelf |
---|
[14072] | 44 | USE isfdynatf , ONLY: isf_dynatf ! ice shelf volume filter correction subroutine |
---|
[4990] | 45 | ! |
---|
[6140] | 46 | USE in_out_manager ! I/O manager |
---|
| 47 | USE iom ! I/O manager library |
---|
| 48 | USE lbclnk ! lateral boundary condition (or mpp link) |
---|
| 49 | USE lib_mpp ! MPP library |
---|
| 50 | USE prtctl ! Print control |
---|
| 51 | USE timing ! Timing |
---|
[13472] | 52 | USE zdfdrg , ONLY : ln_drgice_imp, rCdU_top |
---|
[2528] | 53 | #if defined key_agrif |
---|
[9570] | 54 | USE agrif_oce_interp |
---|
[11050] | 55 | #endif |
---|
[3] | 56 | |
---|
| 57 | IMPLICIT NONE |
---|
| 58 | PRIVATE |
---|
| 59 | |
---|
[11050] | 60 | PUBLIC dyn_atf ! routine called by step.F90 |
---|
[1438] | 61 | |
---|
[14143] | 62 | #if defined key_qco || defined key_linssh |
---|
[13237] | 63 | !!---------------------------------------------------------------------- |
---|
[14143] | 64 | !! 'key_qco' Quasi-Eulerian vertical coordinate |
---|
| 65 | !! OR EMPTY MODULE |
---|
| 66 | !! 'key_linssh' Fix in time vertical coordinate |
---|
[13237] | 67 | !!---------------------------------------------------------------------- |
---|
| 68 | CONTAINS |
---|
| 69 | |
---|
[14143] | 70 | SUBROUTINE dyn_atf( kt, Kbb, Kmm, Kaa, puu, pvv, pe3t, pe3u, pe3v ) |
---|
[13237] | 71 | INTEGER , INTENT(in ) :: kt ! ocean time-step index |
---|
| 72 | INTEGER , INTENT(in ) :: Kbb, Kmm, Kaa ! before and after time level indices |
---|
| 73 | REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv ! velocities to be time filtered |
---|
| 74 | REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: pe3t, pe3u, pe3v ! scale factors to be time filtered |
---|
| 75 | |
---|
| 76 | WRITE(*,*) 'dyn_atf: You should not have seen this print! error?', kt |
---|
| 77 | END SUBROUTINE dyn_atf |
---|
| 78 | |
---|
| 79 | #else |
---|
| 80 | |
---|
[12340] | 81 | !! * Substitutions |
---|
| 82 | # include "do_loop_substitute.h90" |
---|
[2715] | 83 | !!---------------------------------------------------------------------- |
---|
[9598] | 84 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
[14072] | 85 | !! $Id$ |
---|
[10068] | 86 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
[2715] | 87 | !!---------------------------------------------------------------------- |
---|
[3] | 88 | CONTAINS |
---|
| 89 | |
---|
[11050] | 90 | SUBROUTINE dyn_atf ( kt, Kbb, Kmm, Kaa, puu, pvv, pe3t, pe3u, pe3v ) |
---|
[3] | 91 | !!---------------------------------------------------------------------- |
---|
[11050] | 92 | !! *** ROUTINE dyn_atf *** |
---|
[14072] | 93 | !! |
---|
| 94 | !! ** Purpose : Finalize after horizontal velocity. Apply the boundary |
---|
[11475] | 95 | !! condition on the after velocity and apply the Asselin time |
---|
| 96 | !! filter to the now fields. |
---|
[3] | 97 | !! |
---|
[5930] | 98 | !! ** Method : * Ensure after velocities transport matches time splitting |
---|
| 99 | !! estimate (ln_dynspg_ts=T) |
---|
[3] | 100 | !! |
---|
[14072] | 101 | !! * Apply lateral boundary conditions on after velocity |
---|
[1502] | 102 | !! at the local domain boundaries through lbc_lnk call, |
---|
[7646] | 103 | !! at the one-way open boundaries (ln_bdy=T), |
---|
[4990] | 104 | !! at the AGRIF zoom boundaries (lk_agrif=T) |
---|
[3] | 105 | !! |
---|
[11475] | 106 | !! * Apply the Asselin time filter to the now fields |
---|
[1502] | 107 | !! arrays to start the next time step: |
---|
[14072] | 108 | !! (puu(Kmm),pvv(Kmm)) = (puu(Kmm),pvv(Kmm)) |
---|
[12489] | 109 | !! + rn_atfp [ (puu(Kbb),pvv(Kbb)) + (puu(Kaa),pvv(Kaa)) - 2 (puu(Kmm),pvv(Kmm)) ] |
---|
[6140] | 110 | !! Note that with flux form advection and non linear free surface, |
---|
| 111 | !! the time filter is applied on thickness weighted velocity. |
---|
[11475] | 112 | !! As a result, dyn_atf MUST be called after tra_atf. |
---|
[1502] | 113 | !! |
---|
[14072] | 114 | !! ** Action : puu(Kmm),pvv(Kmm) filtered now horizontal velocity |
---|
[3] | 115 | !!---------------------------------------------------------------------- |
---|
[11050] | 116 | INTEGER , INTENT(in ) :: kt ! ocean time-step index |
---|
| 117 | INTEGER , INTENT(in ) :: Kbb, Kmm, Kaa ! before and after time level indices |
---|
| 118 | REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv ! velocities to be time filtered |
---|
| 119 | REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: pe3t, pe3u, pe3v ! scale factors to be time filtered |
---|
[2715] | 120 | ! |
---|
[3] | 121 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
[11050] | 122 | REAL(wp) :: zue3a, zue3n, zue3b, zcoef ! local scalars |
---|
| 123 | REAL(wp) :: zve3a, zve3n, zve3b, z1_2dt ! - - |
---|
[12372] | 124 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zue, zve, zwfld |
---|
[13472] | 125 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zutau, zvtau |
---|
[14072] | 126 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: ze3t_f, ze3u_f, ze3v_f, zua, zva |
---|
[1502] | 127 | !!---------------------------------------------------------------------- |
---|
[3294] | 128 | ! |
---|
[11050] | 129 | IF( ln_timing ) CALL timing_start('dyn_atf') |
---|
[9019] | 130 | IF( ln_dynspg_ts ) ALLOCATE( zue(jpi,jpj) , zve(jpi,jpj) ) |
---|
| 131 | IF( l_trddyn ) ALLOCATE( zua(jpi,jpj,jpk) , zva(jpi,jpj,jpk) ) |
---|
[3294] | 132 | ! |
---|
[3] | 133 | IF( kt == nit000 ) THEN |
---|
| 134 | IF(lwp) WRITE(numout,*) |
---|
[11050] | 135 | IF(lwp) WRITE(numout,*) 'dyn_atf : Asselin time filtering' |
---|
[3] | 136 | IF(lwp) WRITE(numout,*) '~~~~~~~' |
---|
| 137 | ENDIF |
---|
| 138 | |
---|
[5930] | 139 | IF ( ln_dynspg_ts ) THEN |
---|
| 140 | ! Ensure below that barotropic velocities match time splitting estimate |
---|
| 141 | ! Compute actual transport and replace it with ts estimate at "after" time step |
---|
[11050] | 142 | zue(:,:) = pe3u(:,:,1,Kaa) * puu(:,:,1,Kaa) * umask(:,:,1) |
---|
| 143 | zve(:,:) = pe3v(:,:,1,Kaa) * pvv(:,:,1,Kaa) * vmask(:,:,1) |
---|
[5930] | 144 | DO jk = 2, jpkm1 |
---|
[11050] | 145 | zue(:,:) = zue(:,:) + pe3u(:,:,jk,Kaa) * puu(:,:,jk,Kaa) * umask(:,:,jk) |
---|
| 146 | zve(:,:) = zve(:,:) + pe3v(:,:,jk,Kaa) * pvv(:,:,jk,Kaa) * vmask(:,:,jk) |
---|
[1502] | 147 | END DO |
---|
| 148 | DO jk = 1, jpkm1 |
---|
[11050] | 149 | puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kaa) - zue(:,:) * r1_hu(:,:,Kaa) + uu_b(:,:,Kaa) ) * umask(:,:,jk) |
---|
| 150 | pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kaa) - zve(:,:) * r1_hv(:,:,Kaa) + vv_b(:,:,Kaa) ) * vmask(:,:,jk) |
---|
[592] | 151 | END DO |
---|
[6140] | 152 | ! |
---|
| 153 | IF( .NOT.ln_bt_fw ) THEN |
---|
[14072] | 154 | ! Remove advective velocity from "now velocities" |
---|
| 155 | ! prior to asselin filtering |
---|
| 156 | ! In the forward case, this is done below after asselin filtering |
---|
| 157 | ! so that asselin contribution is removed at the same time |
---|
[5930] | 158 | DO jk = 1, jpkm1 |
---|
[11050] | 159 | puu(:,:,jk,Kmm) = ( puu(:,:,jk,Kmm) - un_adv(:,:)*r1_hu(:,:,Kmm) + uu_b(:,:,Kmm) )*umask(:,:,jk) |
---|
| 160 | pvv(:,:,jk,Kmm) = ( pvv(:,:,jk,Kmm) - vn_adv(:,:)*r1_hv(:,:,Kmm) + vv_b(:,:,Kmm) )*vmask(:,:,jk) |
---|
[14072] | 161 | END DO |
---|
[5930] | 162 | ENDIF |
---|
[4292] | 163 | ENDIF |
---|
| 164 | |
---|
[1502] | 165 | ! Update after velocity on domain lateral boundaries |
---|
[14072] | 166 | ! -------------------------------------------------- |
---|
[5930] | 167 | # if defined key_agrif |
---|
| 168 | CALL Agrif_dyn( kt ) !* AGRIF zoom boundaries |
---|
| 169 | # endif |
---|
| 170 | ! |
---|
[13226] | 171 | CALL lbc_lnk_multi( 'dynatf', puu(:,:,:,Kaa), 'U', -1.0_wp, pvv(:,:,:,Kaa), 'V', -1.0_wp ) !* local domain boundaries |
---|
[1502] | 172 | ! |
---|
| 173 | ! !* BDY open boundaries |
---|
[11050] | 174 | IF( ln_bdy .AND. ln_dynspg_exp ) CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa ) |
---|
| 175 | IF( ln_bdy .AND. ln_dynspg_ts ) CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa, dyn3d_only=.true. ) |
---|
[3294] | 176 | |
---|
| 177 | !!$ Do we need a call to bdy_vol here?? |
---|
| 178 | ! |
---|
[4990] | 179 | IF( l_trddyn ) THEN ! prepare the atf trend computation + some diagnostics |
---|
| 180 | ! |
---|
| 181 | ! ! Kinetic energy and Conversion |
---|
[11050] | 182 | IF( ln_KE_trd ) CALL trd_dyn( puu(:,:,:,Kaa), pvv(:,:,:,Kaa), jpdyn_ken, kt, Kmm ) |
---|
[4990] | 183 | ! |
---|
| 184 | IF( ln_dyn_trd ) THEN ! 3D output: total momentum trends |
---|
[12489] | 185 | zua(:,:,:) = ( puu(:,:,:,Kaa) - puu(:,:,:,Kbb) ) * r1_Dt |
---|
| 186 | zva(:,:,:) = ( pvv(:,:,:,Kaa) - pvv(:,:,:,Kbb) ) * r1_Dt |
---|
[4990] | 187 | CALL iom_put( "utrd_tot", zua ) ! total momentum trends, except the asselin time filter |
---|
| 188 | CALL iom_put( "vtrd_tot", zva ) |
---|
| 189 | ENDIF |
---|
| 190 | ! |
---|
[11050] | 191 | zua(:,:,:) = puu(:,:,:,Kmm) ! save the now velocity before the asselin filter |
---|
| 192 | zva(:,:,:) = pvv(:,:,:,Kmm) ! (caution: there will be a shift by 1 timestep in the |
---|
[7753] | 193 | ! ! computation of the asselin filter trends) |
---|
[4990] | 194 | ENDIF |
---|
| 195 | |
---|
[1438] | 196 | ! Time filter and swap of dynamics arrays |
---|
| 197 | ! ------------------------------------------ |
---|
[14072] | 198 | |
---|
| 199 | IF( .NOT. l_1st_euler ) THEN !* Leap-Frog : Asselin time filter |
---|
[2528] | 200 | ! ! =============! |
---|
[6140] | 201 | IF( ln_linssh ) THEN ! Fixed volume ! |
---|
[2528] | 202 | ! ! =============! |
---|
[13295] | 203 | DO_3D( 1, 1, 1, 1, 1, jpkm1 ) |
---|
[12489] | 204 | puu(ji,jj,jk,Kmm) = puu(ji,jj,jk,Kmm) + rn_atfp * ( puu(ji,jj,jk,Kbb) - 2._wp * puu(ji,jj,jk,Kmm) + puu(ji,jj,jk,Kaa) ) |
---|
| 205 | pvv(ji,jj,jk,Kmm) = pvv(ji,jj,jk,Kmm) + rn_atfp * ( pvv(ji,jj,jk,Kbb) - 2._wp * pvv(ji,jj,jk,Kmm) + pvv(ji,jj,jk,Kaa) ) |
---|
[12340] | 206 | END_3D |
---|
[2528] | 207 | ! ! ================! |
---|
| 208 | ELSE ! Variable volume ! |
---|
| 209 | ! ! ================! |
---|
[11050] | 210 | ! Time-filtered scale factor at t-points |
---|
[4292] | 211 | ! ---------------------------------------------------- |
---|
[12372] | 212 | ALLOCATE( ze3t_f(jpi,jpj,jpk), zwfld(jpi,jpj) ) |
---|
[9023] | 213 | DO jk = 1, jpkm1 |
---|
[12489] | 214 | ze3t_f(:,:,jk) = pe3t(:,:,jk,Kmm) + rn_atfp * ( pe3t(:,:,jk,Kbb) - 2._wp * pe3t(:,:,jk,Kmm) + pe3t(:,:,jk,Kaa) ) |
---|
[9023] | 215 | END DO |
---|
| 216 | ! Add volume filter correction: compatibility with tracer advection scheme |
---|
[12372] | 217 | ! => time filter + conservation correction |
---|
[12489] | 218 | zcoef = rn_atfp * rn_Dt * r1_rho0 |
---|
[12372] | 219 | zwfld(:,:) = emp_b(:,:) - emp(:,:) |
---|
| 220 | IF ( ln_rnf ) zwfld(:,:) = zwfld(:,:) - ( rnf_b(:,:) - rnf(:,:) ) |
---|
[13237] | 221 | |
---|
[12372] | 222 | DO jk = 1, jpkm1 |
---|
| 223 | ze3t_f(:,:,jk) = ze3t_f(:,:,jk) - zcoef * zwfld(:,:) * tmask(:,:,jk) & |
---|
[14072] | 224 | & * pe3t(:,:,jk,Kmm) / ( ht(:,:) + 1._wp - ssmask(:,:) ) |
---|
[12372] | 225 | END DO |
---|
[2528] | 226 | ! |
---|
[12150] | 227 | ! ice shelf melting (deal separately as it can be in depth) |
---|
| 228 | ! PM: we could probably define a generic subroutine to do the in depth correction |
---|
| 229 | ! to manage rnf, isf and possibly in the futur icb, tide water glacier (...) |
---|
| 230 | ! ...(kt, coef, ktop, kbot, hz, fwf_b, fwf) |
---|
[12489] | 231 | IF ( ln_isf ) CALL isf_dynatf( kt, Kmm, ze3t_f, rn_atfp * rn_Dt ) |
---|
[12150] | 232 | ! |
---|
[11050] | 233 | pe3t(:,:,1:jpkm1,Kmm) = ze3t_f(:,:,1:jpkm1) ! filtered scale factor at T-points |
---|
| 234 | ! |
---|
[6140] | 235 | IF( ln_dynadv_vec ) THEN ! Asselin filter applied on velocity |
---|
| 236 | ! Before filtered scale factor at (u/v)-points |
---|
[11050] | 237 | CALL dom_vvl_interpol( pe3t(:,:,:,Kmm), pe3u(:,:,:,Kmm), 'U' ) |
---|
| 238 | CALL dom_vvl_interpol( pe3t(:,:,:,Kmm), pe3v(:,:,:,Kmm), 'V' ) |
---|
[13295] | 239 | DO_3D( 1, 1, 1, 1, 1, jpkm1 ) |
---|
[12489] | 240 | puu(ji,jj,jk,Kmm) = puu(ji,jj,jk,Kmm) + rn_atfp * ( puu(ji,jj,jk,Kbb) - 2._wp * puu(ji,jj,jk,Kmm) + puu(ji,jj,jk,Kaa) ) |
---|
| 241 | pvv(ji,jj,jk,Kmm) = pvv(ji,jj,jk,Kmm) + rn_atfp * ( pvv(ji,jj,jk,Kbb) - 2._wp * pvv(ji,jj,jk,Kmm) + pvv(ji,jj,jk,Kaa) ) |
---|
[12340] | 242 | END_3D |
---|
[2528] | 243 | ! |
---|
[6140] | 244 | ELSE ! Asselin filter applied on thickness weighted velocity |
---|
| 245 | ! |
---|
[9019] | 246 | ALLOCATE( ze3u_f(jpi,jpj,jpk) , ze3v_f(jpi,jpj,jpk) ) |
---|
[11050] | 247 | ! Now filtered scale factor at (u/v)-points stored in ze3u_f, ze3v_f |
---|
| 248 | CALL dom_vvl_interpol( pe3t(:,:,:,Kmm), ze3u_f, 'U' ) |
---|
| 249 | CALL dom_vvl_interpol( pe3t(:,:,:,Kmm), ze3v_f, 'V' ) |
---|
[13295] | 250 | DO_3D( 1, 1, 1, 1, 1, jpkm1 ) |
---|
[12340] | 251 | zue3a = pe3u(ji,jj,jk,Kaa) * puu(ji,jj,jk,Kaa) |
---|
| 252 | zve3a = pe3v(ji,jj,jk,Kaa) * pvv(ji,jj,jk,Kaa) |
---|
| 253 | zue3n = pe3u(ji,jj,jk,Kmm) * puu(ji,jj,jk,Kmm) |
---|
| 254 | zve3n = pe3v(ji,jj,jk,Kmm) * pvv(ji,jj,jk,Kmm) |
---|
| 255 | zue3b = pe3u(ji,jj,jk,Kbb) * puu(ji,jj,jk,Kbb) |
---|
| 256 | zve3b = pe3v(ji,jj,jk,Kbb) * pvv(ji,jj,jk,Kbb) |
---|
| 257 | ! |
---|
[12489] | 258 | puu(ji,jj,jk,Kmm) = ( zue3n + rn_atfp * ( zue3b - 2._wp * zue3n + zue3a ) ) / ze3u_f(ji,jj,jk) |
---|
| 259 | pvv(ji,jj,jk,Kmm) = ( zve3n + rn_atfp * ( zve3b - 2._wp * zve3n + zve3a ) ) / ze3v_f(ji,jj,jk) |
---|
[12340] | 260 | END_3D |
---|
[14072] | 261 | pe3u(:,:,1:jpkm1,Kmm) = ze3u_f(:,:,1:jpkm1) |
---|
[11050] | 262 | pe3v(:,:,1:jpkm1,Kmm) = ze3v_f(:,:,1:jpkm1) |
---|
[6140] | 263 | ! |
---|
[9019] | 264 | DEALLOCATE( ze3u_f , ze3v_f ) |
---|
[2528] | 265 | ENDIF |
---|
| 266 | ! |
---|
[12372] | 267 | DEALLOCATE( ze3t_f, zwfld ) |
---|
[3] | 268 | ENDIF |
---|
[2528] | 269 | ! |
---|
[6140] | 270 | IF( ln_dynspg_ts .AND. ln_bt_fw ) THEN |
---|
[11050] | 271 | ! Revert filtered "now" velocities to time split estimate |
---|
[14072] | 272 | ! Doing it here also means that asselin filter contribution is removed |
---|
[11050] | 273 | zue(:,:) = pe3u(:,:,1,Kmm) * puu(:,:,1,Kmm) * umask(:,:,1) |
---|
[14072] | 274 | zve(:,:) = pe3v(:,:,1,Kmm) * pvv(:,:,1,Kmm) * vmask(:,:,1) |
---|
[4990] | 275 | DO jk = 2, jpkm1 |
---|
[11050] | 276 | zue(:,:) = zue(:,:) + pe3u(:,:,jk,Kmm) * puu(:,:,jk,Kmm) * umask(:,:,jk) |
---|
[14072] | 277 | zve(:,:) = zve(:,:) + pe3v(:,:,jk,Kmm) * pvv(:,:,jk,Kmm) * vmask(:,:,jk) |
---|
[4370] | 278 | END DO |
---|
| 279 | DO jk = 1, jpkm1 |
---|
[11050] | 280 | puu(:,:,jk,Kmm) = puu(:,:,jk,Kmm) - (zue(:,:) * r1_hu(:,:,Kmm) - uu_b(:,:,Kmm)) * umask(:,:,jk) |
---|
| 281 | pvv(:,:,jk,Kmm) = pvv(:,:,jk,Kmm) - (zve(:,:) * r1_hv(:,:,Kmm) - vv_b(:,:,Kmm)) * vmask(:,:,jk) |
---|
[4292] | 282 | END DO |
---|
| 283 | ENDIF |
---|
| 284 | ! |
---|
[12489] | 285 | ENDIF ! .NOT. l_1st_euler |
---|
[4354] | 286 | ! |
---|
| 287 | ! Set "now" and "before" barotropic velocities for next time step: |
---|
| 288 | ! JC: Would be more clever to swap variables than to make a full vertical |
---|
| 289 | ! integration |
---|
| 290 | ! |
---|
[6140] | 291 | IF(.NOT.ln_linssh ) THEN |
---|
[11050] | 292 | hu(:,:,Kmm) = pe3u(:,:,1,Kmm ) * umask(:,:,1) |
---|
| 293 | hv(:,:,Kmm) = pe3v(:,:,1,Kmm ) * vmask(:,:,1) |
---|
[6140] | 294 | DO jk = 2, jpkm1 |
---|
[11050] | 295 | hu(:,:,Kmm) = hu(:,:,Kmm) + pe3u(:,:,jk,Kmm ) * umask(:,:,jk) |
---|
| 296 | hv(:,:,Kmm) = hv(:,:,Kmm) + pe3v(:,:,jk,Kmm ) * vmask(:,:,jk) |
---|
[4354] | 297 | END DO |
---|
[11050] | 298 | r1_hu(:,:,Kmm) = ssumask(:,:) / ( hu(:,:,Kmm) + 1._wp - ssumask(:,:) ) |
---|
| 299 | r1_hv(:,:,Kmm) = ssvmask(:,:) / ( hv(:,:,Kmm) + 1._wp - ssvmask(:,:) ) |
---|
[4354] | 300 | ENDIF |
---|
| 301 | ! |
---|
[11050] | 302 | uu_b(:,:,Kaa) = pe3u(:,:,1,Kaa) * puu(:,:,1,Kaa) * umask(:,:,1) |
---|
| 303 | uu_b(:,:,Kmm) = pe3u(:,:,1,Kmm) * puu(:,:,1,Kmm) * umask(:,:,1) |
---|
| 304 | vv_b(:,:,Kaa) = pe3v(:,:,1,Kaa) * pvv(:,:,1,Kaa) * vmask(:,:,1) |
---|
| 305 | vv_b(:,:,Kmm) = pe3v(:,:,1,Kmm) * pvv(:,:,1,Kmm) * vmask(:,:,1) |
---|
[6140] | 306 | DO jk = 2, jpkm1 |
---|
[11050] | 307 | uu_b(:,:,Kaa) = uu_b(:,:,Kaa) + pe3u(:,:,jk,Kaa) * puu(:,:,jk,Kaa) * umask(:,:,jk) |
---|
| 308 | uu_b(:,:,Kmm) = uu_b(:,:,Kmm) + pe3u(:,:,jk,Kmm) * puu(:,:,jk,Kmm) * umask(:,:,jk) |
---|
| 309 | vv_b(:,:,Kaa) = vv_b(:,:,Kaa) + pe3v(:,:,jk,Kaa) * pvv(:,:,jk,Kaa) * vmask(:,:,jk) |
---|
| 310 | vv_b(:,:,Kmm) = vv_b(:,:,Kmm) + pe3v(:,:,jk,Kmm) * pvv(:,:,jk,Kmm) * vmask(:,:,jk) |
---|
[4354] | 311 | END DO |
---|
[11050] | 312 | uu_b(:,:,Kaa) = uu_b(:,:,Kaa) * r1_hu(:,:,Kaa) |
---|
| 313 | vv_b(:,:,Kaa) = vv_b(:,:,Kaa) * r1_hv(:,:,Kaa) |
---|
| 314 | uu_b(:,:,Kmm) = uu_b(:,:,Kmm) * r1_hu(:,:,Kmm) |
---|
| 315 | vv_b(:,:,Kmm) = vv_b(:,:,Kmm) * r1_hv(:,:,Kmm) |
---|
[4354] | 316 | ! |
---|
[6140] | 317 | IF( .NOT.ln_dynspg_ts ) THEN ! output the barotropic currents |
---|
[11050] | 318 | CALL iom_put( "ubar", uu_b(:,:,Kmm) ) |
---|
| 319 | CALL iom_put( "vbar", vv_b(:,:,Kmm) ) |
---|
[6140] | 320 | ENDIF |
---|
[4990] | 321 | IF( l_trddyn ) THEN ! 3D output: asselin filter trends on momentum |
---|
[11050] | 322 | zua(:,:,:) = ( puu(:,:,:,Kmm) - zua(:,:,:) ) * z1_2dt |
---|
| 323 | zva(:,:,:) = ( pvv(:,:,:,Kmm) - zva(:,:,:) ) * z1_2dt |
---|
[10946] | 324 | CALL trd_dyn( zua, zva, jpdyn_atf, kt, Kmm ) |
---|
[4990] | 325 | ENDIF |
---|
| 326 | ! |
---|
[13472] | 327 | IF ( iom_use("utau") ) THEN |
---|
| 328 | IF ( ln_drgice_imp.OR.ln_isfcav ) THEN |
---|
[14072] | 329 | ALLOCATE(zutau(jpi,jpj)) |
---|
[13472] | 330 | DO_2D( 0, 0, 0, 0 ) |
---|
[14072] | 331 | jk = miku(ji,jj) |
---|
[13472] | 332 | zutau(ji,jj) = utau(ji,jj) + 0.5_wp * rho0 * ( rCdU_top(ji+1,jj)+rCdU_top(ji,jj) ) * puu(ji,jj,jk,Kaa) |
---|
| 333 | END_2D |
---|
| 334 | CALL iom_put( "utau", zutau(:,:) ) |
---|
| 335 | DEALLOCATE(zutau) |
---|
| 336 | ELSE |
---|
| 337 | CALL iom_put( "utau", utau(:,:) ) |
---|
| 338 | ENDIF |
---|
| 339 | ENDIF |
---|
| 340 | ! |
---|
| 341 | IF ( iom_use("vtau") ) THEN |
---|
| 342 | IF ( ln_drgice_imp.OR.ln_isfcav ) THEN |
---|
| 343 | ALLOCATE(zvtau(jpi,jpj)) |
---|
| 344 | DO_2D( 0, 0, 0, 0 ) |
---|
| 345 | jk = mikv(ji,jj) |
---|
| 346 | zvtau(ji,jj) = vtau(ji,jj) + 0.5_wp * rho0 * ( rCdU_top(ji,jj+1)+rCdU_top(ji,jj) ) * pvv(ji,jj,jk,Kaa) |
---|
| 347 | END_2D |
---|
| 348 | CALL iom_put( "vtau", zvtau(:,:) ) |
---|
| 349 | DEALLOCATE(zvtau) |
---|
| 350 | ELSE |
---|
| 351 | CALL iom_put( "vtau", vtau(:,:) ) |
---|
| 352 | ENDIF |
---|
| 353 | ENDIF |
---|
| 354 | ! |
---|
[12236] | 355 | IF(sn_cfctl%l_prtctl) CALL prt_ctl( tab3d_1=puu(:,:,:,Kaa), clinfo1=' nxt - puu(:,:,:,Kaa): ', mask1=umask, & |
---|
| 356 | & tab3d_2=pvv(:,:,:,Kaa), clinfo2=' pvv(:,:,:,Kaa): ' , mask2=vmask ) |
---|
[14072] | 357 | ! |
---|
[9019] | 358 | IF( ln_dynspg_ts ) DEALLOCATE( zue, zve ) |
---|
| 359 | IF( l_trddyn ) DEALLOCATE( zua, zva ) |
---|
[11050] | 360 | IF( ln_timing ) CALL timing_stop('dyn_atf') |
---|
[2715] | 361 | ! |
---|
[11050] | 362 | END SUBROUTINE dyn_atf |
---|
[3] | 363 | |
---|
[13237] | 364 | #endif |
---|
| 365 | |
---|
[1502] | 366 | !!========================================================================= |
---|
[11050] | 367 | END MODULE dynatf |
---|