ROMS
Loading...
Searching...
No Matches
tl_t3dmix2_geo.h
Go to the documentation of this file.
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! This subroutine computes tangent linear horizontal harmonic mixing !
11! of tracers along geopotential surfaces. !
12! !
13! BASIC STATE variables needed: diff2, Hz, t, z_r !
14! !
15!=======================================================================
16!
17 implicit none
18!
19 PRIVATE
20 PUBLIC tl_t3dmix2
21!
22 CONTAINS
23!
24!***********************************************************************
25 SUBROUTINE tl_t3dmix2 (ng, tile)
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
92 END SUBROUTINE tl_t3dmix2
93!
94!***********************************************************************
95 SUBROUTINE tl_t3dmix2_geo_tile (ng, tile, &
96 & LBi, UBi, LBj, UBj, &
97 & IminS, ImaxS, JminS, JmaxS, &
98 & nrhs, nstp, nnew, &
99#ifdef MASKING
100 & umask, vmask, &
101#endif
102#ifdef WET_DRY_NOT_YET
103 & umask_wet, vmask_wet, &
104#endif
105 & om_v, on_u, pm, pn, &
106 & Hz, tl_Hz, &
107 & z_r, tl_z_r, &
108#ifdef DIFF_3DCOEF
109 & diff3d_r, &
110#else
111 & diff2, &
112#endif
113#ifdef TS_MIX_CLIMA
114 & tclm, &
115#endif
116#ifdef DIAGNOSTICS_TS
117!! & DiaTwrk, &
118#endif
119 & t, tl_t)
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
608 END SUBROUTINE tl_t3dmix2_geo_tile
609
610 END MODULE tl_t3dmix2_mod
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
real(dp), dimension(:), allocatable dt
logical, dimension(:,:), allocatable ltracerclm
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
integer, dimension(:), allocatable nstp
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, public tl_t3dmix2(ng, tile)
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