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

Functions/Subroutines

subroutine, public tl_t3dmix2 (ng, tile)
 
subroutine tl_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, tl_hz, z_r, tl_z_r, diff3d_r, diff2, tclm, t, tl_t)
 
subroutine tl_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, tl_hz, z_r, tl_z_r, diff3d_r, diff2, pden, tl_pden, tclm, t, tl_t)
 
subroutine tl_t3dmix2_s_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nstp, nnew, umask, vmask, umask_wet, vmask_wet, hz, tl_hz, pmon_u, pnom_v, pm, pn, diff3d_r, diff2, tclm, t, tl_t)
 

Function/Subroutine Documentation

◆ tl_t3dmix2()

subroutine public tl_t3dmix2_mod::tl_t3dmix2 ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 25 of file tl_t3dmix2_geo.h.

26!***********************************************************************
27!
28 USE mod_param
29#ifdef TS_MIX_CLIMA
30 USE mod_clima
31#endif
32#ifdef DIAGNOSTICS_TS
33!! USE mod_diags
34#endif
35 USE mod_grid
36 USE mod_mixing
37 USE mod_ocean
38 USE mod_stepping
39!
40! Imported variable declarations.
41!
42 integer, intent(in) :: ng, tile
43!
44! Local variable declarations.
45!
46 character (len=*), parameter :: MyFile = &
47 & __FILE__
48!
49#include "tile.h"
50!
51#ifdef PROFILE
52 CALL wclock_on (ng, itlm, 25, __line__, myfile)
53#endif
54 CALL tl_t3dmix2_geo_tile (ng, tile, &
55 & lbi, ubi, lbj, ubj, &
56 & imins, imaxs, jmins, jmaxs, &
57 & nrhs(ng), nstp(ng), nnew(ng), &
58#ifdef MASKING
59 & grid(ng) % umask, &
60 & grid(ng) % vmask, &
61#endif
62#ifdef WET_DRY_NOT_YET
63 & grid(ng) % umask_wet, &
64 & grid(ng) % vmask_wet, &
65#endif
66 & grid(ng) % om_v, &
67 & grid(ng) % on_u, &
68 & grid(ng) % pm, &
69 & grid(ng) % pn, &
70 & grid(ng) % Hz, &
71 & grid(ng) % tl_Hz, &
72 & grid(ng) % z_r, &
73 & grid(ng) % tl_z_r, &
74#ifdef DIFF_3DCOEF
75 & mixing(ng) % diff3d_r, &
76#else
77 & mixing(ng) % diff2, &
78#endif
79#ifdef TS_MIX_CLIMA
80 & clima(ng) % tclm, &
81#endif
82#ifdef DIAGNOSTICS_TS
83!! & DIAGS(ng) % DiaTwrk, &
84#endif
85 & ocean(ng) % t, &
86 & ocean(ng) % tl_t)
87#ifdef PROFILE
88 CALL wclock_off (ng, itlm, 25, __line__, myfile)
89#endif
90!
91 RETURN
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
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 itlm
Definition mod_param.F:663
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_grid::grid, mod_param::itlm, mod_mixing::mixing, mod_stepping::nnew, mod_stepping::nrhs, mod_stepping::nstp, mod_ocean::ocean, tl_t3dmix2_geo_tile(), wclock_off(), and wclock_on().

Referenced by tl_rhs3d_mod::tl_rhs3d().

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

◆ tl_t3dmix2_geo_tile()

subroutine tl_t3dmix2_mod::tl_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) tl_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) tl_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),3,nt(ng)), intent(in) t,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),3,nt(ng)), intent(inout) tl_t )
private

Definition at line 95 of file tl_t3dmix2_geo.h.

120!***********************************************************************
121!
122 USE mod_param
123 USE mod_scalars
124!
125! Imported variable declarations.
126!
127 integer, intent(in) :: ng, tile
128 integer, intent(in) :: LBi, UBi, LBj, UBj
129 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
130 integer, intent(in) :: nrhs, nstp, nnew
131
132#ifdef ASSUMED_SHAPE
133# ifdef MASKING
134 real(r8), intent(in) :: umask(LBi:,LBj:)
135 real(r8), intent(in) :: vmask(LBi:,LBj:)
136# endif
137# ifdef WET_DRY_NOT_YET
138 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
139 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
140# endif
141# ifdef DIFF_3DCOEF
142 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
143# else
144 real(r8), intent(in) :: diff2(LBi:,LBj:,:)
145# endif
146 real(r8), intent(in) :: om_v(LBi:,LBj:)
147 real(r8), intent(in) :: on_u(LBi:,LBj:)
148 real(r8), intent(in) :: pm(LBi:,LBj:)
149 real(r8), intent(in) :: pn(LBi:,LBj:)
150 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
151 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
152 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
153# ifdef TS_MIX_CLIMA
154 real(r8), intent(in) :: tclm(LBi:,LBj:,:,:)
155# endif
156 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
157 real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:)
158# ifdef DIAGNOSTICS_TS
159!! real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
160# endif
161
162 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
163#else
164# ifdef MASKING
165 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
166 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
167# endif
168# ifdef WET_DRY
169 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
170 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
171# endif
172# ifdef DIFF_3DCOEF
173 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
174# else
175 real(r8), intent(in) :: diff2(LBi:UBi,LBj:UBj,NT(ng))
176# endif
177 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
178 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
179 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
180 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
181 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
182 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
183 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
184# ifdef TS_MIX_CLIMA
185 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
186# endif
187 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
188 real(r8), intent(in) :: tl_z_r(LBi:UBi,LBj:UBj,N(ng))
189# ifdef DIAGNOSTICS_TS
190!! real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
191!! & NDT)
192# endif
193 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
194#endif
195!
196! Local variable declarations.
197!
198 integer :: i, itrc, j, k, k1, k2
199
200 real(r8) :: cff, cff1, cff2, cff3, cff4
201 real(r8) :: tl_cff, tl_cff1, tl_cff2, tl_cff3, tl_cff4
202
203 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FE
204 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FX
205
206 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdz
207 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdx
208 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTde
209 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZdx
210 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZde
211
212 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_FS
213 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dTdz
214 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dTdx
215 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dTde
216 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dZdx
217 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dZde
218
219#include "set_bounds.h"
220!
221!-----------------------------------------------------------------------
222! Compute horizontal harmonic diffusion along geopotential surfaces.
223!-----------------------------------------------------------------------
224!
225! Compute horizontal and vertical gradients. Notice the recursive
226! blocking sequence. The vertical placement of the gradients is:
227!
228! dTdx,dTde(:,:,k1) k rho-points
229! dTdx,dTde(:,:,k2) k+1 rho-points
230! FS,dTdz(:,:,k1) k-1/2 W-points
231! FS,dTdz(:,:,k2) k+1/2 W-points
232!
233#ifdef TS_MIX_STABILITY
234! In order to increase stability, the biharmonic operator is applied
235! as: 3/4 t(:,:,:,nrhs,:) + 1/4 t(:,:,:,nstp,:).
236!
237#endif
238
239 t_loop : DO itrc=1,nt(ng)
240 k2=1
241 k_loop : DO k=0,n(ng)
242 k1=k2
243 k2=3-k1
244 IF (k.lt.n(ng)) THEN
245 DO j=jstr,jend
246 DO i=istr,iend+1
247 cff=0.5_r8*(pm(i,j)+pm(i-1,j))
248#ifdef MASKING
249 cff=cff*umask(i,j)
250#endif
251#ifdef WET_DRY_NOT_YET
252 cff=cff*umask_wet(i,j)
253#endif
254 dzdx(i,j,k2)=cff*(z_r(i ,j,k+1)- &
255 & z_r(i-1,j,k+1))
256 tl_dzdx(i,j,k2)=cff*(tl_z_r(i ,j,k+1)- &
257 & tl_z_r(i-1,j,k+1))
258#if defined TS_MIX_STABILITY
259 dtdx(i,j,k2)=cff*(0.75_r8*(t(i ,j,k+1,nrhs,itrc)- &
260 & t(i-1,j,k+1,nrhs,itrc))+ &
261 & 0.25_r8*(t(i ,j,k+1,nstp,itrc)- &
262 & t(i-1,j,k+1,nstp,itrc)))
263 tl_dtdx(i,j,k2)=cff* &
264 & (0.75_r8*(tl_t(i ,j,k+1,nrhs,itrc)- &
265 & tl_t(i-1,j,k+1,nrhs,itrc))+ &
266 & 0.25_r8*(tl_t(i ,j,k+1,nstp,itrc)- &
267 & tl_t(i-1,j,k+1,nstp,itrc)))
268#elif defined TS_MIX_CLIMA
269 IF (ltracerclm(itrc,ng)) THEN
270 dtdx(i,j,k2)=cff*((t(i ,j,k+1,nrhs,itrc)- &
271 & tclm(i ,j,k+1,itrc))- &
272 & (t(i-1,j,k+1,nrhs,itrc)- &
273 & tclm(i-1,j,k+1,itrc)))
274 ELSE
275 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
276 & t(i-1,j,k+1,nrhs,itrc))
277 END IF
278 tl_dtdx(i,j,k2)=cff*(tl_t(i ,j,k+1,nrhs,itrc)- &
279 & tl_t(i-1,j,k+1,nrhs,itrc))
280#else
281 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
282 & t(i-1,j,k+1,nrhs,itrc))
283 tl_dtdx(i,j,k2)=cff*(tl_t(i ,j,k+1,nrhs,itrc)- &
284 & tl_t(i-1,j,k+1,nrhs,itrc))
285#endif
286 END DO
287 END DO
288 DO j=jstr,jend+1
289 DO i=istr,iend
290 cff=0.5_r8*(pn(i,j)+pn(i,j-1))
291#ifdef MASKING
292 cff=cff*vmask(i,j)
293#endif
294#ifdef WET_DRY_NOT_YET
295 cff=cff*vmask_wet(i,j)
296#endif
297 dzde(i,j,k2)=cff*(z_r(i,j ,k+1)- &
298 & z_r(i,j-1,k+1))
299 tl_dzde(i,j,k2)=cff*(tl_z_r(i,j ,k+1)- &
300 & tl_z_r(i,j-1,k+1))
301#if defined TS_MIX_STABILITY
302 dtde(i,j,k2)=cff*(0.75_r8*(t(i,j ,k+1,nrhs,itrc)- &
303 & t(i,j-1,k+1,nrhs,itrc))+ &
304 & 0.25_r8*(t(i,j ,k+1,nstp,itrc)- &
305 & t(i,j-1,k+1,nstp,itrc)))
306 tl_dtde(i,j,k2)=cff* &
307 & (0.75_r8*(tl_t(i,j ,k+1,nrhs,itrc)- &
308 & tl_t(i,j-1,k+1,nrhs,itrc))+ &
309 & 0.25_r8*(tl_t(i,j ,k+1,nstp,itrc)- &
310 & tl_t(i,j-1,k+1,nstp,itrc)))
311#elif defined TS_MIX_CLIMA
312 IF (ltracerclm(itrc,ng)) THEN
313 dtde(i,j,k2)=cff*((t(i,j ,k+1,nrhs,itrc)- &
314 & tclm(i,j ,k+1,itrc))- &
315 & (t(i,j-1,k+1,nrhs,itrc)- &
316 & tclm(i,j-1,k+1,itrc)))
317 ELSE
318 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
319 & t(i,j-1,k+1,nrhs,itrc))
320 END IF
321 tl_dtde(i,j,k2)=cff*(tl_t(i,j ,k+1,nrhs,itrc)- &
322 & tl_t(i,j-1,k+1,nrhs,itrc))
323#else
324 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
325 & t(i,j-1,k+1,nrhs,itrc))
326 tl_dtde(i,j,k2)=cff*(tl_t(i,j ,k+1,nrhs,itrc)- &
327 & tl_t(i,j-1,k+1,nrhs,itrc))
328#endif
329 END DO
330 END DO
331 END IF
332 IF ((k.eq.0).or.(k.eq.n(ng))) THEN
333 DO j=jstr-1,jend+1
334 DO i=istr-1,iend+1
335 dtdz(i,j,k2)=0.0_r8
336 tl_dtdz(i,j,k2)=0.0_r8
337!^ FS(i,j,k2)=0.0_r8
338!^
339 tl_fs(i,j,k2)=0.0_r8
340 END DO
341 END DO
342 ELSE
343 DO j=jstr-1,jend+1
344 DO i=istr-1,iend+1
345 cff=1.0_r8/(z_r(i,j,k+1)-z_r(i,j,k))
346 tl_cff=-cff*cff*(tl_z_r(i,j,k+1)-tl_z_r(i,j,k))
347#if defined TS_MIX_STABILITY
348 dtdz(i,j,k2)=cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
349 & t(i,j,k ,nrhs,itrc))+ &
350 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
351 & t(i,j,k ,nstp,itrc)))
352 tl_dtdz(i,j,k2)=tl_cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
353 & t(i,j,k ,nrhs,itrc))+ &
354 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
355 & t(i,j,k ,nstp,itrc)))+&
356 & cff*(0.75_r8*(tl_t(i,j,k+1,nrhs,itrc)- &
357 & tl_t(i,j,k ,nrhs,itrc))+ &
358 & 0.25_r8*(tl_t(i,j,k+1,nstp,itrc)- &
359 & tl_t(i,j,k ,nstp,itrc)))
360#elif defined TS_MIX_CLIMA
361 IF (ltracerclm(itrc,ng)) THEN
362 dtdz(i,j,k2)=cff*((t(i,j,k+1,nrhs,itrc)- &
363 & tclm(i,j,k+1,itrc))- &
364 & (t(i,j,k ,nrhs,itrc)- &
365 & tclm(i,j,k ,itrc)))
366 tl_dtdz(i,j,k2)=tl_cff*((t(i,j,k+1,nrhs,itrc)- &
367 & tclm(i,j,k+1,itrc))- &
368 & (t(i,j,k ,nrhs,itrc)- &
369 & tclm(i,j,k ,itrc)))+ &
370 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
371 & tl_t(i,j,k ,nrhs,itrc))
372
373 ELSE
374 dtdz(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
375 & t(i,j,k ,nrhs,itrc))
376 tl_dtdz(i,j,k2)=tl_cff*(t(i,j,k+1,nrhs,itrc)- &
377 & t(i,j,k ,nrhs,itrc))+ &
378 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
379 & tl_t(i,j,k ,nrhs,itrc))
380 END IF
381#else
382 dtdz(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
383 & t(i,j,k ,nrhs,itrc))
384 tl_dtdz(i,j,k2)=tl_cff*(t(i,j,k+1,nrhs,itrc)- &
385 & t(i,j,k ,nrhs,itrc))+ &
386 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
387 & tl_t(i,j,k ,nrhs,itrc))
388#endif
389 END DO
390 END DO
391 END IF
392!
393! Compute components of the rotated tracer flux (T m3/s) along
394! geopotential surfaces.
395!
396 IF (k.gt.0) THEN
397 DO j=jstr,jend
398 DO i=istr,iend+1
399#ifdef DIFF_3DCOEF
400 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
401 & on_u(i,j)
402#else
403 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
404 & on_u(i,j)
405#endif
406!^ FX(i,j)=cff* &
407!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
408!^ & (dTdx(i,j,k1)- &
409!^ & 0.5_r8*(MIN(dZdx(i,j,k1),0.0_r8)* &
410!^ & (dTdz(i-1,j,k1)+ &
411!^ & dTdz(i ,j,k2))+ &
412!^ & MAX(dZdx(i,j,k1),0.0_r8)* &
413!^ & (dTdz(i-1,j,k2)+ &
414!^ & dTdz(i ,j,k1))))
415!^
416 tl_fx(i,j)=cff* &
417 & (((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
418 & (dtdx(i,j,k1)- &
419 & 0.5_r8*(min(dzdx(i,j,k1),0.0_r8)* &
420 & (dtdz(i-1,j,k1)+ &
421 & dtdz(i ,j,k2))+ &
422 & max(dzdx(i,j,k1),0.0_r8)* &
423 & (dtdz(i-1,j,k2)+ &
424 & dtdz(i ,j,k1)))))+ &
425 & ((hz(i,j,k)+hz(i-1,j,k))* &
426 & (tl_dtdx(i,j,k1)- &
427 & 0.5_r8*(min(dzdx(i,j,k1),0.0_r8)* &
428 & (tl_dtdz(i-1,j,k1)+ &
429 & tl_dtdz(i ,j,k2))+ &
430 & max(dzdx(i,j,k1),0.0_r8)* &
431 & (tl_dtdz(i-1,j,k2)+ &
432 & tl_dtdz(i ,j,k1)))- &
433 & 0.5_r8*((0.5_r8+ &
434 & sign(0.5_r8,-dzdx(i,j,k1)))* &
435 & tl_dzdx(i,j,k1)* &
436 & (dtdz(i-1,j,k1)+dtdz(i,j,k2))+ &
437 & (0.5_r8+ &
438 & sign(0.5_r8, dzdx(i,j,k1)))* &
439 & tl_dzdx(i,j,k1)* &
440 & (dtdz(i-1,j,k2)+dtdz(i,j,k1))))))
441 END DO
442 END DO
443 DO j=jstr,jend+1
444 DO i=istr,iend
445#ifdef DIFF_3DCOEF
446 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
447 & om_v(i,j)
448#else
449 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
450 & om_v(i,j)
451#endif
452!^ FE(i,j)=cff* &
453!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
454!^ & (dTde(i,j,k1)- &
455!^ & 0.5_r8*(MIN(dZde(i,j,k1),0.0_r8)* &
456!^ & (dTdz(i,j-1,k1)+ &
457!^ & dTdz(i,j ,k2))+ &
458!^ & MAX(dZde(i,j,k1),0.0_r8)* &
459!^ & (dTdz(i,j-1,k2)+ &
460!^ & dTdz(i,j ,k1))))
461!^
462 tl_fe(i,j)=cff* &
463 & (((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
464 & (dtde(i,j,k1)- &
465 & 0.5_r8*(min(dzde(i,j,k1),0.0_r8)* &
466 & (dtdz(i,j-1,k1)+ &
467 & dtdz(i,j ,k2))+ &
468 & max(dzde(i,j,k1),0.0_r8)* &
469 & (dtdz(i,j-1,k2)+ &
470 & dtdz(i,j ,k1)))))+ &
471 & ((hz(i,j,k)+hz(i,j-1,k))* &
472 & (tl_dtde(i,j,k1)- &
473 & 0.5_r8*(min(dzde(i,j,k1),0.0_r8)* &
474 & (tl_dtdz(i,j-1,k1)+ &
475 & tl_dtdz(i,j ,k2))+ &
476 & max(dzde(i,j,k1),0.0_r8)* &
477 & (tl_dtdz(i,j-1,k2)+ &
478 & tl_dtdz(i,j ,k1)))- &
479 & 0.5_r8*((0.5_r8+ &
480 & sign(0.5_r8,-dzde(i,j,k1)))* &
481 & tl_dzde(i,j,k1)* &
482 & (dtdz(i,j-1,k1)+dtdz(i,j,k2))+ &
483 & (0.5_r8+ &
484 & sign(0.5_r8, dzde(i,j,k1)))* &
485 & tl_dzde(i,j,k1)* &
486 & (dtdz(i,j-1,k2)+dtdz(i,j,k1))))))
487 END DO
488 END DO
489 IF (k.lt.n(ng)) THEN
490 DO j=jstr,jend
491 DO i=istr,iend
492#ifdef DIFF_3DCOEF
493 cff=0.5_r8*diff3d_r(i,j,k)
494#else
495 cff=0.5_r8*diff2(i,j,itrc)
496#endif
497 cff1=min(dzdx(i ,j,k1),0.0_r8)
498 cff2=min(dzdx(i+1,j,k2),0.0_r8)
499 cff3=max(dzdx(i ,j,k2),0.0_r8)
500 cff4=max(dzdx(i+1,j,k1),0.0_r8)
501 tl_cff1=(0.5_r8+sign(0.5_r8,-dzdx(i ,j,k1)))* &
502 & tl_dzdx(i ,j,k1)
503 tl_cff2=(0.5_r8+sign(0.5_r8,-dzdx(i+1,j,k2)))* &
504 & tl_dzdx(i+1,j,k2)
505 tl_cff3=(0.5_r8+sign(0.5_r8, dzdx(i ,j,k2)))* &
506 & tl_dzdx(i ,j,k2)
507 tl_cff4=(0.5_r8+sign(0.5_r8, dzdx(i+1,j,k1)))* &
508 & tl_dzdx(i+1,j,k1)
509!^ FS(i,j,k2)=cff* &
510!^ & (cff1*(cff1*dTdz(i,j,k2)-dTdx(i ,j,k1))+ &
511!^ & cff2*(cff2*dTdz(i,j,k2)-dTdx(i+1,j,k2))+ &
512!^ & cff3*(cff3*dTdz(i,j,k2)-dTdx(i ,j,k2))+ &
513!^ & cff4*(cff4*dTdz(i,j,k2)-dTdx(i+1,j,k1)))
514!^
515 tl_fs(i,j,k2)=cff* &
516 & (tl_cff1*(cff1*dtdz(i,j,k2)- &
517 & dtdx(i ,j,k1))+ &
518 & tl_cff2*(cff2*dtdz(i,j,k2)- &
519 & dtdx(i+1,j,k2))+ &
520 & tl_cff3*(cff3*dtdz(i,j,k2)- &
521 & dtdx(i ,j,k2))+ &
522 & tl_cff4*(cff4*dtdz(i,j,k2)- &
523 & dtdx(i+1,j,k1))+ &
524 & cff1*(tl_cff1*dtdz(i,j,k2)+ &
525 & cff1*tl_dtdz(i,j,k2)- &
526 & tl_dtdx(i ,j,k1))+ &
527 & cff2*(tl_cff2*dtdz(i,j,k2)+ &
528 & cff2*tl_dtdz(i,j,k2)- &
529 & tl_dtdx(i+1,j,k2))+ &
530 & cff3*(tl_cff3*dtdz(i,j,k2)+ &
531 & cff3*tl_dtdz(i,j,k2)- &
532 & tl_dtdx(i ,j,k2))+ &
533 & cff4*(tl_cff4*dtdz(i,j,k2)+ &
534 & cff4*tl_dtdz(i,j,k2)- &
535 & tl_dtdx(i+1,j,k1)))
536!
537 cff1=min(dzde(i,j ,k1),0.0_r8)
538 cff2=min(dzde(i,j+1,k2),0.0_r8)
539 cff3=max(dzde(i,j ,k2),0.0_r8)
540 cff4=max(dzde(i,j+1,k1),0.0_r8)
541 tl_cff1=(0.5_r8+sign(0.5_r8,-dzde(i,j ,k1)))* &
542 & tl_dzde(i,j ,k1)
543 tl_cff2=(0.5_r8+sign(0.5_r8,-dzde(i,j+1,k2)))* &
544 & tl_dzde(i,j+1,k2)
545 tl_cff3=(0.5_r8+sign(0.5_r8, dzde(i,j ,k2)))* &
546 & tl_dzde(i,j ,k2)
547 tl_cff4=(0.5_r8+sign(0.5_r8, dzde(i,j+1,k1)))* &
548 & tl_dzde(i,j+1,k1)
549!^ FS(i,j,k2)=FS(i,j,k2)+ &
550!^ & cff* &
551!^ & (cff1*(cff1*dTdz(i,j,k2)-dTde(i,j ,k1))+ &
552!^ & cff2*(cff2*dTdz(i,j,k2)-dTde(i,j+1,k2))+ &
553!^ & cff3*(cff3*dTdz(i,j,k2)-dTde(i,j ,k2))+ &
554!^ & cff4*(cff4*dTdz(i,j,k2)-dTde(i,j+1,k1)))
555!^
556 tl_fs(i,j,k2)=tl_fs(i,j,k2)+ &
557 & cff* &
558 & (tl_cff1*(cff1*dtdz(i,j,k2)- &
559 & dtde(i,j ,k1))+ &
560 & tl_cff2*(cff2*dtdz(i,j,k2)- &
561 & dtde(i,j+1,k2))+ &
562 & tl_cff3*(cff3*dtdz(i,j,k2)- &
563 & dtde(i,j ,k2))+ &
564 & tl_cff4*(cff4*dtdz(i,j,k2)- &
565 & dtde(i,j+1,k1))+ &
566 & cff1*(tl_cff1*dtdz(i,j,k2)+ &
567 & cff1*tl_dtdz(i,j,k2)- &
568 & tl_dtde(i,j ,k1))+ &
569 & cff2*(tl_cff2*dtdz(i,j,k2)+ &
570 & cff2*tl_dtdz(i,j,k2)- &
571 & tl_dtde(i,j+1,k2))+ &
572 & cff3*(tl_cff3*dtdz(i,j,k2)+ &
573 & cff3*tl_dtdz(i,j,k2)- &
574 & tl_dtde(i,j ,k2))+ &
575 & cff4*(tl_cff4*dtdz(i,j,k2)+ &
576 & cff4*tl_dtdz(i,j,k2)- &
577 & tl_dtde(i,j+1,k1)))
578 END DO
579 END DO
580 END IF
581!
582! Time-step harmonic, geopotential diffusion term (m Tunits).
583!
584 DO j=jstr,jend
585 DO i=istr,iend
586!^ cff=dt(ng)*pm(i,j)*pn(i,j)* &
587!^ & (FX(i+1,j)-FX(i,j)+ &
588!^ & FE(i,j+1)-FE(i,j))+ &
589!^ & dt(ng)*(FS(i,j,k2)-FS(i,j,k1))
590!^
591 tl_cff=dt(ng)*pm(i,j)*pn(i,j)* &
592 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
593 & tl_fe(i,j+1)-tl_fe(i,j))+ &
594 & dt(ng)*(tl_fs(i,j,k2)-tl_fs(i,j,k1))
595!^ t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff
596!^
597 tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+tl_cff
598#ifdef DIAGNOSTICS_TS
599!! DiaTwrk(i,j,k,itrc,iThdif)=cff
600#endif
601 END DO
602 END DO
603 END IF
604 END DO k_loop
605 END DO t_loop
606!
607 RETURN
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer, dimension(:), allocatable nt
Definition mod_param.F:489
real(dp), dimension(:), allocatable dt
logical, dimension(:,:), allocatable ltracerclm

References mod_scalars::dt, and mod_scalars::ltracerclm.

Referenced by tl_t3dmix2().

Here is the caller graph for this function:

◆ tl_t3dmix2_iso_tile()

subroutine tl_t3dmix2_mod::tl_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) tl_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) tl_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)), intent(in) tl_pden,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),nt(ng)), intent(in) tclm,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),3,nt(ng)), intent(in) t,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),3,nt(ng)), intent(inout) tl_t )
private

Definition at line 97 of file tl_t3dmix2_iso.h.

123!***********************************************************************
124!
125 USE mod_param
126 USE mod_scalars
127!
128! Imported variable declarations.
129!
130 integer, intent(in) :: ng, tile
131 integer, intent(in) :: LBi, UBi, LBj, UBj
132 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
133 integer, intent(in) :: nrhs, nstp, nnew
134
135#ifdef ASSUMED_SHAPE
136# ifdef MASKING
137 real(r8), intent(in) :: umask(LBi:,LBj:)
138 real(r8), intent(in) :: vmask(LBi:,LBj:)
139# endif
140# ifdef WET_DRY_NOT_YET
141 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
142 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
143# endif
144# ifdef DIFF_3DCOEF
145 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
146# else
147 real(r8), intent(in) :: diff2(LBi:,LBj:,:)
148# endif
149 real(r8), intent(in) :: om_v(LBi:,LBj:)
150 real(r8), intent(in) :: on_u(LBi:,LBj:)
151 real(r8), intent(in) :: pm(LBi:,LBj:)
152 real(r8), intent(in) :: pn(LBi:,LBj:)
153 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
154 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
155 real(r8), intent(in) :: pden(LBi:,LBj:,:)
156 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
157# ifdef TS_MIX_CLIMA
158 real(r8), intent(in) :: tclm(LBi:,LBj:,:,:)
159# endif
160 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
161 real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:)
162 real(r8), intent(in) :: tl_pden(LBi:,LBj:,:)
163# ifdef DIAGNOSTICS_TS
164!! real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
165# endif
166
167 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
168#else
169# ifdef MASKING
170 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
171 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
172# endif
173# ifdef WET_DRY_NOT_YET
174 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
175 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
176# endif
177# ifdef DIFF_3DCOEF
178 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
179# else
180 real(r8), intent(in) :: diff2(LBi:UBi,LBj:UBj,NT(ng))
181# endif
182 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
183 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
184 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
185 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
186 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
187 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
188 real(r8), intent(in) :: pden(LBi:UBi,LBj:UBj,N(ng))
189 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
190# ifdef TS_MIX_CLIMA
191 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
192# endif
193 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
194 real(r8), intent(in) :: tl_z_r(LBi:UBi,LBj:UBj,N(ng))
195 real(r8), intent(in) :: tl_pden(LBi:UBi,LBj:UBj,N(ng))
196# ifdef DIAGNOSTICS_TS
197!! real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
198!! & NDT)
199# endif
200 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
201#endif
202!
203! Local variable declarations.
204!
205 integer :: i, itrc, j, k, k1, k2
206
207 real(r8), parameter :: eps = 0.5_r8
208 real(r8), parameter :: small = 1.0e-14_r8
209 real(r8), parameter :: slope_max = 0.0001_r8
210 real(r8), parameter :: strat_min = 0.1_r8
211
212 real(r8) :: cff, cff1, cff2, cff3, cff4
213 real(r8) :: tl_cff, tl_cff1, tl_cff2, tl_cff3, tl_cff4
214
215 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FE
216 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FX
217
218 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: FS
219 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dRde
220 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dRdx
221 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTde
222 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdr
223 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dTdx
224
225 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_FS
226 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dRde
227 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dRdx
228 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dTde
229 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dTdr
230 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: tl_dTdx
231
232#include "set_bounds.h"
233!
234!-----------------------------------------------------------------------
235! Compute horizontal harmonic diffusion along isopycnic surfaces.
236!-----------------------------------------------------------------------
237!
238! Compute horizontal and density gradients. Notice the recursive
239! blocking sequence. The vertical placement of the gradients is:
240!
241! dTdx,dTde(:,:,k1) k rho-points
242! dTdx,dTde(:,:,k2) k+1 rho-points
243! FS,dTdr(:,:,k1) k-1/2 W-points
244! FS,dTdr(:,:,k2) k+1/2 W-points
245!
246 t_loop : DO itrc=1,nt(ng)
247 k2=1
248 k_loop : DO k=0,n(ng)
249 k1=k2
250 k2=3-k1
251 IF (k.lt.n(ng)) THEN
252 DO j=jstr,jend
253 DO i=istr,iend+1
254 cff=0.5_r8*(pm(i,j)+pm(i-1,j))
255#ifdef MASKING
256 cff=cff*umask(i,j)
257#endif
258#ifdef WET_DRY_NOT_YET
259 cff=cff*umask_wet(i,j)
260#endif
261 drdx(i,j,k2)=cff*(pden(i ,j,k+1)- &
262 & pden(i-1,j,k+1))
263 tl_drdx(i,j,k2)=cff*(tl_pden(i ,j,k+1)- &
264 & tl_pden(i-1,j,k+1))
265#if defined TS_MIX_STABILITY
266 dtdx(i,j,k2)=cff*(0.75_r8*(t(i ,j,k+1,nrhs,itrc)- &
267 & t(i-1,j,k+1,nrhs,itrc))+ &
268 & 0.25_r8*(t(i ,j,k+1,nstp,itrc)- &
269 & t(i-1,j,k+1,nstp,itrc)))
270 tl_dtdx(i,j,k2)=cff* &
271 & (0.75_r8*(tl_t(i ,j,k+1,nrhs,itrc)- &
272 & tl_t(i-1,j,k+1,nrhs,itrc))+ &
273 & 0.25_r8*(tl_t(i ,j,k+1,nstp,itrc)- &
274 & tl_t(i-1,j,k+1,nstp,itrc)))
275#elif defined TS_MIX_CLIMA
276 IF (ltracerclm(itrc,ng)) THEN
277 dtdx(i,j,k2)=cff*((t(i ,j,k+1,nrhs,itrc)- &
278 & tclm(i ,j,k+1,itrc))- &
279 & (t(i-1,j,k+1,nrhs,itrc)- &
280 & tclm(i-1,j,k+1,itrc)))
281 ELSE
282 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
283 & t(i-1,j,k+1,nrhs,itrc))
284 END IF
285 tl_dtdx(i,j,k2)=cff*(tl_t(i ,j,k+1,nrhs,itrc)- &
286 & tl_t(i-1,j,k+1,nrhs,itrc))
287#else
288 dtdx(i,j,k2)=cff*(t(i ,j,k+1,nrhs,itrc)- &
289 & t(i-1,j,k+1,nrhs,itrc))
290 tl_dtdx(i,j,k2)=cff*(tl_t(i ,j,k+1,nrhs,itrc)- &
291 & tl_t(i-1,j,k+1,nrhs,itrc))
292#endif
293 END DO
294 END DO
295 DO j=jstr,jend+1
296 DO i=istr,iend
297 cff=0.5_r8*(pn(i,j)+pn(i,j-1))
298#ifdef MASKING
299 cff=cff*vmask(i,j)
300#endif
301#ifdef WET_DRY_NOT_YET
302 cff=cff*vmask_wet(i,j)
303#endif
304 drde(i,j,k2)=cff*(pden(i,j ,k+1)- &
305 & pden(i,j-1,k+1))
306 tl_drde(i,j,k2)=cff*(tl_pden(i,j ,k+1)- &
307 & tl_pden(i,j-1,k+1))
308#if defined TS_MIX_STABILITY
309 dtde(i,j,k2)=cff*(0.75_r8*(t(i,j ,k+1,nrhs,itrc)- &
310 & t(i,j-1,k+1,nrhs,itrc))+ &
311 & 0.25_r8*(t(i,j ,k+1,nstp,itrc)- &
312 & t(i,j-1,k+1,nstp,itrc)))
313 tl_dtde(i,j,k2)=cff* &
314 & (0.75_r8*(tl_t(i,j ,k+1,nrhs,itrc)- &
315 & tl_t(i,j-1,k+1,nrhs,itrc))+ &
316 & 0.25_r8*(tl_t(i,j ,k+1,nstp,itrc)- &
317 & tl_t(i,j-1,k+1,nstp,itrc)))
318#elif defined TS_MIX_CLIMA
319 IF (ltracerclm(itrc,ng)) THEN
320 dtde(i,j,k2)=cff*((t(i,j ,k+1,nrhs,itrc)- &
321 & tclm(i,j ,k+1,itrc))- &
322 & (t(i,j-1,k+1,nrhs,itrc)- &
323 & tclm(i,j-1,k+1,itrc)))
324 ELSE
325 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
326 & t(i,j-1,k+1,nrhs,itrc))
327 END IF
328 tl_dtde(i,j,k2)=cff*(tl_t(i,j ,k+1,nrhs,itrc)- &
329 & tl_t(i,j-1,k+1,nrhs,itrc))
330#else
331 dtde(i,j,k2)=cff*(t(i,j ,k+1,nrhs,itrc)- &
332 & t(i,j-1,k+1,nrhs,itrc))
333 tl_dtde(i,j,k2)=cff*(tl_t(i,j ,k+1,nrhs,itrc)- &
334 & tl_t(i,j-1,k+1,nrhs,itrc))
335#endif
336 END DO
337 END DO
338 END IF
339 IF ((k.eq.0).or.(k.eq.n(ng))) THEN
340 DO j=jstr-1,jend+1
341 DO i=istr-1,iend+1
342 dtdr(i,j,k2)=0.0_r8
343 tl_dtdr(i,j,k2)=0.0_r8
344 fs(i,j,k2)=0.0_r8
345 tl_fs(i,j,k2)=0.0_r8
346 END DO
347 END DO
348 ELSE
349 DO j=jstr-1,jend+1
350 DO i=istr-1,iend+1
351#if defined TS_MIX_MAX_SLOPE
352 cff1=sqrt(drdx(i,j,k2)**2+drdx(i+1,j,k2)**2+ &
353 & drdx(i,j,k1)**2+drdx(i+1,j,k1)**2+ &
354 & drde(i,j,k2)**2+drde(i,j+1,k2)**2+ &
355 & drde(i,j,k1)**2+drde(i,j+1,k1)**2)
356 IF (cff1.ne.0.0_r8) THEN
357 tl_cff1=(drdx(i ,j,k2)*tl_drdx(i ,j,k2)+ &
358 & drdx(i+1,j,k2)*tl_drdx(i+1,j,k2)+ &
359 & drdx(i ,j,k1)*tl_drdx(i ,j,k1)+ &
360 & drdx(i+1,j,k1)*tl_drdx(i+1,j,k1)+ &
361 & drde(i,j ,k2)*tl_drde(i,j ,k2)+ &
362 & drde(i,j+1,k2)*tl_drde(i,j+1,k2)+ &
363 & drde(i,j ,k1)*tl_drde(i,j ,k1)+ &
364 & drde(i,j+1,k1)*tl_drde(i,j+1,k1))/cff1
365 ELSE
366 tl_cff1=0.0_r8
367 END IF
368 cff2=0.25_r8*slope_max* &
369 & (z_r(i,j,k+1)-z_r(i,j,k))*cff1
370 tl_cff2=0.25_r8*slope_max* &
371 & ((tl_z_r(i,j,k+1)-tl_z_r(i,j,k))*cff1+ &
372 & (z_r(i,j,k+1)-z_r(i,j,k))*tl_cff1)
373 cff3=max(pden(i,j,k)-pden(i,j,k+1),small)
374 tl_cff3=(0.5_r8+sign(0.5_r8,pden(i,j,k)-pden(i,j,k+1)- &
375 & small))* &
376 & (tl_pden(i,j,k)-tl_pden(i,j,k+1))
377 cff4=max(cff2,cff3)
378 tl_cff4=(0.5_r8+sign(0.5_r8,cff2-cff3))*tl_cff2+ &
379 & (0.5_r8-sign(0.5_r8,cff2-cff3))*tl_cff3
380 cff=-1.0_r8/cff4
381 tl_cff=cff*cff*tl_cff4
382#elif defined TS_MIX_MIN_STRAT
383 cff1=max(pden(i,j,k)-pden(i,j,k+1), &
384 & strat_min*(z_r(i,j,k+1)-z_r(i,j,k)))
385 tl_cff1=(0.5_r8+sign(0.5_r8, &
386 & pden(i,j,k)-pden(i,j,k+1)- &
387 & strat_min*(z_r(i,j,k+1)- &
388 & z_r(i,j,k ))))* &
389 & (tl_pden(i,j,k)-tl_pden(i,j,k+1))+ &
390 & (0.5_r8-sign(0.5_r8, &
391 & pden(i,j,k)-pden(i,j,k+1)- &
392 & strat_min*(z_r(i,j,k+1)- &
393 & z_r(i,j,k ))))* &
394 & (strat_min*(tl_z_r(i,j,k+1)-tl_z_r(i,j,k )))
395 cff=-1.0_r8/cff1
396 tl_cff=cff*cff*tl_cff1
397#else
398 cff1=max(pden(i,j,k)-pden(i,j,k+1),eps)
399 tl_cff1=(0.5_r8+sign(0.5_r8, &
400 & pden(i,j,k)-pden(i,j,k+1)-eps))* &
401 & (tl_pden(i,j,k)-tl_pden(i,j,k+1))
402 cff=-1.0_r8/cff1
403 tl_cff=cff*cff*tl_cff1
404#endif
405#if defined TS_MIX_STABILITY
406 dtdr(i,j,k2)=cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
407 & t(i,j,k ,nrhs,itrc))+ &
408 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
409 & t(i,j,k ,nstp,itrc)))
410 tl_dtdr(i,j,k2)=tl_cff* &
411 & (0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
412 & t(i,j,k ,nrhs,itrc))+ &
413 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
414 & t(i,j,k ,nstp,itrc)))+ &
415 & cff* &
416 & (0.75_r8*(tl_t(i,j,k+1,nrhs,itrc)- &
417 & tl_t(i,j,k ,nrhs,itrc))+ &
418 & 0.25_r8*(tl_t(i,j,k+1,nstp,itrc)- &
419 & tl_t(i,j,k ,nstp,itrc)))
420#elif defined TS_MIX_CLIMA
421 IF (ltracerclm(itrc,ng)) THEN
422 dtdr(i,j,k2)=cff*((t(i,j,k+1,nrhs,itrc)- &
423 & tclm(i,j,k+1,itrc))- &
424 & (t(i,j,k ,nrhs,itrc)- &
425 & tclm(i,j,k ,itrc)))
426 tl_dtdr(i,j,k2)=tl_cff*((t(i,j,k+1,nrhs,itrc)- &
427 & tclm(i,j,k+1,itrc))- &
428 & (t(i,j,k ,nrhs,itrc)- &
429 & tclm(i,j,k ,itrc)))+ &
430 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
431 & tl_t(i,j,k ,nrhs,itrc))
432 ELSE
433 dtdr(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
434 & t(i,j,k ,nrhs,itrc))
435 tl_dtdr(i,j,k2)=tl_cff*(t(i,j,k+1,nrhs,itrc)- &
436 & t(i,j,k ,nrhs,itrc))+ &
437 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
438 & tl_t(i,j,k ,nrhs,itrc))
439 END IF
440#else
441 dtdr(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
442 & t(i,j,k ,nrhs,itrc))
443 tl_dtdr(i,j,k2)=tl_cff*(t(i,j,k+1,nrhs,itrc)- &
444 & t(i,j,k ,nrhs,itrc))+ &
445 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
446 & tl_t(i,j,k ,nrhs,itrc))
447#endif
448 fs(i,j,k2)=cff*(z_r(i,j,k+1)-z_r(i,j,k))
449 tl_fs(i,j,k2)=tl_cff*(z_r(i,j,k+1)-z_r(i,j,k))+ &
450 & cff*(tl_z_r(i,j,k+1)-tl_z_r(i,j,k))
451 END DO
452 END DO
453 END IF
454!
455! Compute components of the rotated tracer flux (T m4/s) along
456! isopycnic surfaces.
457!
458 IF (k.gt.0) THEN
459 DO j=jstr,jend
460 DO i=istr,iend+1
461#ifdef DIFF_3DCOEF
462 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
463 & on_u(i,j)
464#else
465 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
466 & on_u(i,j)
467#endif
468!^ FX(i,j)=cff* &
469!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
470!^ & (dTdx(i,j,k1)- &
471!^ & 0.5_r8*(MAX(dRdx(i,j,k1),0.0_r8)* &
472!^ & (dTdr(i-1,j,k1)+ &
473!^ & dTdr(i ,j,k2))+ &
474!^ & MIN(dRdx(i,j,k1),0.0_r8)* &
475!^ & (dTdr(i-1,j,k2)+ &
476!^ & dTdr(i,j,k1))))
477!^
478 tl_fx(i,j)=cff* &
479 & (((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
480 & (dtdx(i,j,k1)- &
481 & 0.5_r8*(max(drdx(i,j,k1),0.0_r8)* &
482 & (dtdr(i-1,j,k1)+ &
483 & dtdr(i ,j,k2))+ &
484 & min(drdx(i,j,k1),0.0_r8)* &
485 & (dtdr(i-1,j,k2)+ &
486 & dtdr(i ,j,k1)))))+ &
487 & ((hz(i,j,k)+hz(i-1,j,k))* &
488 & (tl_dtdx(i,j,k1)- &
489 & 0.5_r8*(max(drdx(i,j,k1),0.0_r8)* &
490 & (tl_dtdr(i-1,j,k1)+ &
491 & tl_dtdr(i ,j,k2))+ &
492 & min(drdx(i,j,k1),0.0_r8)* &
493 & (tl_dtdr(i-1,j,k2)+ &
494 & tl_dtdr(i ,j,k1)))- &
495 & 0.5_r8*((0.5_r8+ &
496 & sign(0.5_r8, drdx(i,j,k1)))* &
497 & tl_drdx(i,j,k1)* &
498 & (dtdr(i-1,j,k1)+dtdr(i,j,k2))+ &
499 & (0.5_r8+ &
500 & sign(0.5_r8,-drdx(i,j,k1)))* &
501 & tl_drdx(i,j,k1)* &
502 & (dtdr(i-1,j,k2)+dtdr(i,j,k1))))))
503 END DO
504 END DO
505 DO j=jstr,jend+1
506 DO i=istr,iend
507#ifdef DIFF_3DCOEF
508 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
509 & om_v(i,j)
510#else
511 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
512 & om_v(i,j)
513#endif
514!^ FE(i,j)=cff* &
515!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
516!^ & (dTde(i,j,k1)- &
517!^ & 0.5_r8*(MAX(dRde(i,j,k1),0.0_r8)* &
518!^ & (dTdr(i,j-1,k1)+ &
519!^ & dTdr(i,j ,k2))+ &
520!^ & MIN(dRde(i,j,k1),0.0_r8)* &
521!^ & (dTdr(i,j-1,k2)+ &
522!^ & dTdr(i,j ,k1))))
523!^
524 tl_fe(i,j)=cff* &
525 & (((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
526 & (dtde(i,j,k1)- &
527 & 0.5_r8*(max(drde(i,j,k1),0.0_r8)* &
528 & (dtdr(i,j-1,k1)+ &
529 & dtdr(i,j ,k2))+ &
530 & min(drde(i,j,k1),0.0_r8)* &
531 & (dtdr(i,j-1,k2)+ &
532 & dtdr(i,j ,k1)))))+ &
533 & ((hz(i,j,k)+hz(i,j-1,k))* &
534 & (tl_dtde(i,j,k1)- &
535 & 0.5_r8*(max(drde(i,j,k1),0.0_r8)* &
536 & (tl_dtdr(i,j-1,k1)+ &
537 & tl_dtdr(i,j ,k2))+ &
538 & min(drde(i,j,k1),0.0_r8)* &
539 & (tl_dtdr(i,j-1,k2)+ &
540 & tl_dtdr(i,j ,k1)))- &
541 & 0.5_r8*((0.5_r8+ &
542 & sign(0.5_r8, drde(i,j,k1)))* &
543 & tl_drde(i,j,k1)* &
544 & (dtdr(i,j-1,k1)+dtdr(i,j,k2))+ &
545 & (0.5_r8+ &
546 & sign(0.5_r8,-drde(i,j,k1)))* &
547 & tl_drde(i,j,k1)* &
548 & (dtdr(i,j-1,k2)+dtdr(i,j,k1))))))
549 END DO
550 END DO
551 IF (k.lt.n(ng)) THEN
552 DO j=jstr,jend
553 DO i=istr,iend
554 cff1=max(drdx(i ,j,k1),0.0_r8)
555 cff2=max(drdx(i+1,j,k2),0.0_r8)
556 cff3=min(drdx(i ,j,k2),0.0_r8)
557 cff4=min(drdx(i+1,j,k1),0.0_r8)
558 tl_cff1=(0.5_r8+sign(0.5_r8, drdx(i ,j,k1)))* &
559 & tl_drdx(i ,j,k1)
560 tl_cff2=(0.5_r8+sign(0.5_r8, drdx(i+1,j,k2)))* &
561 & tl_drdx(i+1,j,k2)
562 tl_cff3=(0.5_r8+sign(0.5_r8,-drdx(i ,j,k2)))* &
563 & tl_drdx(i ,j,k2)
564 tl_cff4=(0.5_r8+sign(0.5_r8,-drdx(i+1,j,k1)))* &
565 & tl_drdx(i+1,j,k1)
566 cff=cff1*(cff1*dtdr(i,j,k2)-dtdx(i ,j,k1))+ &
567 & cff2*(cff2*dtdr(i,j,k2)-dtdx(i+1,j,k2))+ &
568 & cff3*(cff3*dtdr(i,j,k2)-dtdx(i ,j,k2))+ &
569 & cff4*(cff4*dtdr(i,j,k2)-dtdx(i+1,j,k1))
570 tl_cff=tl_cff1*(cff1*dtdr(i ,j,k2)- &
571 & dtdx(i ,j,k1))+ &
572 & tl_cff2*(cff2*dtdr(i,j,k2)- &
573 & dtdx(i+1,j,k2))+ &
574 & tl_cff3*(cff3*dtdr(i,j,k2)- &
575 & dtdx(i ,j,k2))+ &
576 & tl_cff4*(cff4*dtdr(i,j,k2)- &
577 & dtdx(i+1,j,k1))+ &
578 & cff1*(tl_cff1*dtdr(i,j,k2)+ &
579 & cff1*tl_dtdr(i,j,k2)- &
580 & tl_dtdx(i ,j,k1))+ &
581 & cff2*(tl_cff2*dtdr(i,j,k2)+ &
582 & cff2*tl_dtdr(i,j,k2)- &
583 & tl_dtdx(i+1,j,k2))+ &
584 & cff3*(tl_cff3*dtdr(i,j,k2)+ &
585 & cff3*tl_dtdr(i,j,k2)- &
586 & tl_dtdx(i ,j,k2))+ &
587 & cff4*(tl_cff4*dtdr(i,j,k2)+ &
588 & cff4*tl_dtdr(i,j,k2)- &
589 & tl_dtdx(i+1,j,k1))
590 cff1=max(drde(i,j ,k1),0.0_r8)
591 cff2=max(drde(i,j+1,k2),0.0_r8)
592 cff3=min(drde(i,j ,k2),0.0_r8)
593 cff4=min(drde(i,j+1,k1),0.0_r8)
594 tl_cff1=(0.5_r8+sign(0.5_r8, drde(i,j ,k1)))* &
595 & tl_drde(i,j ,k1)
596 tl_cff2=(0.5_r8+sign(0.5_r8, drde(i,j+1,k2)))* &
597 & tl_drde(i,j+1,k2)
598 tl_cff3=(0.5_r8+sign(0.5_r8,-drde(i,j ,k2)))* &
599 & tl_drde(i,j ,k2)
600 tl_cff4=(0.5_r8+sign(0.5_r8,-drde(i,j+1,k1)))* &
601 & tl_drde(i,j+1,k1)
602 cff=cff+ &
603 & cff1*(cff1*dtdr(i,j,k2)-dtde(i,j ,k1))+ &
604 & cff2*(cff2*dtdr(i,j,k2)-dtde(i,j+1,k2))+ &
605 & cff3*(cff3*dtdr(i,j,k2)-dtde(i,j ,k2))+ &
606 & cff4*(cff4*dtdr(i,j,k2)-dtde(i,j+1,k1))
607 tl_cff=tl_cff+ &
608 & tl_cff1*(cff1*dtdr(i,j,k2)- &
609 & dtde(i,j ,k1))+ &
610 & tl_cff2*(cff2*dtdr(i,j,k2)- &
611 & dtde(i,j+1,k2))+ &
612 & tl_cff3*(cff3*dtdr(i,j,k2)- &
613 & dtde(i,j ,k2))+ &
614 & tl_cff4*(cff4*dtdr(i,j,k2)- &
615 & dtde(i,j+1,k1))+ &
616 & cff1*(tl_cff1*dtdr(i,j,k2)+ &
617 & cff1*tl_dtdr(i,j,k2)- &
618 & tl_dtde(i,j ,k1))+ &
619 & cff2*(tl_cff2*dtdr(i,j,k2)+ &
620 & cff2*tl_dtdr(i,j,k2)- &
621 & tl_dtde(i,j+1,k2))+ &
622 & cff3*(tl_cff3*dtdr(i,j,k2)+ &
623 & cff3*tl_dtdr(i,j,k2)- &
624 & tl_dtde(i,j ,k2))+ &
625 & cff4*(tl_cff4*dtdr(i,j,k2)+ &
626 & cff4*tl_dtdr(i,j,k2)- &
627 & tl_dtde(i,j+1,k1))
628#ifdef DIFF_3DCOEF
629!^ FS(i,j,k2)=0.5_r8*cff*diff3d_r(i,j,k)*FS(i,j,k2)
630!^
631 tl_fs(i,j,k2)=0.5_r8*diff3d_r(i,j,k)* &
632 & (tl_cff*fs(i,j,k2)+ &
633 & cff*tl_fs(i,j,k2))
634#else
635!^ FS(i,j,k2)=0.5_r8*cff*diff2(i,j,itrc)*FS(i,j,k2)
636!^
637 tl_fs(i,j,k2)=0.5_r8*diff2(i,j,itrc)* &
638 & (tl_cff*fs(i,j,k2)+ &
639 & cff*tl_fs(i,j,k2))
640#endif
641 END DO
642 END DO
643 END IF
644!
645! Time-step harmonic, isopycnic diffusion term (m Tunits).
646!
647 DO j=jstr,jend
648 DO i=istr,iend
649!^ cff=dt(ng)*pm(i,j)*pn(i,j)* &
650!^ & (FX(i+1,j)-FX(i,j)+ &
651!^ & FE(i,j+1)-FE(i,j))+ &
652!^ & dt(ng)*(FS(i,j,k2)-FS(i,j,k1))
653!^
654 tl_cff=dt(ng)*pm(i,j)*pn(i,j)* &
655 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
656 & tl_fe(i,j+1)-tl_fe(i,j))+ &
657 & dt(ng)*(tl_fs(i,j,k2)-tl_fs(i,j,k1))
658!^ t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff
659!^
660 tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+tl_cff
661#ifdef DIAGNOSTICS_TS
662!! DiaTwrk(i,j,k,itrc,iThdif)=cff
663#endif
664 END DO
665 END DO
666 END IF
667 END DO k_loop
668 END DO t_loop
669!
670 RETURN

References mod_scalars::dt, and mod_scalars::ltracerclm.

◆ tl_t3dmix2_s_tile()

subroutine tl_t3dmix2_mod::tl_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,n(ng)), intent(in) tl_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),3,nt(ng)), intent(in) t,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),3,nt(ng)), intent(inout) tl_t )
private

Definition at line 93 of file tl_t3dmix2_s.h.

117!***********************************************************************
118!
119 USE mod_param
120 USE mod_scalars
121!
122! Imported variable declarations.
123!
124 integer, intent(in) :: ng, tile
125 integer, intent(in) :: LBi, UBi, LBj, UBj
126 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
127 integer, intent(in) :: nrhs, nstp, nnew
128
129#ifdef ASSUMED_SHAPE
130# ifdef MASKING
131 real(r8), intent(in) :: umask(LBi:,LBj:)
132 real(r8), intent(in) :: vmask(LBi:,LBj:)
133# endif
134# ifdef WET_DRY_NOT_YET
135 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
136 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
137# endif
138# ifdef DIFF_3DCOEF
139 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
140# else
141 real(r8), intent(in) :: diff2(LBi:,LBj:,:)
142# endif
143 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
144 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
145 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
146 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
147 real(r8), intent(in) :: pm(LBi:,LBj:)
148 real(r8), intent(in) :: pn(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(in) :: t(LBi:,LBj:,:,:,:)
156
157 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
158#else
159# ifdef MASKING
160 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
161 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
162# endif
163# ifdef WET_DRY_NOT_YET
164 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
165 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
166# endif
167# ifdef DIFF_3DCOEF
168 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
169# else
170 real(r8), intent(in) :: diff2(LBi:UBi,LBj:UBj,NT(ng))
171# endif
172 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
173 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
174 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
175 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
176 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
177 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
178# ifdef TS_MIX_CLIMA
179 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
180# endif
181# ifdef DIAGNOSTICS_TS
182!! real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
183!! & NDT)
184# endif
185 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
186
187 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
188#endif
189!
190! Local variable declarations.
191!
192 integer :: i, itrc, j, k
193
194 real(r8) :: cff, cff1, tl_cff, tl_cff1
195
196 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FE
197 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FX
198
199#include "set_bounds.h"
200!
201!-----------------------------------------------------------------------
202! Compute tangent linear horizontal harmonic diffusion along constant
203! S-surfaces.
204!-----------------------------------------------------------------------
205!
206 DO itrc=1,nt(ng)
207 DO k=1,n(ng)
208!
209! Compute XI- and ETA-components of diffusive tracer flux (T m3/s).
210!
211 DO j=jstr,jend
212 DO i=istr,iend+1
213#ifdef DIFF_3DCOEF
214 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
215 & pmon_u(i,j)
216#else
217 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
218 & pmon_u(i,j)
219#endif
220#if defined TS_MIX_STABILITY
221!^ FX(i,j)=cff* &
222!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
223!^ & (0.75_r8*(t(i ,j,k,nrhs,itrc)- &
224!^ & t(i-1,j,k,nrhs,itrc))+ &
225!^ & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
226!^ & t(i-1,j,k,nstp,itrc)))
227!^
228 tl_fx(i,j)=cff* &
229 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
230 & (0.75_r8*(t(i ,j,k,nrhs,itrc)- &
231 & t(i-1,j,k,nrhs,itrc))+ &
232 & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
233 & t(i-1,j,k,nstp,itrc)))+ &
234 & (hz(i,j,k)+hz(i-1,j,k))* &
235 & (0.75_r8*(tl_t(i ,j,k,nrhs,itrc)- &
236 & tl_t(i-1,j,k,nrhs,itrc))+ &
237 & 0.25_r8*(tl_t(i ,j,k,nstp,itrc)- &
238 & tl_t(i-1,j,k,nstp,itrc))))
239#elif defined TS_MIX_CLIMA
240 IF (ltracerclm(itrc,ng)) THEN
241!^ FX(i,j)=cff* &
242!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
243!^ & ((t(i ,j,k,nrhs,itrc)-tclm(i ,j,k,itrc))- &
244!^ & (t(i-1,j,k,nrhs,itrc)-tclm(i-1,j,k,itrc)))
245!^
246 tl_fx(i,j)=cff* &
247 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
248 & ((t(i ,j,k,nrhs,itrc)- &
249 & tclm(i ,j,k,itrc))- &
250 & (t(i-1,j,k,nrhs,itrc)- &
251 & tclm(i-1,j,k,itrc)))+ &
252 & (hz(i,j,k)+hz(i-1,j,k))* &
253 & (tl_t(i ,j,k,nrhs,itrc)- &
254 & tl_t(i-1,j,k,nrhs,itrc)))
255 ELSE
256!^ FX(i,j)=cff* &
257!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
258!^ & (t(i,j,k,nrhs,itrc)-t(i-1,j,k,nrhs,itrc))
259!^
260 tl_fx(i,j)=cff* &
261 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
262 & (t(i ,j,k,nrhs,itrc)- &
263 & t(i-1,j,k,nrhs,itrc))+ &
264 & (hz(i,j,k)+hz(i-1,j,k))* &
265 & (tl_t(i ,j,k,nrhs,itrc)- &
266 & tl_t(i-1,j,k,nrhs,itrc)))
267 END IF
268#else
269!^ FX(i,j)=cff* &
270!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
271!^ & (t(i,j,k,nrhs,itrc)-t(i-1,j,k,nrhs,itrc))
272!^
273 tl_fx(i,j)=cff* &
274 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
275 & (t(i ,j,k,nrhs,itrc)- &
276 & t(i-1,j,k,nrhs,itrc))+ &
277 & (hz(i,j,k)+hz(i-1,j,k))* &
278 & (tl_t(i ,j,k,nrhs,itrc)- &
279 & tl_t(i-1,j,k,nrhs,itrc)))
280#endif
281#ifdef MASKING
282!^ FX(i,j)=FX(i,j)*umask(i,j)
283!^
284 tl_fx(i,j)=tl_fx(i,j)*umask(i,j)
285#endif
286#ifdef WET_DRY_NOT_YET
287 fx(i,j)=fx(i,j)*umask_wet(i,j)
288#endif
289 END DO
290 END DO
291 DO j=jstr,jend+1
292 DO i=istr,iend
293#ifdef DIFF_3DCOEF
294 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
295 & pnom_v(i,j)
296#else
297 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
298 & pnom_v(i,j)
299#endif
300#if defined TS_MIX_STABILITY
301!^ FE(i,j)=cff* &
302!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
303!^ & (0.75_r8*(t(i,j ,k,nrhs,itrc)- &
304!^ & t(i,j-1,k,nrhs,itrc))+ &
305!^ & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
306!^ & t(i,j-1,k,nstp,itrc)))
307!^
308 tl_fe(i,j)=cff* &
309 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
310 & (0.75_r8*(t(i,j ,k,nrhs,itrc)- &
311 & t(i,j-1,k,nrhs,itrc))+ &
312 & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
313 & t(i,j-1,k,nstp,itrc)))+ &
314 & (hz(i,j,k)+hz(i,j-1,k))* &
315 & (0.75_r8*(tl_t(i,j ,k,nrhs,itrc)- &
316 & tl_t(i,j-1,k,nrhs,itrc))+ &
317 & 0.25_r8*(tl_t(i,j ,k,nstp,itrc)- &
318 & tl_t(i,j-1,k,nstp,itrc))))
319#elif defined TS_MIX_CLIMA
320 IF (ltracerclm(itrc,ng)) THEN
321!^ FE(i,j)=cff* &
322!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
323!^ & ((t(i,j ,k,nrhs,itrc)-tclm(i,j ,k,itrc))- &
324!^ & (t(i,j-1,k,nrhs,itrc)-tclm(i,j-1,k,itrc)))
325!^
326 tl_fe(i,j)=cff* &
327 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
328 & ((t(i,j ,k,nrhs,itrc)- &
329 & tclm(i,j ,k,itrc))- &
330 & (t(i,j-1,k,nrhs,itrc)- &
331 & tclm(i,j-1,k,itrc)))+ &
332 & (hz(i,j,k)+hz(i,j-1,k))* &
333 & (tl_t(i,j ,k,nrhs,itrc)- &
334 & tl_t(i,j-1,k,nrhs,itrc)))
335 ELSE
336!^ FE(i,j)=cff* &
337!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
338!^ & (t(i,j,k,nrhs,itrc)-t(i,j-1,k,nrhs,itrc))
339!^
340 tl_fe(i,j)=cff* &
341 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
342 & (t(i,j ,k,nrhs,itrc)- &
343 & t(i,j-1,k,nrhs,itrc))+ &
344 & (hz(i,j,k)+hz(i,j-1,k))* &
345 & (tl_t(i,j ,k,nrhs,itrc)- &
346 & tl_t(i,j-1,k,nrhs,itrc)))
347 END IF
348#else
349!^ FE(i,j)=cff* &
350!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
351!^ & (t(i,j,k,nrhs,itrc)-t(i,j-1,k,nrhs,itrc))
352!^
353 tl_fe(i,j)=cff* &
354 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
355 & (t(i,j ,k,nrhs,itrc)- &
356 & t(i,j-1,k,nrhs,itrc))+ &
357 & (hz(i,j,k)+hz(i,j-1,k))* &
358 & (tl_t(i,j ,k,nrhs,itrc)- &
359 & tl_t(i,j-1,k,nrhs,itrc)))
360#endif
361#ifdef MASKING
362!^ FE(i,j)=FE(i,j)*vmask(i,j)
363!^
364 tl_fe(i,j)=tl_fe(i,j)*vmask(i,j)
365#endif
366#ifdef WET_DRY_NOT_YET
367 fe(i,j)=fe(i,j)*vmask_wet(i,j)
368#endif
369 END DO
370 END DO
371!
372! Time-step harmonic, S-surfaces diffusion term (m Tunits).
373!
374 DO j=jstr,jend
375 DO i=istr,iend
376!^ cff=dt(ng)*pm(i,j)*pn(i,j)* &
377!^ & (FX(i+1,j)-FX(i,j)+ &
378!^ & FE(i,j+1)-FE(i,j))
379!^
380 tl_cff=dt(ng)*pm(i,j)*pn(i,j)* &
381 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
382 & tl_fe(i,j+1)-tl_fe(i,j))
383!^ t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff
384!^
385 tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+tl_cff
386#ifdef DIAGNOSTICS_TS
387!! DiaTwrk(i,j,k,itrc,iThdif)=cff
388#endif
389 END DO
390 END DO
391 END DO
392 END DO
393!
394 RETURN

References mod_scalars::dt, and mod_scalars::ltracerclm.