ROMS
Loading...
Searching...
No Matches
rp_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 representers tangent linear horizontal !
11! harmonic mixing 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 rp_t3dmix2
21!
22 CONTAINS
23!
24!***********************************************************************
25 SUBROUTINE rp_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, irpm, 25, __line__, myfile)
53#endif
54 CALL rp_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, irpm, 25, __line__, myfile)
89#endif
90!
91 RETURN
92 END SUBROUTINE rp_t3dmix2
93!
94!***********************************************************************
95 SUBROUTINE rp_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#ifdef TL_IOMS
348 & 2.0_r8*cff
349#endif
350#if defined TS_MIX_STABILITY
351 dtdz(i,j,k2)=cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
352 & t(i,j,k ,nrhs,itrc))+ &
353 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
354 & t(i,j,k ,nstp,itrc)))
355 tl_dtdz(i,j,k2)=tl_cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
356 & t(i,j,k ,nrhs,itrc))+ &
357 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
358 & t(i,j,k ,nstp,itrc)))+&
359 & cff*(0.75_r8*(tl_t(i,j,k+1,nrhs,itrc)- &
360 & tl_t(i,j,k ,nrhs,itrc))+ &
361 & 0.25_r8*(tl_t(i,j,k+1,nstp,itrc)- &
362 & tl_t(i,j,k ,nstp,itrc)))-&
363# ifdef TL_IOMS
364 & dtdz(i,j,k2)
365# endif
366#elif defined TS_MIX_CLIMA
367 IF (ltracerclm(itrc,ng)) THEN
368 dtdz(i,j,k2)=cff*((t(i,j,k+1,nrhs,itrc)- &
369 & tclm(i,j,k+1,itrc))- &
370 & (t(i,j,k ,nrhs,itrc)- &
371 & tclm(i,j,k ,itrc)))
372 tl_dtdz(i,j,k2)=tl_cff*((t(i,j,k+1,nrhs,itrc)- &
373 & tclm(i,j,k+1,itrc))- &
374 & (t(i,j,k ,nrhs,itrc)- &
375 & tclm(i,j,k ,itrc)))+ &
376 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
377 & tl_t(i,j,k ,nrhs,itrc))- &
378# ifdef TL_IOMS
379 & dtdz(i,j,k2)
380# endif
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# ifdef TL_IOMS
389 & dtdz(i,j,k2)
390# endif
391 END IF
392#else
393 dtdz(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
394 & t(i,j,k ,nrhs,itrc))
395 tl_dtdz(i,j,k2)=tl_cff*(t(i,j,k+1,nrhs,itrc)- &
396 & t(i,j,k ,nrhs,itrc))+ &
397 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
398 & tl_t(i,j,k ,nrhs,itrc))- &
399# ifdef TL_IOMS
400 & dtdz(i,j,k2)
401# endif
402#endif
403 END DO
404 END DO
405 END IF
406!
407! Compute components of the rotated tracer flux (T m3/s) along
408! geopotential surfaces.
409!
410 IF (k.gt.0) THEN
411 DO j=jstr,jend
412 DO i=istr,iend+1
413#ifdef DIFF_3DCOEF
414 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
415 & on_u(i,j)
416#else
417 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
418 & on_u(i,j)
419#endif
420!^ FX(i,j)=cff* &
421!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
422!^ & (dTdx(i,j,k1)- &
423!^ & 0.5_r8*(MIN(dZdx(i,j,k1),0.0_r8)* &
424!^ & (dTdz(i-1,j,k1)+ &
425!^ & dTdz(i ,j,k2))+ &
426!^ & MAX(dZdx(i,j,k1),0.0_r8)* &
427!^ & (dTdz(i-1,j,k2)+ &
428!^ & dTdz(i ,j,k1))))
429!^
430 tl_fx(i,j)=cff* &
431 & (((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
432 & (dtdx(i,j,k1)- &
433 & 0.5_r8*(min(dzdx(i,j,k1),0.0_r8)* &
434 & (dtdz(i-1,j,k1)+ &
435 & dtdz(i ,j,k2))+ &
436 & max(dzdx(i,j,k1),0.0_r8)* &
437 & (dtdz(i-1,j,k2)+ &
438 & dtdz(i ,j,k1)))))+ &
439 & ((hz(i,j,k)+hz(i-1,j,k))* &
440 & (tl_dtdx(i,j,k1)- &
441 & 0.5_r8*(min(dzdx(i,j,k1),0.0_r8)* &
442 & (tl_dtdz(i-1,j,k1)+ &
443 & tl_dtdz(i ,j,k2))+ &
444 & max(dzdx(i,j,k1),0.0_r8)* &
445 & (tl_dtdz(i-1,j,k2)+ &
446 & tl_dtdz(i ,j,k1)))- &
447 & 0.5_r8*((0.5_r8+ &
448 & sign(0.5_r8,-dzdx(i,j,k1)))* &
449 & tl_dzdx(i,j,k1)* &
450 & (dtdz(i-1,j,k1)+dtdz(i,j,k2))+ &
451 & (0.5_r8+ &
452 & sign(0.5_r8, dzdx(i,j,k1)))* &
453 & tl_dzdx(i,j,k1)* &
454 & (dtdz(i-1,j,k2)+dtdz(i,j,k1))))))-&
455#ifdef TL_IOMS
456 & cff* &
457 & (hz(i,j,k)+hz(i-1,j,k))* &
458 & (dtdx(i,j,k1)- &
459 & (min(dzdx(i,j,k1),0.0_r8)* &
460 & (dtdz(i-1,j,k1)+ &
461 & dtdz(i ,j,k2))+ &
462 & max(dzdx(i,j,k1),0.0_r8)* &
463 & (dtdz(i-1,j,k2)+ &
464 & dtdz(i ,j,k1))))
465#endif
466 END DO
467 END DO
468 DO j=jstr,jend+1
469 DO i=istr,iend
470#ifdef DIFF_3DCOEF
471 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
472 & om_v(i,j)
473#else
474 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
475 & om_v(i,j)
476#endif
477!^ FE(i,j)=cff* &
478!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
479!^ & (dTde(i,j,k1)- &
480!^ & 0.5_r8*(MIN(dZde(i,j,k1),0.0_r8)* &
481!^ & (dTdz(i,j-1,k1)+ &
482!^ & dTdz(i,j ,k2))+ &
483!^ & MAX(dZde(i,j,k1),0.0_r8)* &
484!^ & (dTdz(i,j-1,k2)+ &
485!^ & dTdz(i,j ,k1))))
486!^
487 tl_fe(i,j)=cff* &
488 & (((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
489 & (dtde(i,j,k1)- &
490 & 0.5_r8*(min(dzde(i,j,k1),0.0_r8)* &
491 & (dtdz(i,j-1,k1)+ &
492 & dtdz(i,j ,k2))+ &
493 & max(dzde(i,j,k1),0.0_r8)* &
494 & (dtdz(i,j-1,k2)+ &
495 & dtdz(i,j ,k1)))))+ &
496 & ((hz(i,j,k)+hz(i,j-1,k))* &
497 & (tl_dtde(i,j,k1)- &
498 & 0.5_r8*(min(dzde(i,j,k1),0.0_r8)* &
499 & (tl_dtdz(i,j-1,k1)+ &
500 & tl_dtdz(i,j ,k2))+ &
501 & max(dzde(i,j,k1),0.0_r8)* &
502 & (tl_dtdz(i,j-1,k2)+ &
503 & tl_dtdz(i,j ,k1)))- &
504 & 0.5_r8*((0.5_r8+ &
505 & sign(0.5_r8,-dzde(i,j,k1)))* &
506 & tl_dzde(i,j,k1)* &
507 & (dtdz(i,j-1,k1)+dtdz(i,j,k2))+ &
508 & (0.5_r8+ &
509 & sign(0.5_r8, dzde(i,j,k1)))* &
510 & tl_dzde(i,j,k1)* &
511 & (dtdz(i,j-1,k2)+dtdz(i,j,k1))))))-&
512#ifdef TL_IOMS
513 & cff* &
514 & (hz(i,j,k)+hz(i,j-1,k))* &
515 & (dtde(i,j,k1)- &
516 & (min(dzde(i,j,k1),0.0_r8)* &
517 & (dtdz(i,j-1,k1)+ &
518 & dtdz(i,j ,k2))+ &
519 & max(dzde(i,j,k1),0.0_r8)* &
520 & (dtdz(i,j-1,k2)+ &
521 & dtdz(i,j ,k1))))
522#endif
523 END DO
524 END DO
525 IF (k.lt.n(ng)) THEN
526 DO j=jstr,jend
527 DO i=istr,iend
528#ifdef DIFF_3DCOEF
529 cff=0.5_r8*diff3d_r(i,j,k)
530#else
531 cff=0.5_r8*diff2(i,j,itrc)
532#endif
533 cff1=min(dzdx(i ,j,k1),0.0_r8)
534 cff2=min(dzdx(i+1,j,k2),0.0_r8)
535 cff3=max(dzdx(i ,j,k2),0.0_r8)
536 cff4=max(dzdx(i+1,j,k1),0.0_r8)
537 tl_cff1=(0.5_r8+sign(0.5_r8,-dzdx(i ,j,k1)))* &
538 & tl_dzdx(i ,j,k1)
539 tl_cff2=(0.5_r8+sign(0.5_r8,-dzdx(i+1,j,k2)))* &
540 & tl_dzdx(i+1,j,k2)
541 tl_cff3=(0.5_r8+sign(0.5_r8, dzdx(i ,j,k2)))* &
542 & tl_dzdx(i ,j,k2)
543 tl_cff4=(0.5_r8+sign(0.5_r8, dzdx(i+1,j,k1)))* &
544 & tl_dzdx(i+1,j,k1)
545!^ FS(i,j,k2)=cff* &
546!^ & (cff1*(cff1*dTdz(i,j,k2)-dTdx(i ,j,k1))+ &
547!^ & cff2*(cff2*dTdz(i,j,k2)-dTdx(i+1,j,k2))+ &
548!^ & cff3*(cff3*dTdz(i,j,k2)-dTdx(i ,j,k2))+ &
549!^ & cff4*(cff4*dTdz(i,j,k2)-dTdx(i+1,j,k1)))
550!^
551 tl_fs(i,j,k2)=cff* &
552 & (tl_cff1*(cff1*dtdz(i,j,k2)- &
553 & dtdx(i ,j,k1))+ &
554 & tl_cff2*(cff2*dtdz(i,j,k2)- &
555 & dtdx(i+1,j,k2))+ &
556 & tl_cff3*(cff3*dtdz(i,j,k2)- &
557 & dtdx(i ,j,k2))+ &
558 & tl_cff4*(cff4*dtdz(i,j,k2)- &
559 & dtdx(i+1,j,k1))+ &
560 & cff1*(tl_cff1*dtdz(i,j,k2)+ &
561 & cff1*tl_dtdz(i,j,k2)- &
562 & tl_dtdx(i ,j,k1))+ &
563 & cff2*(tl_cff2*dtdz(i,j,k2)+ &
564 & cff2*tl_dtdz(i,j,k2)- &
565 & tl_dtdx(i+1,j,k2))+ &
566 & cff3*(tl_cff3*dtdz(i,j,k2)+ &
567 & cff3*tl_dtdz(i,j,k2)- &
568 & tl_dtdx(i ,j,k2))+ &
569 & cff4*(tl_cff4*dtdz(i,j,k2)+ &
570 & cff4*tl_dtdz(i,j,k2)- &
571 & tl_dtdx(i+1,j,k1)))- &
572#ifdef TL_IOMS
573 & cff* &
574 & (cff1*(2.0_r8*cff1*dtdz(i,j,k2)- &
575 & dtdx(i,j,k1))+ &
576 & cff2*(2.0_r8*cff2*dtdz(i ,j,k2)- &
577 & dtdx(i+1,j,k2))+ &
578 & cff3*(2.0_r8*cff3*dtdz(i,j,k2)- &
579 & dtdx(i,j,k2))+ &
580 & cff4*(2.0_r8*cff4*dtdz(i ,j,k2)- &
581 & dtdx(i+1,j,k1)))
582#endif
583!
584 cff1=min(dzde(i,j ,k1),0.0_r8)
585 cff2=min(dzde(i,j+1,k2),0.0_r8)
586 cff3=max(dzde(i,j ,k2),0.0_r8)
587 cff4=max(dzde(i,j+1,k1),0.0_r8)
588 tl_cff1=(0.5_r8+sign(0.5_r8,-dzde(i,j ,k1)))* &
589 & tl_dzde(i,j ,k1)
590 tl_cff2=(0.5_r8+sign(0.5_r8,-dzde(i,j+1,k2)))* &
591 & tl_dzde(i,j+1,k2)
592 tl_cff3=(0.5_r8+sign(0.5_r8, dzde(i,j ,k2)))* &
593 & tl_dzde(i,j ,k2)
594 tl_cff4=(0.5_r8+sign(0.5_r8, dzde(i,j+1,k1)))* &
595 & tl_dzde(i,j+1,k1)
596!^ FS(i,j,k2)=FS(i,j,k2)+ &
597!^ & cff* &
598!^ & (cff1*(cff1*dTdz(i,j,k2)-dTde(i,j ,k1))+ &
599!^ & cff2*(cff2*dTdz(i,j,k2)-dTde(i,j+1,k2))+ &
600!^ & cff3*(cff3*dTdz(i,j,k2)-dTde(i,j ,k2))+ &
601!^ & cff4*(cff4*dTdz(i,j,k2)-dTde(i,j+1,k1)))
602!^
603 tl_fs(i,j,k2)=tl_fs(i,j,k2)+ &
604 & cff* &
605 & (tl_cff1*(cff1*dtdz(i,j,k2)- &
606 & dtde(i,j ,k1))+ &
607 & tl_cff2*(cff2*dtdz(i,j,k2)- &
608 & dtde(i,j+1,k2))+ &
609 & tl_cff3*(cff3*dtdz(i,j,k2)- &
610 & dtde(i,j ,k2))+ &
611 & tl_cff4*(cff4*dtdz(i,j,k2)- &
612 & dtde(i,j+1,k1))+ &
613 & cff1*(tl_cff1*dtdz(i,j,k2)+ &
614 & cff1*tl_dtdz(i,j,k2)- &
615 & tl_dtde(i,j ,k1))+ &
616 & cff2*(tl_cff2*dtdz(i,j,k2)+ &
617 & cff2*tl_dtdz(i,j,k2)- &
618 & tl_dtde(i,j+1,k2))+ &
619 & cff3*(tl_cff3*dtdz(i,j,k2)+ &
620 & cff3*tl_dtdz(i,j,k2)- &
621 & tl_dtde(i,j ,k2))+ &
622 & cff4*(tl_cff4*dtdz(i,j,k2)+ &
623 & cff4*tl_dtdz(i,j,k2)- &
624 & tl_dtde(i,j+1,k1)))- &
625#ifdef TL_IOMS
626 & cff* &
627 & (cff1*(2.0_r8*cff1*dtdz(i,j,k2)- &
628 & dtde(i,j,k1))+ &
629 & cff2*(2.0_r8*cff2*dtdz(i,j ,k2)- &
630 & dtde(i,j+1,k2))+ &
631 & cff3*(2.0_r8*cff3*dtdz(i,j,k2)- &
632 & dtde(i,j,k2))+ &
633 & cff4*(2.0_r8*cff4*dtdz(i,j ,k2)- &
634 & dtde(i,j+1,k1)))
635#endif
636 END DO
637 END DO
638 END IF
639!
640! Time-step harmonic, geopotential diffusion term (m Tunits).
641!
642 DO j=jstr,jend
643 DO i=istr,iend
644!^ cff=dt(ng)*pm(i,j)*pn(i,j)* &
645!^ & (FX(i+1,j)-FX(i,j)+ &
646!^ & FE(i,j+1)-FE(i,j))+ &
647!^ & dt(ng)*(FS(i,j,k2)-FS(i,j,k1))
648!^
649 tl_cff=dt(ng)*pm(i,j)*pn(i,j)* &
650 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
651 & tl_fe(i,j+1)-tl_fe(i,j))+ &
652 & dt(ng)*(tl_fs(i,j,k2)-tl_fs(i,j,k1))
653!^ t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff
654!^
655 tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+tl_cff
656#ifdef DIAGNOSTICS_TS
657!! DiaTwrk(i,j,k,itrc,iThdif)=cff
658#endif
659 END DO
660 END DO
661 END IF
662 END DO k_loop
663 END DO t_loop
664!
665 RETURN
666 END SUBROUTINE rp_t3dmix2_geo_tile
667
668 END MODULE rp_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 irpm
Definition mod_param.F:664
real(dp), dimension(:), allocatable dt
logical, dimension(:,:), allocatable ltracerclm
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
integer, dimension(:), allocatable nstp
subroutine rp_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 rp_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