ROMS
Loading...
Searching...
No Matches
t3dmix2_mod Module Reference

Functions/Subroutines

subroutine, public t3dmix2 (ng, tile)
 
subroutine t3dmix2_geo_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nstp, nnew, umask, vmask, umask_wet, vmask_wet, om_v, on_u, pm, pn, hz, z_r, diff3d_r, diff2, tclm, diatwrk, t)
 
subroutine t3dmix2_iso_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nstp, nnew, umask, vmask, umask_wet, vmask_wet, om_v, on_u, pm, pn, hz, z_r, diff3d_r, diff2, pden, tclm, diatwrk, t)
 
subroutine t3dmix2_s_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nstp, nnew, umask, vmask, umask_wet, vmask_wet, hz, pmon_u, pnom_v, pm, pn, diff3d_r, diff2, tclm, diatwrk, t)
 

Function/Subroutine Documentation

◆ t3dmix2()

subroutine public t3dmix2_mod::t3dmix2 ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 23 of file t3dmix2_geo.h.

24!***********************************************************************
25!
26 USE mod_param
27#ifdef TS_MIX_CLIMA
28 USE mod_clima
29#endif
30#ifdef DIAGNOSTICS_TS
31 USE mod_diags
32#endif
33 USE mod_grid
34 USE mod_mixing
35 USE mod_ocean
36 USE mod_stepping
37!
38! Imported variable declarations.
39!
40 integer, intent(in) :: ng, tile
41!
42! Local variable declarations.
43!
44 character (len=*), parameter :: MyFile = &
45 & __FILE__
46!
47#include "tile.h"
48!
49#ifdef PROFILE
50 CALL wclock_on (ng, inlm, 25, __line__, myfile)
51#endif
52 CALL t3dmix2_geo_tile (ng, tile, &
53 & lbi, ubi, lbj, ubj, &
54 & imins, imaxs, jmins, jmaxs, &
55 & nrhs(ng), nstp(ng), nnew(ng), &
56#ifdef MASKING
57 & grid(ng) % umask, &
58 & grid(ng) % vmask, &
59#endif
60#ifdef WET_DRY
61 & grid(ng) % umask_wet, &
62 & grid(ng) % vmask_wet, &
63#endif
64 & grid(ng) % om_v, &
65 & grid(ng) % on_u, &
66 & grid(ng) % pm, &
67 & grid(ng) % pn, &
68 & grid(ng) % Hz, &
69 & grid(ng) % z_r, &
70#ifdef DIFF_3DCOEF
71 & mixing(ng) % diff3d_r, &
72#else
73 & mixing(ng) % diff2, &
74#endif
75#ifdef TS_MIX_CLIMA
76 & clima(ng) % tclm, &
77#endif
78#ifdef DIAGNOSTICS_TS
79 & diags(ng) % DiaTwrk, &
80#endif
81 & ocean(ng) % t)
82#ifdef PROFILE
83 CALL wclock_off (ng, inlm, 25, __line__, myfile)
84#endif
85!
86 RETURN
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
type(t_diags), dimension(:), allocatable diags
Definition mod_diags.F:100
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_mixing), dimension(:), allocatable mixing
Definition mod_mixing.F:399
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
integer, dimension(:), allocatable nstp
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3

References mod_clima::clima, mod_diags::diags, mod_grid::grid, mod_param::inlm, mod_mixing::mixing, mod_stepping::nnew, mod_stepping::nrhs, mod_stepping::nstp, mod_ocean::ocean, t3dmix2_geo_tile(), wclock_off(), and wclock_on().

Referenced by rhs3d_mod::rhs3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ t3dmix2_geo_tile()

subroutine t3dmix2_mod::t3dmix2_geo_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nrhs,
integer, intent(in) nstp,
integer, intent(in) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) hz,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) z_r,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) diff3d_r,
real(r8), dimension(lbi:ubi,lbj:ubj,nt(ng)), intent(in) diff2,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),nt(ng)), intent(in) tclm,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),nt(ng), ndt), intent(inout) diatwrk,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),3,nt(ng)), intent(inout) t )
private

Definition at line 90 of file t3dmix2_geo.h.

114!***********************************************************************
115!
116 USE mod_param
117 USE mod_scalars
118!
119! Imported variable declarations.
120!
121 integer, intent(in) :: ng, tile
122 integer, intent(in) :: LBi, UBi, LBj, UBj
123 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
124 integer, intent(in) :: nrhs, nstp, nnew
125
126#ifdef ASSUMED_SHAPE
127# ifdef MASKING
128 real(r8), intent(in) :: umask(LBi:,LBj:)
129 real(r8), intent(in) :: vmask(LBi:,LBj:)
130# endif
131# ifdef WET_DRY
132 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
133 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
134# endif
135# ifdef DIFF_3DCOEF
136 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
137# else
138 real(r8), intent(in) :: diff2(LBi:,LBj:,:)
139# endif
140 real(r8), intent(in) :: om_v(LBi:,LBj:)
141 real(r8), intent(in) :: on_u(LBi:,LBj:)
142 real(r8), intent(in) :: pm(LBi:,LBj:)
143 real(r8), intent(in) :: pn(LBi:,LBj:)
144 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
145 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
146# ifdef TS_MIX_CLIMA
147 real(r8), intent(in) :: tclm(LBi:,LBj:,:,:)
148# endif
149# ifdef DIAGNOSTICS_TS
150 real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
151# endif
152 real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
153#else
154# ifdef MASKING
155 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
156 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
157# endif
158# ifdef WET_DRY
159 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
160 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
161# endif
162# ifdef DIFF_3DCOEF
163 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
164# else
165 real(r8), intent(in) :: diff2(LBi:UBi,LBj:UBj,NT(ng))
166# endif
167 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
168 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
169 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
170 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
171 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
172 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
173# ifdef TS_MIX_CLIMA
174 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
175# endif
176# ifdef DIAGNOSTICS_TS
177 real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
178 & NDT)
179# endif
180 real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
181#endif
182!
183! Local variable declarations.
184!
185 integer :: i, itrc, j, k, k1, k2
186
187 real(r8) :: cff, cff1, cff2, cff3, cff4
188
189 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE
190 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX
191
192 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: FS
193 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdz
194 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdx
195 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTde
196 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZdx
197 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZde
198
199#include "set_bounds.h"
200!
201!-----------------------------------------------------------------------
202! Compute horizontal harmonic diffusion along geopotential surfaces.
203!-----------------------------------------------------------------------
204!
205! Compute horizontal and vertical gradients. Notice the recursive
206! blocking sequence. The vertical placement of the gradients is:
207!
208! dTdx,dTde(:,:,k1) k rho-points
209! dTdx,dTde(:,:,k2) k+1 rho-points
210! FS,dTdz(:,:,k1) k-1/2 W-points
211! FS,dTdz(:,:,k2) k+1/2 W-points
212!
213#ifdef TS_MIX_STABILITY
214! In order to increase stability, the biharmonic operator is applied
215! as: 3/4 t(:,:,:,nrhs,:) + 1/4 t(:,:,:,nstp,:).
216!
217#endif
218
219 t_loop : DO itrc=1,nt(ng)
220 k2=1
221 k_loop : DO k=0,n(ng)
222 k1=k2
223 k2=3-k1
224 IF (k.lt.n(ng)) THEN
225 DO j=jstr,jend
226 DO i=istr,iend+1
227 cff=0.5_r8*(pm(i,j)+pm(i-1,j))
228#ifdef MASKING
229 cff=cff*umask(i,j)
230#endif
231#ifdef WET_DRY
232 cff=cff*umask_wet(i,j)
233#endif
234 dzdx(i,j,k2)=cff*(z_r(i ,j,k+1)- &
235 & z_r(i-1,j,k+1))
236#if defined TS_MIX_STABILITY
237 dtdx(i,j,k2)=cff*(0.75_r8*(t(i ,j,k+1,nrhs,itrc)- &
238 & t(i-1,j,k+1,nrhs,itrc))+ &
239 & 0.25_r8*(t(i ,j,k+1,nstp,itrc)- &
240 & t(i-1,j,k+1,nstp,itrc)))
241#elif defined TS_MIX_CLIMA
242 IF (ltracerclm(itrc,ng)) THEN
243 dtdx(i,j,k2)=cff*((t(i ,j,k+1,nrhs,itrc)- &
244 & tclm(i ,j,k+1,itrc))- &
245 & (t(i-1,j,k+1,nrhs,itrc)- &
246 & tclm(i-1,j,k+1,itrc)))
247 ELSE
248 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
249 & t(i-1,j,k+1,nrhs,itrc))
250 END IF
251#else
252 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
253 & t(i-1,j,k+1,nrhs,itrc))
254#endif
255 END DO
256 END DO
257 DO j=jstr,jend+1
258 DO i=istr,iend
259 cff=0.5_r8*(pn(i,j)+pn(i,j-1))
260#ifdef MASKING
261 cff=cff*vmask(i,j)
262#endif
263#ifdef WET_DRY
264 cff=cff*vmask_wet(i,j)
265#endif
266 dzde(i,j,k2)=cff*(z_r(i,j ,k+1)- &
267 & z_r(i,j-1,k+1))
268#if defined TS_MIX_STABILITY
269 dtde(i,j,k2)=cff*(0.75_r8*(t(i,j ,k+1,nrhs,itrc)- &
270 & t(i,j-1,k+1,nrhs,itrc))+ &
271 & 0.25_r8*(t(i,j ,k+1,nstp,itrc)- &
272 & t(i,j-1,k+1,nstp,itrc)))
273#elif defined TS_MIX_CLIMA
274 IF (ltracerclm(itrc,ng)) THEN
275 dtde(i,j,k2)=cff*((t(i,j ,k+1,nrhs,itrc)- &
276 & tclm(i,j ,k+1,itrc))- &
277 & (t(i,j-1,k+1,nrhs,itrc)- &
278 & tclm(i,j-1,k+1,itrc)))
279 ELSE
280 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
281 & t(i,j-1,k+1,nrhs,itrc))
282 END IF
283#else
284 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
285 & t(i,j-1,k+1,nrhs,itrc))
286#endif
287 END DO
288 END DO
289 END IF
290 IF ((k.eq.0).or.(k.eq.n(ng))) THEN
291 DO j=jstr-1,jend+1
292 DO i=istr-1,iend+1
293 dtdz(i,j,k2)=0.0_r8
294 fs(i,j,k2)=0.0_r8
295 END DO
296 END DO
297 ELSE
298 DO j=jstr-1,jend+1
299 DO i=istr-1,iend+1
300 cff=1.0_r8/(z_r(i,j,k+1)-z_r(i,j,k))
301#if defined TS_MIX_STABILITY
302 dtdz(i,j,k2)=cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
303 & t(i,j,k ,nrhs,itrc))+ &
304 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
305 & t(i,j,k ,nstp,itrc)))
306#elif defined TS_MIX_CLIMA
307 IF (ltracerclm(itrc,ng)) THEN
308 dtdz(i,j,k2)=cff*((t(i,j,k+1,nrhs,itrc)- &
309 & tclm(i,j,k+1,itrc))- &
310 & (t(i,j,k ,nrhs,itrc)- &
311 & tclm(i,j,k ,itrc)))
312 ELSE
313 dtdz(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
314 & t(i,j,k ,nrhs,itrc))
315 END IF
316#else
317 dtdz(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
318 & t(i,j,k ,nrhs,itrc))
319#endif
320 END DO
321 END DO
322 END IF
323!
324! Compute components of the rotated tracer flux (T m3/s) along
325! geopotential surfaces.
326!
327 IF (k.gt.0) THEN
328 DO j=jstr,jend
329 DO i=istr,iend+1
330#ifdef DIFF_3DCOEF
331 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
332 & on_u(i,j)
333#else
334 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
335 & on_u(i,j)
336#endif
337 fx(i,j)=cff* &
338 & (hz(i,j,k)+hz(i-1,j,k))* &
339 & (dtdx(i,j,k1)- &
340 & 0.5_r8*(min(dzdx(i,j,k1),0.0_r8)* &
341 & (dtdz(i-1,j,k1)+ &
342 & dtdz(i ,j,k2))+ &
343 & max(dzdx(i,j,k1),0.0_r8)* &
344 & (dtdz(i-1,j,k2)+ &
345 & dtdz(i ,j,k1))))
346 END DO
347 END DO
348 DO j=jstr,jend+1
349 DO i=istr,iend
350#ifdef DIFF_3DCOEF
351 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
352 & om_v(i,j)
353#else
354 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
355 & om_v(i,j)
356#endif
357 fe(i,j)=cff* &
358 & (hz(i,j,k)+hz(i,j-1,k))* &
359 & (dtde(i,j,k1)- &
360 & 0.5_r8*(min(dzde(i,j,k1),0.0_r8)* &
361 & (dtdz(i,j-1,k1)+ &
362 & dtdz(i,j ,k2))+ &
363 & max(dzde(i,j,k1),0.0_r8)* &
364 & (dtdz(i,j-1,k2)+ &
365 & dtdz(i,j ,k1))))
366 END DO
367 END DO
368 IF (k.lt.n(ng)) THEN
369 DO j=jstr,jend
370 DO i=istr,iend
371#ifdef DIFF_3DCOEF
372 cff=0.5_r8*diff3d_r(i,j,k)
373#else
374 cff=0.5_r8*diff2(i,j,itrc)
375#endif
376 cff1=min(dzdx(i ,j,k1),0.0_r8)
377 cff2=min(dzdx(i+1,j,k2),0.0_r8)
378 cff3=max(dzdx(i ,j,k2),0.0_r8)
379 cff4=max(dzdx(i+1,j,k1),0.0_r8)
380 fs(i,j,k2)=cff* &
381 & (cff1*(cff1*dtdz(i,j,k2)-dtdx(i ,j,k1))+ &
382 & cff2*(cff2*dtdz(i,j,k2)-dtdx(i+1,j,k2))+ &
383 & cff3*(cff3*dtdz(i,j,k2)-dtdx(i ,j,k2))+ &
384 & cff4*(cff4*dtdz(i,j,k2)-dtdx(i+1,j,k1)))
385 cff1=min(dzde(i,j ,k1),0.0_r8)
386 cff2=min(dzde(i,j+1,k2),0.0_r8)
387 cff3=max(dzde(i,j ,k2),0.0_r8)
388 cff4=max(dzde(i,j+1,k1),0.0_r8)
389 fs(i,j,k2)=fs(i,j,k2)+ &
390 & cff* &
391 & (cff1*(cff1*dtdz(i,j,k2)-dtde(i,j ,k1))+ &
392 & cff2*(cff2*dtdz(i,j,k2)-dtde(i,j+1,k2))+ &
393 & cff3*(cff3*dtdz(i,j,k2)-dtde(i,j ,k2))+ &
394 & cff4*(cff4*dtdz(i,j,k2)-dtde(i,j+1,k1)))
395 END DO
396 END DO
397 END IF
398!
399! Time-step harmonic, geopotential diffusion term (m Tunits).
400!
401 DO j=jstr,jend
402 DO i=istr,iend
403 cff=dt(ng)*pm(i,j)*pn(i,j)
404 cff1=cff*(fx(i+1,j )-fx(i,j))
405 cff2=cff*(fe(i ,j+1)-fe(i,j))
406 cff3=dt(ng)*(fs(i,j,k2)-fs(i,j,k1))
407 cff4=cff1+cff2+cff3
408 t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff4
409#ifdef DIAGNOSTICS_TS
410 diatwrk(i,j,k,itrc,itxdif)=cff1
411 diatwrk(i,j,k,itrc,itydif)=cff2
412 diatwrk(i,j,k,itrc,itsdif)=cff3
413 diatwrk(i,j,k,itrc,ithdif)=cff4
414#endif
415 END DO
416 END DO
417 END IF
418 END DO k_loop
419 END DO t_loop
420!
421 RETURN
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer, dimension(:), allocatable nt
Definition mod_param.F:489
real(dp), dimension(:), allocatable dt
integer itxdif
integer itsdif
integer itydif
logical, dimension(:,:), allocatable ltracerclm
integer ithdif

References mod_scalars::dt, mod_scalars::ithdif, mod_scalars::itsdif, mod_scalars::itxdif, mod_scalars::itydif, and mod_scalars::ltracerclm.

Referenced by t3dmix2().

Here is the caller graph for this function:

◆ t3dmix2_iso_tile()

subroutine t3dmix2_mod::t3dmix2_iso_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nrhs,
integer, intent(in) nstp,
integer, intent(in) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) hz,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) z_r,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) diff3d_r,
real(r8), dimension(lbi:ubi,lbj:ubj,nt(ng)), intent(in) diff2,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) pden,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),nt(ng)), intent(in) tclm,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),nt(ng), ndt), intent(inout) diatwrk,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),3,nt(ng)), intent(inout) t )
private

Definition at line 91 of file t3dmix2_iso.h.

116!***********************************************************************
117!
118 USE mod_param
119 USE mod_scalars
120!
121! Imported variable declarations.
122!
123 integer, intent(in) :: ng, tile
124 integer, intent(in) :: LBi, UBi, LBj, UBj
125 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
126 integer, intent(in) :: nrhs, nstp, nnew
127
128#ifdef ASSUMED_SHAPE
129# ifdef MASKING
130 real(r8), intent(in) :: umask(LBi:,LBj:)
131 real(r8), intent(in) :: vmask(LBi:,LBj:)
132# endif
133# ifdef WET_DRY
134 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
135 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
136# endif
137# ifdef DIFF_3DCOEF
138 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
139# else
140 real(r8), intent(in) :: diff2(LBi:,LBj:,:)
141# endif
142 real(r8), intent(in) :: om_v(LBi:,LBj:)
143 real(r8), intent(in) :: on_u(LBi:,LBj:)
144 real(r8), intent(in) :: pm(LBi:,LBj:)
145 real(r8), intent(in) :: pn(LBi:,LBj:)
146 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
147 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
148 real(r8), intent(in) :: pden(LBi:,LBj:,:)
149# ifdef TS_MIX_CLIMA
150 real(r8), intent(in) :: tclm(LBi:,LBj:,:,:)
151# endif
152# ifdef DIAGNOSTICS_TS
153 real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
154# endif
155 real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
156#else
157# ifdef MASKING
158 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
159 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
160# endif
161# ifdef WET_DRY
162 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
163 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
164# endif
165# ifdef DIFF_3DCOEF
166 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
167# else
168 real(r8), intent(in) :: diff2(LBi:UBi,LBj:UBj,NT(ng))
169# endif
170 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
171 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
172 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
173 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
174 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
175 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
176 real(r8), intent(in) :: pden(LBi:UBi,LBj:UBj,N(ng))
177# ifdef TS_MIX_CLIMA
178 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
179# endif
180# ifdef DIAGNOSTICS_TS
181 real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
182 & NDT)
183# endif
184 real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
185#endif
186!
187! Local variable declarations.
188!
189 integer :: i, itrc, j, k, k1, k2
190
191 real(r8), parameter :: eps = 0.5_r8
192 real(r8), parameter :: small = 1.0e-14_r8
193 real(r8), parameter :: slope_max = 0.0001_r8
194 real(r8), parameter :: strat_min = 0.1_r8
195
196 real(r8) :: cff, cff1, cff2, cff3, cff4
197
198 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE
199 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX
200
201 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: FS
202 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dRde
203 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dRdx
204 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTde
205 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdr
206 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdx
207
208#include "set_bounds.h"
209!
210!-----------------------------------------------------------------------
211! Compute horizontal harmonic diffusion along isopycnic surfaces.
212!-----------------------------------------------------------------------
213!
214! Compute horizontal and density gradients. Notice the recursive
215! blocking sequence. The vertical placement of the gradients is:
216!
217! dTdx,dTde(:,:,k1) k rho-points
218! dTdx,dTde(:,:,k2) k+1 rho-points
219! FS,dTdr(:,:,k1) k-1/2 W-points
220! FS,dTdr(:,:,k2) k+1/2 W-points
221!
222 t_loop : DO itrc=1,nt(ng)
223 k2=1
224 k_loop : DO k=0,n(ng)
225 k1=k2
226 k2=3-k1
227 IF (k.lt.n(ng)) THEN
228 DO j=jstr,jend
229 DO i=istr,iend+1
230 cff=0.5_r8*(pm(i,j)+pm(i-1,j))
231#ifdef MASKING
232 cff=cff*umask(i,j)
233#endif
234#ifdef WET_DRY
235 cff=cff*umask_wet(i,j)
236#endif
237 drdx(i,j,k2)=cff*(pden(i ,j,k+1)- &
238 & pden(i-1,j,k+1))
239#if defined TS_MIX_STABILITY
240 dtdx(i,j,k2)=cff*(0.75_r8*(t(i ,j,k+1,nrhs,itrc)- &
241 & t(i-1,j,k+1,nrhs,itrc))+ &
242 & 0.25_r8*(t(i ,j,k+1,nstp,itrc)- &
243 & t(i-1,j,k+1,nstp,itrc)))
244#elif defined TS_MIX_CLIMA
245 IF (ltracerclm(itrc,ng)) THEN
246 dtdx(i,j,k2)=cff*((t(i ,j,k+1,nrhs,itrc)- &
247 & tclm(i ,j,k+1,itrc))- &
248 & (t(i-1,j,k+1,nrhs,itrc)- &
249 & tclm(i-1,j,k+1,itrc)))
250 ELSE
251 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
252 & t(i-1,j,k+1,nrhs,itrc))
253 END IF
254#else
255 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
256 & t(i-1,j,k+1,nrhs,itrc))
257#endif
258 END DO
259 END DO
260 DO j=jstr,jend+1
261 DO i=istr,iend
262 cff=0.5_r8*(pn(i,j)+pn(i,j-1))
263#ifdef MASKING
264 cff=cff*vmask(i,j)
265#endif
266#ifdef WET_DRY
267 cff=cff*vmask_wet(i,j)
268#endif
269 drde(i,j,k2)=cff*(pden(i,j ,k+1)- &
270 & pden(i,j-1,k+1))
271#if defined TS_MIX_STABILITY
272 dtde(i,j,k2)=cff*(0.75_r8*(t(i,j ,k+1,nrhs,itrc)- &
273 & t(i,j-1,k+1,nrhs,itrc))+ &
274 & 0.25_r8*(t(i,j ,k+1,nstp,itrc)- &
275 & t(i,j-1,k+1,nstp,itrc)))
276#elif defined TS_MIX_CLIMA
277 IF (ltracerclm(itrc,ng)) THEN
278 dtde(i,j,k2)=cff*((t(i,j ,k+1,nrhs,itrc)- &
279 & tclm(i,j ,k+1,itrc))- &
280 & (t(i,j-1,k+1,nrhs,itrc)- &
281 & tclm(i,j-1,k+1,itrc)))
282 ELSE
283 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
284 & t(i,j-1,k+1,nrhs,itrc))
285 END IF
286#else
287 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
288 & t(i,j-1,k+1,nrhs,itrc))
289#endif
290 END DO
291 END DO
292 END IF
293 IF ((k.eq.0).or.(k.eq.n(ng))) THEN
294 DO j=jstr-1,jend+1
295 DO i=istr-1,iend+1
296 dtdr(i,j,k2)=0.0_r8
297 fs(i,j,k2)=0.0_r8
298 END DO
299 END DO
300 ELSE
301 DO j=jstr-1,jend+1
302 DO i=istr-1,iend+1
303#if defined TS_MIX_MAX_SLOPE
304 cff1=sqrt(drdx(i,j,k2)**2+drdx(i+1,j,k2)**2+ &
305 & drdx(i,j,k1)**2+drdx(i+1,j,k1)**2+ &
306 & drde(i,j,k2)**2+drde(i,j+1,k2)**2+ &
307 & drde(i,j,k1)**2+drde(i,j+1,k1)**2)
308 cff2=0.25_r8*slope_max* &
309 & (z_r(i,j,k+1)-z_r(i,j,k))*cff1
310 cff3=max(pden(i,j,k)-pden(i,j,k+1),small)
311 cff4=max(cff2,cff3)
312 cff=-1.0_r8/cff4
313#elif defined TS_MIX_MIN_STRAT
314 cff1=max(pden(i,j,k)-pden(i,j,k+1), &
315 & strat_min*(z_r(i,j,k+1)-z_r(i,j,k)))
316 cff=-1.0_r8/cff1
317#else
318 cff1=max(pden(i,j,k)-pden(i,j,k+1),eps)
319 cff=-1.0_r8/cff1
320#endif
321#if defined TS_MIX_STABILITY
322 dtdr(i,j,k2)=cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
323 & t(i,j,k ,nrhs,itrc))+ &
324 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
325 & t(i,j,k ,nstp,itrc)))
326#elif defined TS_MIX_CLIMA
327 IF (ltracerclm(itrc,ng)) THEN
328 dtdr(i,j,k2)=cff*((t(i,j,k+1,nrhs,itrc)- &
329 & tclm(i,j,k+1,itrc))- &
330 & (t(i,j,k ,nrhs,itrc)- &
331 & tclm(i,j,k ,itrc)))
332 ELSE
333 dtdr(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
334 & t(i,j,k ,nrhs,itrc))
335 END IF
336#else
337 dtdr(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
338 & t(i,j,k ,nrhs,itrc))
339#endif
340 fs(i,j,k2)=cff*(z_r(i,j,k+1)-z_r(i,j,k))
341 END DO
342 END DO
343 END IF
344!
345! Compute components of the rotated tracer flux (T m4/s) along
346! isopycnic surfaces.
347!
348 IF (k.gt.0) THEN
349 DO j=jstr,jend
350 DO i=istr,iend+1
351#ifdef DIFF_3DCOEF
352 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
353 & on_u(i,j)
354#else
355 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
356 & on_u(i,j)
357#endif
358 fx(i,j)=cff* &
359 & (hz(i,j,k)+hz(i-1,j,k))* &
360 & (dtdx(i,j,k1)- &
361 & 0.5_r8*(max(drdx(i,j,k1),0.0_r8)* &
362 & (dtdr(i-1,j,k1)+ &
363 & dtdr(i ,j,k2))+ &
364 & min(drdx(i,j,k1),0.0_r8)* &
365 & (dtdr(i-1,j,k2)+ &
366 & dtdr(i ,j,k1))))
367 END DO
368 END DO
369 DO j=jstr,jend+1
370 DO i=istr,iend
371#ifdef DIFF_3DCOEF
372 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
373 & om_v(i,j)
374#else
375 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
376 & om_v(i,j)
377#endif
378 fe(i,j)=cff* &
379 & (hz(i,j,k)+hz(i,j-1,k))* &
380 & (dtde(i,j,k1)- &
381 & 0.5_r8*(max(drde(i,j,k1),0.0_r8)* &
382 & (dtdr(i,j-1,k1)+ &
383 & dtdr(i,j ,k2))+ &
384 & min(drde(i,j,k1),0.0_r8)* &
385 & (dtdr(i,j-1,k2)+ &
386 & dtdr(i,j ,k1))))
387 END DO
388 END DO
389 IF (k.lt.n(ng)) THEN
390 DO j=jstr,jend
391 DO i=istr,iend
392 cff1=max(drdx(i ,j,k1),0.0_r8)
393 cff2=max(drdx(i+1,j,k2),0.0_r8)
394 cff3=min(drdx(i ,j,k2),0.0_r8)
395 cff4=min(drdx(i+1,j,k1),0.0_r8)
396 cff=cff1*(cff1*dtdr(i,j,k2)-dtdx(i ,j,k1))+ &
397 & cff2*(cff2*dtdr(i,j,k2)-dtdx(i+1,j,k2))+ &
398 & cff3*(cff3*dtdr(i,j,k2)-dtdx(i ,j,k2))+ &
399 & cff4*(cff4*dtdr(i,j,k2)-dtdx(i+1,j,k1))
400 cff1=max(drde(i,j ,k1),0.0_r8)
401 cff2=max(drde(i,j+1,k2),0.0_r8)
402 cff3=min(drde(i,j ,k2),0.0_r8)
403 cff4=min(drde(i,j+1,k1),0.0_r8)
404 cff=cff+ &
405 & cff1*(cff1*dtdr(i,j,k2)-dtde(i,j ,k1))+ &
406 & cff2*(cff2*dtdr(i,j,k2)-dtde(i,j+1,k2))+ &
407 & cff3*(cff3*dtdr(i,j,k2)-dtde(i,j ,k2))+ &
408 & cff4*(cff4*dtdr(i,j,k2)-dtde(i,j+1,k1))
409#ifdef DIFF_3DCOEF
410 fs(i,j,k2)=0.5_r8*cff*diff3d_r(i,j,k)*fs(i,j,k2)
411#else
412 fs(i,j,k2)=0.5_r8*cff*diff2(i,j,itrc)*fs(i,j,k2)
413#endif
414 END DO
415 END DO
416 END IF
417!
418! Time-step harmonic, isopycnic diffusion term (m Tunits).
419!
420 DO j=jstr,jend
421 DO i=istr,iend
422 cff=dt(ng)*pm(i,j)*pn(i,j)
423 cff1=cff*(fx(i+1,j )-fx(i,j))
424 cff2=cff*(fe(i ,j+1)-fe(i,j))
425 cff3=dt(ng)*(fs(i,j,k2)-fs(i,j,k1))
426 cff4=cff1+cff2+cff3
427 t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff4
428#ifdef DIAGNOSTICS_TS
429 diatwrk(i,j,k,itrc,itxdif)=cff1
430 diatwrk(i,j,k,itrc,itydif)=cff2
431 diatwrk(i,j,k,itrc,itsdif)=cff3
432 diatwrk(i,j,k,itrc,ithdif)=cff4
433#endif
434 END DO
435 END DO
436 END IF
437 END DO k_loop
438 END DO t_loop
439!
440 RETURN

References mod_scalars::dt, mod_scalars::ithdif, mod_scalars::itsdif, mod_scalars::itxdif, mod_scalars::itydif, and mod_scalars::ltracerclm.

◆ t3dmix2_s_tile()

subroutine t3dmix2_mod::t3dmix2_s_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nrhs,
integer, intent(in) nstp,
integer, intent(in) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask_wet,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) hz,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) diff3d_r,
real(r8), dimension(lbi:ubi,lbj:ubj,nt(ng)), intent(in) diff2,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),nt(ng)), intent(in) tclm,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),nt(ng), ndt), intent(inout) diatwrk,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),3,nt(ng)), intent(inout) t )
private

Definition at line 89 of file t3dmix2_s.h.

112!***********************************************************************
113!
114 USE mod_param
115 USE mod_scalars
116!
117! Imported variable declarations.
118!
119 integer, intent(in) :: ng, tile
120 integer, intent(in) :: LBi, UBi, LBj, UBj
121 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
122 integer, intent(in) :: nrhs, nstp, nnew
123
124#ifdef ASSUMED_SHAPE
125# ifdef MASKING
126 real(r8), intent(in) :: umask(LBi:,LBj:)
127 real(r8), intent(in) :: vmask(LBi:,LBj:)
128# endif
129# ifdef WET_DRY
130 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
131 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
132# endif
133# ifdef DIFF_3DCOEF
134 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
135# else
136 real(r8), intent(in) :: diff2(LBi:,LBj:,:)
137# endif
138 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
139 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
140 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
141 real(r8), intent(in) :: pm(LBi:,LBj:)
142 real(r8), intent(in) :: pn(LBi:,LBj:)
143# ifdef TS_MIX_CLIMA
144 real(r8), intent(in) :: tclm(LBi:,LBj:,:,:)
145# endif
146# ifdef DIAGNOSTICS_TS
147 real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
148# endif
149 real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
150#else
151# ifdef MASKING
152 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
153 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
154# endif
155# ifdef WET_DRY
156 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
157 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
158# endif
159# ifdef DIFF_3DCOEF
160 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
161# else
162 real(r8), intent(in) :: diff2(LBi:UBi,LBj:UBj,NT(ng))
163# endif
164 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
165 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
166 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
167 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
168 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
169# ifdef TS_MIX_CLIMA
170 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
171# endif
172# ifdef DIAGNOSTICS_TS
173 real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
174 & NDT)
175# endif
176 real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
177#endif
178!
179! Local variable declarations.
180!
181 integer :: i, itrc, j, k
182
183 real(r8) :: cff, cff1, cff2, cff3
184
185 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE
186 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX
187
188#include "set_bounds.h"
189!
190!-----------------------------------------------------------------------
191! Compute horizontal harmonic diffusion along constant S-surfaces.
192#ifdef TS_MIX_STABILITY
193! In order to increase stability, the harmonic operator is applied
194! as: 3/4 t(:,:,:,nrhs,:) + 1/4 t(:,:,:,nstp,:).
195#endif
196!-----------------------------------------------------------------------
197!
198 DO itrc=1,nt(ng)
199 DO k=1,n(ng)
200!
201! Compute XI- and ETA-components of diffusive tracer flux (T m3/s).
202!
203 DO j=jstr,jend
204 DO i=istr,iend+1
205#ifdef DIFF_3DCOEF
206 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
207 & pmon_u(i,j)
208#else
209 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
210 & pmon_u(i,j)
211#endif
212#if defined TS_MIX_STABILITY
213 fx(i,j)=cff* &
214 & (hz(i,j,k)+hz(i-1,j,k))* &
215 & (0.75_r8*(t(i ,j,k,nrhs,itrc)- &
216 & t(i-1,j,k,nrhs,itrc))+ &
217 & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
218 & t(i-1,j,k,nstp,itrc)))
219#elif defined TS_MIX_CLIMA
220 IF (ltracerclm(itrc,ng)) THEN
221 fx(i,j)=cff* &
222 & (hz(i,j,k)+hz(i-1,j,k))* &
223 & ((t(i ,j,k,nrhs,itrc)-tclm(i ,j,k,itrc))- &
224 & (t(i-1,j,k,nrhs,itrc)-tclm(i-1,j,k,itrc)))
225 ELSE
226 fx(i,j)=cff* &
227 & (hz(i,j,k)+hz(i-1,j,k))* &
228 & (t(i,j,k,nrhs,itrc)-t(i-1,j,k,nrhs,itrc))
229 END IF
230#else
231 fx(i,j)=cff* &
232 & (hz(i,j,k)+hz(i-1,j,k))* &
233 & (t(i,j,k,nrhs,itrc)-t(i-1,j,k,nrhs,itrc))
234#endif
235#ifdef MASKING
236 fx(i,j)=fx(i,j)*umask(i,j)
237#endif
238#ifdef WET_DRY
239 fx(i,j)=fx(i,j)*umask_wet(i,j)
240#endif
241 END DO
242 END DO
243 DO j=jstr,jend+1
244 DO i=istr,iend
245#ifdef DIFF_3DCOEF
246 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
247 & pnom_v(i,j)
248#else
249 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
250 & pnom_v(i,j)
251#endif
252#if defined TS_MIX_STABILITY
253 fe(i,j)=cff* &
254 & (hz(i,j,k)+hz(i,j-1,k))* &
255 & (0.75_r8*(t(i,j ,k,nrhs,itrc)- &
256 & t(i,j-1,k,nrhs,itrc))+ &
257 & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
258 & t(i,j-1,k,nstp,itrc)))
259#elif defined TS_MIX_CLIMA
260 IF (ltracerclm(itrc,ng)) THEN
261 fe(i,j)=cff* &
262 & (hz(i,j,k)+hz(i,j-1,k))* &
263 & ((t(i,j ,k,nrhs,itrc)-tclm(i,j ,k,itrc))- &
264 & (t(i,j-1,k,nrhs,itrc)-tclm(i,j-1,k,itrc)))
265 ELSE
266 fe(i,j)=cff* &
267 & (hz(i,j,k)+hz(i,j-1,k))* &
268 & (t(i,j,k,nrhs,itrc)-t(i,j-1,k,nrhs,itrc))
269 END IF
270#else
271 fe(i,j)=cff* &
272 & (hz(i,j,k)+hz(i,j-1,k))* &
273 & (t(i,j,k,nrhs,itrc)-t(i,j-1,k,nrhs,itrc))
274#endif
275#ifdef MASKING
276 fe(i,j)=fe(i,j)*vmask(i,j)
277#endif
278#ifdef WET_DRY
279 fe(i,j)=fe(i,j)*vmask_wet(i,j)
280#endif
281 END DO
282 END DO
283!
284! Time-step harmonic, S-surfaces diffusion term (m Tunits).
285!
286 DO j=jstr,jend
287 DO i=istr,iend
288 cff=dt(ng)*pm(i,j)*pn(i,j)
289 cff1=cff*(fx(i+1,j )-fx(i,j))
290 cff2=cff*(fe(i ,j+1)-fe(i,j))
291 cff3=cff1+cff2
292 t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff3
293#ifdef DIAGNOSTICS_TS
294 diatwrk(i,j,k,itrc,itxdif)=cff1
295 diatwrk(i,j,k,itrc,itydif)=cff2
296 diatwrk(i,j,k,itrc,ithdif)=cff3
297#endif
298 END DO
299 END DO
300 END DO
301 END DO
302!
303 RETURN

References mod_scalars::dt, mod_scalars::ithdif, mod_scalars::itxdif, mod_scalars::itydif, and mod_scalars::ltracerclm.