Changeset 791 for branches/dev_001_GM/NEMO/OPA_SRC/TRA/traadv_ubs.F90
- Timestamp:
- 2008-01-12T21:33:34+01:00 (16 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/dev_001_GM/NEMO/OPA_SRC/TRA/traadv_ubs.F90
r786 r791 4 4 !! Ocean tracers: horizontal & vertical advective trend 5 5 !!============================================================================== 6 !! History : 1.0 ! 06-08 (L. Debreu, R. Benshila) Original code7 !! 2.4 ! 08-01 (G. Madec) Merge TRA-TRC6 !! History : 1.0 ! 2006-08 (L. Debreu, R. Benshila) Original code 7 !! 2.4 ! 2008-01 (G. Madec) merge TRC-TRA + switch from velocity to transport 8 8 !!---------------------------------------------------------------------- 9 9 … … 27 27 PUBLIC tra_adv_ubs ! routine called by traadv module 28 28 29 REAL(wp), DIMENSION(jpi,jpj) :: e1e2tr ! = 1/(e1t * e2t)30 31 29 !! * Substitutions 32 30 # include "domzgr_substitute.h90" … … 34 32 !!---------------------------------------------------------------------- 35 33 !! NEMO/OPA & TRP 2.4 , LOCEAN-IPSL (2008) 36 !! $Id :$34 !! $Id$ 37 35 !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) 38 36 !!---------------------------------------------------------------------- … … 71 69 !! 72 70 !! Reference : Shchepetkin, A. F., J. C. McWilliams, 2005, Ocean Modelling, 9, 347-404. 73 !! Farrow, D.E., Stevens, D.P., 1995, J. Phys. Ocean. 25, 1731 –1741.71 !! Farrow, D.E., Stevens, D.P., 1995, J. Phys. Ocean. 25, 1731-1741. 74 72 !!---------------------------------------------------------------------- 75 73 INTEGER , INTENT(in ) :: kt ! ocean time-step index … … 80 78 REAL(wp) , INTENT(inout), DIMENSION(jpi,jpj,jpk) :: pta ! tracer trend 81 79 !! 82 INTEGER :: ji, jj, jk ! dummy loop indices 83 REAL(wp) :: zta, zbtr, zcoef ! temporary scalars 84 REAL(wp) :: zfui, zfp_ui, zfm_ui, zcenut ! " " 85 REAL(wp) :: zfvj, zfp_vj, zfm_vj, zcenvt ! " " 86 REAL(wp) :: z_hdivn ! " " 87 REAL(wp), DIMENSION(jpi,jpj) :: zeeu, zeev ! temporary 2D workspace 88 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zwx, zwy ! temporary 3D workspace 89 REAL(wp), DIMENSION(jpi,jpj,jpk) :: ztu , ztv , zltu , zltv, ztrdt ! " " 90 !!---------------------------------------------------------------------- 91 80 INTEGER :: ji, jj, jk ! dummy loop indices 81 REAL(wp) :: zbtr, zcoef, z_hdivn ! temporary scalars 82 REAL(wp) :: zeeu, zfp_ui, zfm_ui, zcenut ! " " 83 REAL(wp) :: zeev, zfp_vj, zfm_vj, zcenvt ! " " 84 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zwx, ztu, zltu ! temporary 3D workspace 85 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zwy, ztv, zltv ! " " 86 !!---------------------------------------------------------------------- 87 88 !!gm this can be optimized: we don't need zero everywhere 92 89 zltu(:,:,:) = 0.e0 93 90 zltv(:,:,:) = 0.e0 … … 97 94 IF(lwp) WRITE(numout,*) 'tra_adv_ubs : horizontal UBS advection scheme' 98 95 IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 99 ! 100 e1e2tr(:,:) = 1. / ( e1t(:,:) * e2t(:,:) ) 101 ENDIF 102 103 ! store pta trends 104 ztrdt(:,:,:) = pta(:,:,:) 105 106 zcoef = 1./6. 107 ! ! =============== 108 DO jk = 1, jpkm1 ! Horizontal slab 109 ! ! =============== 110 111 ! Initialization of metric arrays (for z- or s-coordinates) 96 ENDIF 97 98 ! Laplacian 99 DO jk = 1, jpkm1 100 DO jj = 1, jpjm1 ! First derivative (gradient) 101 DO ji = 1, fs_jpim1 ! vector opt. 102 zeeu = e2u(ji,jj) * fse3u(ji,jj,jk) / e1u(ji,jj) * umask(ji,jj,jk) 103 zeev = e1v(ji,jj) * fse3v(ji,jj,jk) / e2v(ji,jj) * vmask(ji,jj,jk) 104 ztu(ji,jj,jk) = zeeu * ( ptb(ji+1,jj ,jk) - ptb(ji,jj,jk) ) 105 ztv(ji,jj,jk) = zeev * ( ptb(ji ,jj+1,jk) - ptb(ji,jj,jk) ) 106 END DO 107 END DO 108 DO jj = 2, jpjm1 ! Second derivative (divergence) 109 DO ji = fs_2, fs_jpim1 ! vector opt. 110 zcoef = 1. / ( 6. * fse3t(ji,jj,jk) ) 111 zltu(ji,jj,jk) = ( ztu(ji,jj,jk) - ztu(ji-1,jj,jk) ) * zcoef 112 zltv(ji,jj,jk) = ( ztv(ji,jj,jk) - ztv(ji,jj-1,jk) ) * zcoef 113 END DO 114 END DO 115 END DO 116 117 ! Lateral boundary conditions on the laplacian (zlt,zls) (unchanged sgn) 118 CALL lbc_lnk( zltu, 'T', 1. ) ; CALL lbc_lnk( zltv, 'T', 1. ) 119 120 121 ! Horizontal advective fluxes 122 DO jk = 1, jpkm1 112 123 DO jj = 1, jpjm1 113 124 DO ji = 1, fs_jpim1 ! vector opt. 114 #if defined key_zco 115 ! z-coordinates, no vertical scale factors 116 zeeu(ji,jj) = e2u(ji,jj) / e1u(ji,jj) * umask(ji,jj,jk) 117 zeev(ji,jj) = e1v(ji,jj) / e2v(ji,jj) * vmask(ji,jj,jk) 118 #else 119 ! s-coordinates, vertical scale factor are used 120 zeeu(ji,jj) = e2u(ji,jj) * fse3u(ji,jj,jk) / e1u(ji,jj) * umask(ji,jj,jk) 121 zeev(ji,jj) = e1v(ji,jj) * fse3v(ji,jj,jk) / e2v(ji,jj) * vmask(ji,jj,jk) 122 #endif 123 END DO 124 END DO 125 126 ! Laplacian 127 ! First derivative (gradient) 128 DO jj = 1, jpjm1 129 DO ji = 1, fs_jpim1 ! vector opt. 130 ztu(ji,jj,jk) = zeeu(ji,jj) * ( ptb(ji+1,jj ,jk) - ptb(ji,jj,jk) ) 131 ztv(ji,jj,jk) = zeev(ji,jj) * ( ptb(ji ,jj+1,jk) - ptb(ji,jj,jk) ) 132 END DO 133 END DO 134 ! Second derivative (divergence) 135 DO jj = 2, jpjm1 136 DO ji = fs_2, fs_jpim1 ! vector opt. 137 #if ! defined key_zco 138 zcoef = 1. / ( 6. * fse3t(ji,jj,jk) ) 139 #endif 140 zltu(ji,jj,jk) = ( ztu(ji,jj,jk) - ztu(ji-1,jj,jk) ) * zcoef 141 zltv(ji,jj,jk) = ( ztv(ji,jj,jk) - ztv(ji,jj-1,jk) ) * zcoef 142 END DO 143 END DO 144 ! ! ================= 145 END DO ! End of slab 146 ! ! ================= 147 148 ! Lateral boundary conditions on the laplacian (zlt,zls) (unchanged sgn) 149 CALL lbc_lnk( zltu, 'T', 1. ) ; CALL lbc_lnk( zltv, 'T', 1. ) 150 151 ! ! =============== 152 DO jk = 1, jpkm1 ! Horizontal slab 153 ! ! =============== 154 ! Horizontal advective fluxes 155 DO jj = 1, jpjm1 156 DO ji = 1, fs_jpim1 ! vector opt. 157 ! volume fluxes * 1/2 158 #if defined key_zco 159 zfui = 0.5 * e2u(ji,jj) * pun(ji,jj,jk) 160 zfvj = 0.5 * e1v(ji,jj) * pvn(ji,jj,jk) 161 #else 162 zfui = 0.5 * e2u(ji,jj) * fse3u(ji,jj,jk) * pun(ji,jj,jk) 163 zfvj = 0.5 * e1v(ji,jj) * fse3v(ji,jj,jk) * pvn(ji,jj,jk) 164 #endif 165 ! upstream scheme 166 zfp_ui = zfui + ABS( zfui ) 167 zfp_vj = zfvj + ABS( zfvj ) 168 zfm_ui = zfui - ABS( zfui ) 169 zfm_vj = zfvj - ABS( zfvj ) 125 ! upstream transport 126 zfp_ui = pun(ji,jj,jk) + ABS( pun(ji,jj,jk) ) ! = 2.pun if pun>0, =0 if pun<0 127 zfm_ui = pun(ji,jj,jk) - ABS( pun(ji,jj,jk) ) ! = 0 if pun>0, = 2.pun if pun<0 128 zfp_vj = pvn(ji,jj,jk) + ABS( pvn(ji,jj,jk) ) ! idem on pvn 129 zfm_vj = pvn(ji,jj,jk) - ABS( pvn(ji,jj,jk) ) ! 170 130 ! centered scheme 171 zcenut = zfui * ( ptn(ji,jj,jk) + ptn(ji+1,jj ,jk) ) 172 zcenvt = zfvj * ( ptn(ji,jj,jk) + ptn(ji ,jj+1,jk) ) 173 ! mixed centered / upstream scheme 174 zwx(ji,jj,jk) = zcenut - zfp_ui * zltu(ji,jj,jk) -zfm_ui * zltu(ji+1,jj,jk) 175 zwy(ji,jj,jk) = zcenvt - zfp_vj * zltv(ji,jj,jk) -zfm_vj * zltv(ji,jj+1,jk) 176 END DO 177 END DO 178 179 ! Tracer flux divergence at t-point added to the general trend 180 DO jj = 2, jpjm1 181 DO ji = fs_2, fs_jpim1 ! vector opt. 182 ! horizontal advective trends 183 #if defined key_zco 184 zbtr = e1e2tr(ji,jj) 185 #else 186 zbtr = e1e2tr(ji,jj) / fse3t(ji,jj,jk) 187 #endif 188 zta = - zbtr * ( zwx(ji,jj,jk) - zwx(ji-1,jj ,jk) & 189 & + zwy(ji,jj,jk) - zwy(ji ,jj-1,jk) ) 190 ! add it to the general tracer trends 191 pta(ji,jj,jk) = pta(ji,jj,jk) + zta 192 END DO 193 END DO 194 ! ! =============== 195 END DO ! End of slab 196 ! ! =============== 197 198 ! Horizontal trend used in tra_adv_ztvd subroutine 199 zltu(:,:,:) = pta(:,:,:) - ztrdt(:,:,:) 200 201 ! Save the horizontal advective trends for diagnostic 202 ! ----------------------------------------------------- 203 IF( l_trdtra ) THEN 131 zcenut = pun(ji,jj,jk) * ( ptn(ji,jj,jk) + ptn(ji+1,jj ,jk) ) 132 zcenvt = pvn(ji,jj,jk) * ( ptn(ji,jj,jk) + ptn(ji ,jj+1,jk) ) 133 ! UBS scheme 134 zwx(ji,jj,jk) = 0.5 * ( zcenut - zfp_ui * zltu(ji,jj,jk) -zfm_ui * zltu(ji+1,jj,jk) ) 135 zwy(ji,jj,jk) = 0.5 * ( zcenvt - zfp_vj * zltv(ji,jj,jk) -zfm_vj * zltv(ji,jj+1,jk) ) 136 END DO 137 END DO 138 END DO 139 140 zltu(:,:,:) = pta(:,:,:) ! store pta trends 141 142 ! Horizontal advective trends 143 DO jk = 1, jpkm1 144 DO jj = 2, jpjm1 145 DO ji = fs_2, fs_jpim1 ! vector opt. 146 zbtr = 1.e0 / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 147 ! horizontal advective trends added to the general tracer trends 148 pta(ji,jj,jk) = pta(ji,jj,jk) - zbtr * ( zwx(ji,jj,jk) - zwx(ji-1,jj ,jk) & 149 & + zwy(ji,jj,jk) - zwy(ji ,jj-1,jk) ) 150 END DO 151 END DO 152 END DO 153 154 zltu(:,:,:) = pta(:,:,:) - zltu(:,:,:) ! store the Horizontal advective trend (used in tra_adv_ztvd subroutine) 155 156 157 IF( l_trdtra ) THEN ! trend diagnostics (from advective fluxes) 204 158 CALL trd_tra_adv( kt, ktra, jpt_trd_xad, cdtype, zwx, pun, ptn ) 205 159 CALL trd_tra_adv( kt, ktra, jpt_trd_yad, cdtype, zwy, pvn, ptn ) 206 160 ENDIF 207 208 ! "Poleward" heat or salt transport 161 ! ! "Poleward" heat or salt transports 209 162 IF( cdtype == 'TRA' .AND. ln_diaptr .AND. ( MOD( kt, nf_ptr ) == 0 ) ) THEN 210 IF( lk_zco ) THEN 211 DO jk = 1, jpkm1 212 DO jj = 2, jpjm1 213 DO ji = fs_2, fs_jpim1 ! vector opt. 214 zwy(ji,jj,jk) = zwy(ji,jj,jk) * fse3v(ji,jj,jk) 215 END DO 216 END DO 217 END DO 218 ENDIF 219 IF( ktra == jp_tem) pht_adv(:) = ptr_vj( zwy(:,:,:) ) 220 IF( ktra == jp_sal) pst_adv(:) = ptr_vj( zwy(:,:,:) ) 221 ENDIF 222 163 IF( ktra == jp_tem ) pht_adv(:) = ptr_vj( zwy(:,:,:) ) 164 IF( ktra == jp_sal ) pst_adv(:) = ptr_vj( zwy(:,:,:) ) 165 ENDIF 166 ! ! control print 223 167 IF(ln_ctl) CALL prt_ctl( tab3d_1=pta, clinfo1=' ubs - had: ', mask1=tmask, clinfo3=cdtype ) 224 168 … … 226 170 ! II. Vertical advection 227 171 ! ---------------------- 228 IF( l_trdtra ) z trdt(:,:,:) = pta(:,:,:) ! Save ta and sa trends172 IF( l_trdtra ) zltv(:,:,:) = pta(:,:,:) ! store pta if trend diag. 229 173 230 174 ! TVD scheme the vertical direction 231 175 CALL tra_adv_ztvd( kt, pwn, zltu, ptb, ptn, pta ) 232 176 233 IF( l_trdtra ) THEN !vertical advective trend diagnostics234 DO jk = 1, jpkm1 177 IF( l_trdtra ) THEN ! vertical advective trend diagnostics 178 DO jk = 1, jpkm1 ! (compute -w.dk[ptn]= -dk[w.ptn] + ptn.dk[w]) 235 179 DO jj = 2, jpjm1 236 180 DO ji = fs_2, fs_jpim1 ! vector opt. 237 z_hdivn = ( pwn(ji,jj,jk) - pwn(ji,jj,jk+1) ) / fse3t(ji,jj,jk) 238 ztrdt(ji,jj,jk) = pta(ji,jj,jk) - ztrdt(ji,jj,jk+1) + ptn(ji,jj,jk) * z_hdivn 181 zbtr = 1.e0 / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 182 z_hdivn = ( pwn(ji,jj,jk) - pwn(ji,jj,jk+1) ) * zbtr 183 zltv(ji,jj,jk) = pta(ji,jj,jk) - zltv(ji,jj,jk+1) + ptn(ji,jj,jk) * z_hdivn 239 184 END DO 240 185 END DO 241 186 END DO 242 CALL trd_tra( kt, ktra, jpt_trd_zad, cdtype, ptrd3d=ztrdt ) 243 ! 244 ENDIF 245 187 CALL trd_tra( kt, ktra, jpt_trd_zad, cdtype, ptrd3d=zltv ) 188 ENDIF 189 ! ! control print 246 190 IF(ln_ctl) CALL prt_ctl( tab3d_1=pta, clinfo1=' ubs - zad: ', mask1=tmask, clinfo3=cdtype ) 247 191 ! … … 249 193 250 194 251 SUBROUTINE tra_adv_ztvd( kt, pwn, zttrd, ptb, ptn, pta )195 SUBROUTINE tra_adv_ztvd( kt, pwn, phtrd, ptb, ptn, pta ) 252 196 !!---------------------------------------------------------------------- 253 197 !! *** ROUTINE tra_adv_ztvd *** 254 198 !! 255 !! ** Purpose : Compute the now trend due to total advection of256 !! tracers and add it to the general trend of tracer equations199 !! ** Purpose : Compute the vertical advective trend of a tracer 200 !! and add it to its general trend 257 201 !! 258 202 !! ** Method : TVD scheme, i.e. 2nd order centered scheme with 259 203 !! corrected flux (monotonic correction) 260 !! note: - this advection scheme needs a leap-frog time scheme 261 !! 262 !! ** Action : - update (ta,sa) with the now advective tracer trends 263 !! - save the trends in (ztrdt,ztrds) ('key_trdtra') 264 !!---------------------------------------------------------------------- 265 INTEGER , INTENT(in ) :: kt ! ocean time-step 266 REAL(wp), INTENT(in ), DIMENSION(jpi,jpj,jpk) :: pwn ! verical effective velocity 267 REAL(wp), INTENT(in ), DIMENSION(jpi,jpj,jpk) :: zttrd ! lateral advective trends on T & S 268 REAL(wp), INTENT(inout), DIMENSION(jpi,jpj,jpk) :: ptb, ptn ! before and now tracer fields 269 REAL(wp), INTENT(inout), DIMENSION(jpi,jpj,jpk) :: pta ! tracer trend 204 !! 205 !! ** Action : - update pta with the vertical advective trend 206 !! - trend diagnostic (lk_trdtra=T) 207 !!---------------------------------------------------------------------- 208 INTEGER , INTENT(in ) :: kt ! ocean time-step 209 REAL(wp), INTENT(in ), DIMENSION(jpi,jpj,jpk) :: pwn ! verical effective velocity 210 REAL(wp), INTENT(in ), DIMENSION(jpi,jpj,jpk) :: phtrd ! lateral advective trends on T & S 211 REAL(wp), INTENT(inout), DIMENSION(jpi,jpj,jpk) :: ptb, ptn ! before and now tracer fields 212 REAL(wp), INTENT(inout), DIMENSION(jpi,jpj,jpk) :: pta ! tracer trend 270 213 !! 271 214 INTEGER :: ji, jj, jk ! dummy loop indices 272 REAL(wp) :: z2dtt, zbtr, z ew, z2! temporary scalar273 REAL(wp) :: ztak, zfp_wk, zfm_wk 215 REAL(wp) :: z2dtt, zbtr, z2 ! temporary scalar 216 REAL(wp) :: ztak, zfp_wk, zfm_wk ! " " 274 217 REAL(wp), DIMENSION (jpi,jpj,jpk) :: zti, ztw ! temporary 3D workspace 275 218 !!---------------------------------------------------------------------- … … 285 228 ENDIF 286 229 287 ! Bottom value : flux set to zero288 ! --------------289 ztw(:,:,jpk) = 0.e0 ; zti (:,:,:) = 0.e0290 291 230 292 231 ! upstream advection with initial mass fluxes & intermediate update 293 ! ------------------------------------------------------------------- 294 ! Surface value 295 IF( lk_dynspg_rl .OR. lk_vvl ) THEN ! rigid lid : flux set to zero 296 ztw(:,:,1) = 0.e0 297 ELSE ! free surface 298 ztw(:,:,1) = e1t(:,:) * e2t(:,:) * pwn(:,:,1) * ptb(:,:,1) 299 ENDIF 300 301 ! Interior value 302 DO jk = 2, jpkm1 232 ztw(:,:,jpk) = 0.e0 ! Bottom value : flux set to zero 233 zti(:,:,jpk) = 0.e0 234 ! ! Surface value : 235 IF( lk_dynspg_rl .OR. lk_vvl ) THEN ; ztw(:,:, 1 ) = 0.e0 ! rigid lid or non-linear fs 236 ELSE ; ztw(:,:, 1 ) = pwn(:,:,1) * ptb(:,:,1) ! linear free surface 237 ENDIF 238 ! 239 DO jk = 2, jpkm1 ! interior values 303 240 DO jj = 1, jpj 304 241 DO ji = 1, jpi 305 zew = 0.5 * e1t(ji,jj) * e2t(ji,jj) * pwn(ji,jj,jk) 306 zfp_wk = zew + ABS( zew ) 307 zfm_wk = zew - ABS( zew ) 308 ztw(ji,jj,jk) = zfp_wk * ptb(ji,jj,jk) + zfm_wk * ptb(ji,jj,jk-1) 242 zfp_wk = pwn(ji,jj,jk) + ABS( pwn(ji,jj,jk) ) 243 zfm_wk = pwn(ji,jj,jk) - ABS( pwn(ji,jj,jk) ) 244 ztw(ji,jj,jk) = 0.5 * ( zfp_wk * ptb(ji,jj,jk) + zfm_wk * ptb(ji,jj,jk-1) ) 309 245 END DO 310 246 END DO … … 318 254 zbtr = 1./ ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 319 255 ztak = - ( ztw(ji,jj,jk) - ztw(ji,jj,jk+1) ) * zbtr 320 pta(ji,jj,jk) = pta(ji,jj,jk) + ztak 321 zti (ji,jj,jk) = ( ptb(ji,jj,jk) + z2dtt * ( ztak + zttrd(ji,jj,jk) ) ) * tmask(ji,jj,jk) 322 END DO 323 END DO 324 END DO 325 326 ! Lateral boundary conditions on zti, zsi (unchanged sign) 327 CALL lbc_lnk( zti, 'T', 1. ) 256 pta(ji,jj,jk) = pta(ji,jj,jk) + ztak 257 zti(ji,jj,jk) = ( ptb(ji,jj,jk) + z2dtt * ( ztak + phtrd(ji,jj,jk) ) ) * tmask(ji,jj,jk) 258 END DO 259 END DO 260 END DO 261 ! 262 CALL lbc_lnk( zti, 'T', 1. ) ! Lateral boundary conditions on zti, zsi (unchanged sign) 328 263 329 264 330 265 ! antidiffusive flux : high order minus low order 331 ! -------------------------------------------------332 266 ztw(:,:,1) = 0.e0 ! Surface value 333 334 267 DO jk = 2, jpkm1 ! Interior value 335 268 DO jj = 1, jpj 336 269 DO ji = 1, jpi 337 zew = 0.5 * e1t(ji,jj) * e2t(ji,jj) * pwn(ji,jj,jk) 338 ztw(ji,jj,jk) = zew * ( ptn(ji,jj,jk) + ptn(ji,jj,jk-1) ) - ztw(ji,jj,jk) 339 END DO 340 END DO 341 END DO 342 343 ! monotonicity algorithm 344 ! ------------------------ 345 CALL nonosc_z( ptb, ztw, zti, z2 ) 270 ztw(ji,jj,jk) = 0.5 * pwn(ji,jj,jk) * ( ptn(ji,jj,jk) + ptn(ji,jj,jk-1) ) - ztw(ji,jj,jk) 271 END DO 272 END DO 273 END DO 274 ! 275 CALL nonosc_z( ptb, ztw, zti, z2 ) ! monotonicity algorithm 346 276 347 277 348 278 ! final trend with corrected fluxes 349 ! -----------------------------------350 279 DO jk = 1, jpkm1 351 280 DO jj = 2, jpjm1 352 281 DO ji = fs_2, fs_jpim1 ! vector opt. 353 282 zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) 354 ! k- vertical advective trends 355 ztak = - ( ztw(ji,jj,jk) - ztw(ji,jj,jk+1) ) * zbtr 356 ! add them to the general tracer trends 357 pta(ji,jj,jk) = pta(ji,jj,jk) + ztak 283 ! k- vertical advective trends added to the general tracer trends 284 pta(ji,jj,jk) = pta(ji,jj,jk) - ( ztw(ji,jj,jk) - ztw(ji,jj,jk+1) ) * zbtr 358 285 END DO 359 286 END DO
Note: See TracChangeset
for help on using the changeset viewer.