ROMS
Loading...
Searching...
No Matches
rp_t3dmix2_iso.h
Go to the documentation of this file.
1 MODULE rp_t3dmix2_mod
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 isopycnic surfaces. !
12! !
13! BASIC STATE variables needed: diff2, Hz, rho, 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, 26, __line__, myfile)
53#endif
54 CALL rp_t3dmix2_iso_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 & ocean(ng) % pden, &
80 & ocean(ng) % tl_pden, &
81#ifdef TS_MIX_CLIMA
82 & clima(ng) % tclm, &
83#endif
84#ifdef DIAGNOSTICS_TS
85!! & DIAGS(ng) % DiaTwrk, &
86#endif
87 & ocean(ng) % t, &
88 & ocean(ng) % tl_t)
89#ifdef PROFILE
90 CALL wclock_off (ng, irpm, 26, __line__, myfile)
91#endif
92!
93 RETURN
94 END SUBROUTINE rp_t3dmix2
95!
96!***********************************************************************
97 SUBROUTINE rp_t3dmix2_iso_tile (ng, tile, &
98 & LBi, UBi, LBj, UBj, &
99 & IminS, ImaxS, JminS, JmaxS, &
100 & nrhs, nstp, nnew, &
101#ifdef MASKING
102 & umask, vmask, &
103#endif
104#ifdef WET_DRY_NOT_YET
105 & umask_wet, vmask_wet, &
106#endif
107 & om_v, on_u, pm, pn, &
108 & Hz, tl_Hz, &
109 & z_r, tl_z_r, &
110#ifdef DIFF_3DCOEF
111 & diff3d_r, &
112#else
113 & diff2, &
114#endif
115 & pden, tl_pden, &
116#ifdef TS_MIX_CLIMA
117 & tclm, &
118#endif
119#ifdef DIAGNOSTICS_TS
120!! & DiaTwrk, &
121#endif
122 & t, tl_t)
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) :: rho(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_rho(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# ifdef TL_IOMS
374 & cff2
375# endif
376 cff3=max(pden(i,j,k)-pden(i,j,k+1),small)
377 tl_cff3=(0.5_r8+sign(0.5_r8,pden(i,j,k)-pden(i,j,k+1)- &
378 & small))* &
379 & (tl_pden(i,j,k)-tl_pden(i,j,k+1))+ &
380# ifdef TL_IOMS
381 & (0.5_r8-sign(0.5_r8, &
382 & pden(i,j,k)-pden(i,j,k+1)-small))* &
383 & small
384# endif
385 cff4=max(cff2,cff3)
386 tl_cff4=(0.5_r8+sign(0.5_r8,cff2-cff3))*tl_cff2+ &
387 & (0.5_r8-sign(0.5_r8,cff2-cff3))*tl_cff3
388 cff=-1.0_r8/cff4
389 tl_cff=cff*cff*tl_cff4+ &
390# ifdef TL_IOMS
391 & 2.0_r8*cff
392# endif
393#elif defined TS_MIX_MIN_STRAT
394 cff1=max(pden(i,j,k)-pden(i,j,k+1), &
395 & strat_min*(z_r(i,j,k+1)-z_r(i,j,k)))
396 tl_cff1=(0.5_r8+sign(0.5_r8, &
397 & pden(i,j,k)-pden(i,j,k+1)- &
398 & strat_min*(z_r(i,j,k+1)- &
399 & z_r(i,j,k ))))* &
400 & (tl_pden(i,j,k)-tl_pden(i,j,k+1))+ &
401 & (0.5_r8-sign(0.5_r8, &
402 & pden(i,j,k)-pden(i,j,k+1)- &
403 & strat_min*(z_r(i,j,k+1)- &
404 & z_r(i,j,k ))))* &
405 & (strat_min*(tl_z_r(i,j,k+1)-tl_z_r(i,j,k )))
406 cff=-1.0_r8/cff1
407 tl_cff=cff*cff*tl_cff1+ &
408# ifdef TL_IOMS
409 & 2.0_r8*cff
410# endif
411#else
412 cff1=max(pden(i,j,k)-pden(i,j,k+1),eps)
413 tl_cff1=(0.5_r8+sign(0.5_r8, &
414 & pden(i,j,k)-pden(i,j,k+1)-eps))* &
415 & (tl_pden(i,j,k)-tl_pden(i,j,k+1))+ &
416# ifdef TL_IOMS
417 & (0.5_r8- &
418 & sign(0.5_r8,pden(i,j,k)-pden(i,j,k+1)-eps))*eps
419# endif
420 cff=-1.0_r8/cff1
421 tl_cff=cff*cff*tl_cff1+ &
422# ifdef TL_IOMS
423 & 2.0_r8*cff
424# endif
425#endif
426#if defined TS_MIX_STABILITY
427 dtdr(i,j,k2)=cff*(0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
428 & t(i,j,k ,nrhs,itrc))+ &
429 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
430 & t(i,j,k ,nstp,itrc)))
431 tl_dtdr(i,j,k2)=tl_cff* &
432 & (0.75_r8*(t(i,j,k+1,nrhs,itrc)- &
433 & t(i,j,k ,nrhs,itrc))+ &
434 & 0.25_r8*(t(i,j,k+1,nstp,itrc)- &
435 & t(i,j,k ,nstp,itrc)))+ &
436 & cff* &
437 & (0.75_r8*(tl_t(i,j,k+1,nrhs,itrc)- &
438 & tl_t(i,j,k ,nrhs,itrc))+ &
439 & 0.25_r8*(tl_t(i,j,k+1,nstp,itrc)- &
440 & tl_t(i,j,k ,nstp,itrc)))- &
441# ifdef TL_IOMS
442 & dtdr(i,j,k2)
443# endif
444#elif defined TS_MIX_CLIMA
445 IF (ltracerclm(itrc,ng)) THEN
446 dtdr(i,j,k2)=cff*((t(i,j,k+1,nrhs,itrc)- &
447 & tclm(i,j,k+1,itrc))- &
448 & (t(i,j,k ,nrhs,itrc)- &
449 & tclm(i,j,k ,itrc)))
450 tl_dtdr(i,j,k2)=tl_cff*((t(i,j,k+1,nrhs,itrc)- &
451 & tclm(i,j,k+1,itrc))- &
452 & (t(i,j,k ,nrhs,itrc)- &
453 & tclm(i,j,k ,itrc)))+ &
454 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
455 & tl_t(i,j,k ,nrhs,itrc))- &
456# ifdef TL_IOMS
457 & dtdr(i,j,k2)
458# endif
459 ELSE
460 dtdr(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
461 & t(i,j,k ,nrhs,itrc))
462 tl_dtdr(i,j,k2)=tl_cff*(t(i,j,k+1,nrhs,itrc)- &
463 & t(i,j,k ,nrhs,itrc))+ &
464 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
465 & tl_t(i,j,k ,nrhs,itrc))- &
466# ifdef TL_IOMS
467 & dtdr(i,j,k2)
468# endif
469 END IF
470#else
471 dtdr(i,j,k2)=cff*(t(i,j,k+1,nrhs,itrc)- &
472 & t(i,j,k ,nrhs,itrc))
473 tl_dtdr(i,j,k2)=tl_cff*(t(i,j,k+1,nrhs,itrc)- &
474 & t(i,j,k ,nrhs,itrc))+ &
475 & cff*(tl_t(i,j,k+1,nrhs,itrc)- &
476 & tl_t(i,j,k ,nrhs,itrc))- &
477# ifdef TL_IOMS
478 & dtdr(i,j,k2)
479# endif
480#endif
481 fs(i,j,k2)=cff*(z_r(i,j,k+1)-z_r(i,j,k))
482 tl_fs(i,j,k2)=tl_cff*(z_r(i,j,k+1)-z_r(i,j,k))+ &
483 & cff*(tl_z_r(i,j,k+1)-tl_z_r(i,j,k))- &
484#ifdef TL_IOMS
485 & fs(i,j,k2)
486#endif
487 END DO
488 END DO
489 END IF
490!
491! Compute components of the rotated tracer flux (T m4/s) along
492! isopycnic surfaces.
493!
494 IF (k.gt.0) THEN
495 DO j=jstr,jend
496 DO i=istr,iend+1
497#ifdef DIFF_3DCOEF
498 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
499 & on_u(i,j)
500#else
501 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i-1,j,itrc))* &
502 & on_u(i,j)
503#endif
504!^ FX(i,j)=cff* &
505!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
506!^ & (dTdx(i,j,k1)- &
507!^ & 0.5_r8*(MAX(dRdx(i,j,k1),0.0_r8)* &
508!^ & (dTdr(i-1,j,k1)+ &
509!^ & dTdr(i ,j,k2))+ &
510!^ & MIN(dRdx(i,j,k1),0.0_r8)* &
511!^ & (dTdr(i-1,j,k2)+ &
512!^ & dTdr(i,j,k1))))
513!^
514 tl_fx(i,j)=cff* &
515 & (((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
516 & (dtdx(i,j,k1)- &
517 & 0.5_r8*(max(drdx(i,j,k1),0.0_r8)* &
518 & (dtdr(i-1,j,k1)+ &
519 & dtdr(i ,j,k2))+ &
520 & min(drdx(i,j,k1),0.0_r8)* &
521 & (dtdr(i-1,j,k2)+ &
522 & dtdr(i ,j,k1)))))+ &
523 & ((hz(i,j,k)+hz(i-1,j,k))* &
524 & (tl_dtdx(i,j,k1)- &
525 & 0.5_r8*(max(drdx(i,j,k1),0.0_r8)* &
526 & (tl_dtdr(i-1,j,k1)+ &
527 & tl_dtdr(i ,j,k2))+ &
528 & min(drdx(i,j,k1),0.0_r8)* &
529 & (tl_dtdr(i-1,j,k2)+ &
530 & tl_dtdr(i ,j,k1)))- &
531 & 0.5_r8*((0.5_r8+ &
532 & sign(0.5_r8, drdx(i,j,k1)))* &
533 & tl_drdx(i,j,k1)* &
534 & (dtdr(i-1,j,k1)+dtdr(i,j,k2))+ &
535 & (0.5_r8+ &
536 & sign(0.5_r8,-drdx(i,j,k1)))* &
537 & tl_drdx(i,j,k1)* &
538 & (dtdr(i-1,j,k2)+dtdr(i,j,k1))))))-&
539#ifdef TL_IOMS
540 & cff* &
541 & (hz(i,j,k)+hz(i-1,j,k))* &
542 & (dtdx(i,j,k1)- &
543 & (max(drdx(i,j,k1),0.0_r8)* &
544 & (dtdr(i-1,j,k1)+ &
545 & dtdr(i ,j,k2))+ &
546 & min(drdx(i,j,k1),0.0_r8)* &
547 & (dtdr(i-1,j,k2)+ &
548 & dtdr(i,j,k1))))
549#endif
550 END DO
551 END DO
552 DO j=jstr,jend+1
553 DO i=istr,iend
554#ifdef DIFF_3DCOEF
555 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
556 & om_v(i,j)
557#else
558 cff=0.25_r8*(diff2(i,j,itrc)+diff2(i,j-1,itrc))* &
559 & om_v(i,j)
560#endif
561!^ FE(i,j)=cff* &
562!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
563!^ & (dTde(i,j,k1)- &
564!^ & 0.5_r8*(MAX(dRde(i,j,k1),0.0_r8)* &
565!^ & (dTdr(i,j-1,k1)+ &
566!^ & dTdr(i,j ,k2))+ &
567!^ & MIN(dRde(i,j,k1),0.0_r8)* &
568!^ & (dTdr(i,j-1,k2)+ &
569!^ & dTdr(i,j ,k1))))
570!^
571 tl_fe(i,j)=cff* &
572 & (((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
573 & (dtde(i,j,k1)- &
574 & 0.5_r8*(max(drde(i,j,k1),0.0_r8)* &
575 & (dtdr(i,j-1,k1)+ &
576 & dtdr(i,j ,k2))+ &
577 & min(drde(i,j,k1),0.0_r8)* &
578 & (dtdr(i,j-1,k2)+ &
579 & dtdr(i,j ,k1)))))+ &
580 & ((hz(i,j,k)+hz(i,j-1,k))* &
581 & (tl_dtde(i,j,k1)- &
582 & 0.5_r8*(max(drde(i,j,k1),0.0_r8)* &
583 & (tl_dtdr(i,j-1,k1)+ &
584 & tl_dtdr(i,j ,k2))+ &
585 & min(drde(i,j,k1),0.0_r8)* &
586 & (tl_dtdr(i,j-1,k2)+ &
587 & tl_dtdr(i,j ,k1)))- &
588 & 0.5_r8*((0.5_r8+ &
589 & sign(0.5_r8, drde(i,j,k1)))* &
590 & tl_drde(i,j,k1)* &
591 & (dtdr(i,j-1,k1)+dtdr(i,j,k2))+ &
592 & (0.5_r8+ &
593 & sign(0.5_r8,-drde(i,j,k1)))* &
594 & tl_drde(i,j,k1)* &
595 & (dtdr(i,j-1,k2)+dtdr(i,j,k1))))))-&
596#ifdef TL_IOMS
597 & cff* &
598 & (hz(i,j,k)+hz(i,j-1,k))* &
599 & (dtde(i,j,k1)- &
600 & (max(drde(i,j,k1),0.0_r8)* &
601 & (dtdr(i,j-1,k1)+ &
602 & dtdr(i,j ,k2))+ &
603 & min(drde(i,j,k1),0.0_r8)* &
604 & (dtdr(i,j-1,k2)+ &
605 & dtdr(i,j ,k1))))
606#endif
607 END DO
608 END DO
609 IF (k.lt.n(ng)) THEN
610 DO j=jstr,jend
611 DO i=istr,iend
612 cff1=max(drdx(i ,j,k1),0.0_r8)
613 cff2=max(drdx(i+1,j,k2),0.0_r8)
614 cff3=min(drdx(i ,j,k2),0.0_r8)
615 cff4=min(drdx(i+1,j,k1),0.0_r8)
616 tl_cff1=(0.5_r8+sign(0.5_r8, drdx(i ,j,k1)))* &
617 & tl_drdx(i ,j,k1)
618 tl_cff2=(0.5_r8+sign(0.5_r8, drdx(i+1,j,k2)))* &
619 & tl_drdx(i+1,j,k2)
620 tl_cff3=(0.5_r8+sign(0.5_r8,-drdx(i ,j,k2)))* &
621 & tl_drdx(i ,j,k2)
622 tl_cff4=(0.5_r8+sign(0.5_r8,-drdx(i+1,j,k1)))* &
623 & tl_drdx(i+1,j,k1)
624 cff=cff1*(cff1*dtdr(i,j,k2)-dtdx(i ,j,k1))+ &
625 & cff2*(cff2*dtdr(i,j,k2)-dtdx(i+1,j,k2))+ &
626 & cff3*(cff3*dtdr(i,j,k2)-dtdx(i ,j,k2))+ &
627 & cff4*(cff4*dtdr(i,j,k2)-dtdx(i+1,j,k1))
628 tl_cff=tl_cff1*(cff1*dtdr(i ,j,k2)- &
629 & dtdx(i ,j,k1))+ &
630 & tl_cff2*(cff2*dtdr(i,j,k2)- &
631 & dtdx(i+1,j,k2))+ &
632 & tl_cff3*(cff3*dtdr(i,j,k2)- &
633 & dtdx(i ,j,k2))+ &
634 & tl_cff4*(cff4*dtdr(i,j,k2)- &
635 & dtdx(i+1,j,k1))+ &
636 & cff1*(tl_cff1*dtdr(i,j,k2)+ &
637 & cff1*tl_dtdr(i,j,k2)- &
638 & tl_dtdx(i ,j,k1))+ &
639 & cff2*(tl_cff2*dtdr(i,j,k2)+ &
640 & cff2*tl_dtdr(i,j,k2)- &
641 & tl_dtdx(i+1,j,k2))+ &
642 & cff3*(tl_cff3*dtdr(i,j,k2)+ &
643 & cff3*tl_dtdr(i,j,k2)- &
644 & tl_dtdx(i ,j,k2))+ &
645 & cff4*(tl_cff4*dtdr(i,j,k2)+ &
646 & cff4*tl_dtdr(i,j,k2)- &
647 & tl_dtdx(i+1,j,k1))- &
648#ifdef TL_IOMS
649 & cff1*(2.0_r8*cff1*dtdr(i,j,k2)- &
650 & dtdx(i,j,k1))- &
651 & cff2*(2.0_r8*cff2*dtdr(i ,j,k2)- &
652 & dtdx(i+1,j,k2))- &
653 & cff3*(2.0_r8*cff3*dtdr(i,j,k2)- &
654 & dtdx(i,j,k2))- &
655 & cff4*(2.0_r8*cff4*dtdr(i ,j,k2)- &
656 & dtdx(i+1,j,k1))
657#endif
658 cff1=max(drde(i,j ,k1),0.0_r8)
659 cff2=max(drde(i,j+1,k2),0.0_r8)
660 cff3=min(drde(i,j ,k2),0.0_r8)
661 cff4=min(drde(i,j+1,k1),0.0_r8)
662 tl_cff1=(0.5_r8+sign(0.5_r8, drde(i,j ,k1)))* &
663 & tl_drde(i,j ,k1)
664 tl_cff2=(0.5_r8+sign(0.5_r8, drde(i,j+1,k2)))* &
665 & tl_drde(i,j+1,k2)
666 tl_cff3=(0.5_r8+sign(0.5_r8,-drde(i,j ,k2)))* &
667 & tl_drde(i,j ,k2)
668 tl_cff4=(0.5_r8+sign(0.5_r8,-drde(i,j+1,k1)))* &
669 & tl_drde(i,j+1,k1)
670 cff=cff+ &
671 & cff1*(cff1*dtdr(i,j,k2)-dtde(i,j ,k1))+ &
672 & cff2*(cff2*dtdr(i,j,k2)-dtde(i,j+1,k2))+ &
673 & cff3*(cff3*dtdr(i,j,k2)-dtde(i,j ,k2))+ &
674 & cff4*(cff4*dtdr(i,j,k2)-dtde(i,j+1,k1))
675 tl_cff=tl_cff+ &
676 & tl_cff1*(cff1*dtdr(i,j,k2)- &
677 & dtde(i,j ,k1))+ &
678 & tl_cff2*(cff2*dtdr(i,j,k2)- &
679 & dtde(i,j+1,k2))+ &
680 & tl_cff3*(cff3*dtdr(i,j,k2)- &
681 & dtde(i,j ,k2))+ &
682 & tl_cff4*(cff4*dtdr(i,j,k2)- &
683 & dtde(i,j+1,k1))+ &
684 & cff1*(tl_cff1*dtdr(i,j,k2)+ &
685 & cff1*tl_dtdr(i,j,k2)- &
686 & tl_dtde(i,j ,k1))+ &
687 & cff2*(tl_cff2*dtdr(i,j,k2)+ &
688 & cff2*tl_dtdr(i,j,k2)- &
689 & tl_dtde(i,j+1,k2))+ &
690 & cff3*(tl_cff3*dtdr(i,j,k2)+ &
691 & cff3*tl_dtdr(i,j,k2)- &
692 & tl_dtde(i,j ,k2))+ &
693 & cff4*(tl_cff4*dtdr(i,j,k2)+ &
694 & cff4*tl_dtdr(i,j,k2)- &
695 & tl_dtde(i,j+1,k1))- &
696#ifdef TL_IOMS
697 & cff1*(2.0_r8*cff1*dtdr(i,j,k2)- &
698 & dtde(i,j,k1))- &
699 & cff2*(2.0_r8*cff2*dtdr(i,j ,k2)- &
700 & dtde(i,j+1,k2))- &
701 & cff3*(2.0_r8*cff3*dtdr(i,j,k2)- &
702 & dtde(i,j,k2))- &
703 & cff4*(2.0_r8*cff4*dtdr(i,j ,k2)- &
704 & dtde(i,j+1,k1))
705#endif
706#ifdef DIFF_3DCOEF
707!^ FS(i,j,k2)=0.5_r8*cff*diff3d_r(i,j,k)*FS(i,j,k2)
708!^
709 tl_fs(i,j,k2)=0.5_r8*diff3d_r(i,j,k)* &
710 & (tl_cff*fs(i,j,k2)+ &
711 & cff*tl_fs(i,j,k2))- &
712# ifdef TL_IOMS
713 & 0.5_r8*diff3d_r(i,j,k)*cff*fs(i,j,k2)
714# endif
715#else
716!^ FS(i,j,k2)=0.5_r8*cff*diff2(i,j,itrc)*FS(i,j,k2)
717!^
718 tl_fs(i,j,k2)=0.5_r8*diff2(i,j,itrc)* &
719 & (tl_cff*fs(i,j,k2)+ &
720 & cff*tl_fs(i,j,k2))- &
721# ifdef TL_IOMS
722 & 0.5_r8*diff2(i,j,itrc)*cff*fs(i,j,k2)
723# endif
724#endif
725 END DO
726 END DO
727 END IF
728!
729! Time-step harmonic, isopycnic diffusion term (m Tunits).
730!
731 DO j=jstr,jend
732 DO i=istr,iend
733!^ cff=dt(ng)*pm(i,j)*pn(i,j)* &
734!^ & (FX(i+1,j)-FX(i,j)+ &
735!^ & FE(i,j+1)-FE(i,j))+ &
736!^ & dt(ng)*(FS(i,j,k2)-FS(i,j,k1))
737!^
738 tl_cff=dt(ng)*pm(i,j)*pn(i,j)* &
739 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
740 & tl_fe(i,j+1)-tl_fe(i,j))+ &
741 & dt(ng)*(tl_fs(i,j,k2)-tl_fs(i,j,k1))
742!^ t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)+cff
743!^
744 tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+tl_cff
745#ifdef DIAGNOSTICS_TS
746!! DiaTwrk(i,j,k,itrc,iThdif)=cff
747#endif
748 END DO
749 END DO
750 END IF
751 END DO k_loop
752 END DO t_loop
753!
754 RETURN
755 END SUBROUTINE rp_t3dmix2_iso_tile
756
757 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_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, 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