ROMS
Loading...
Searching...
No Matches
rp_t3dmix4_s.h
Go to the documentation of this file.
1 MODULE rp_t3dmix4_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! biharmonic mixing of tracers along S-coordinate levels surfaces. !
12! !
13!=======================================================================
14!
15 implicit none
16!
17 PRIVATE
18 PUBLIC rp_t3dmix4
19!
20 CONTAINS
21!
22!***********************************************************************
23 SUBROUTINE rp_t3dmix4 (ng, tile)
24!***********************************************************************
25!
26 USE mod_param
27#ifdef TS_MIX_CLIMA
28 USE mod_clima
29#endif
30#ifdef DIAGNOSTICS_TS
31!! USE mod_diags
32#endif
33 USE mod_grid
34 USE mod_mixing
35 USE mod_ocean
36 USE mod_stepping
37!
38! Imported variable declarations.
39!
40 integer, intent(in) :: ng, tile
41!
42! Local variable declarations.
43!
44 character (len=*), parameter :: MyFile = &
45 & __FILE__
46!
47#include "tile.h"
48!
49#ifdef PROFILE
50 CALL wclock_on (ng, irpm, 27, __line__, myfile)
51#endif
52 CALL rp_t3dmix4_s_tile (ng, tile, &
53 & lbi, ubi, lbj, ubj, &
54 & imins, imaxs, jmins, jmaxs, &
55 & nrhs(ng), nstp(ng), nnew(ng), &
56#ifdef MASKING
57 & grid(ng) % umask, &
58 & grid(ng) % vmask, &
59#endif
60#ifdef WET_DRY_NOT_YET
61 & grid(ng) % umask_wet, &
62 & grid(ng) % vmask_wet, &
63#endif
64 & grid(ng) % Hz, &
65 & grid(ng) % tl_Hz, &
66 & grid(ng) % pmon_u, &
67 & grid(ng) % pnom_v, &
68 & grid(ng) % pm, &
69 & grid(ng) % pn, &
70#ifdef DIFF_3DCOEF
71# ifdef TS_U3ADV_SPLIT_NOT_YET
72 & mixing(ng) % diff3d_u, &
73 & mixing(ng) % diff3d_v, &
74# else
75 & mixing(ng) % diff3d_r, &
76# endif
77#else
78 & mixing(ng) % diff4, &
79#endif
80#ifdef TS_MIX_CLIMA
81 & clima(ng) % tclm, &
82#endif
83#ifdef DIAGNOSTICS_TS
84!! & DIAGS(ng) % DiaTwrk, &
85#endif
86 & ocean(ng) % t, &
87 & ocean(ng) % tl_t)
88#ifdef PROFILE
89 CALL wclock_off (ng, irpm, 27, __line__, myfile)
90#endif
91!
92 RETURN
93 END SUBROUTINE rp_t3dmix4
94!
95!***********************************************************************
96 SUBROUTINE rp_t3dmix4_s_tile (ng, tile, &
97 & LBi, UBi, LBj, UBj, &
98 & IminS, ImaxS, JminS, JmaxS, &
99 & nrhs, nstp, nnew, &
100#ifdef MASKING
101 & umask, vmask, &
102#endif
103#ifdef WET_DRY_NOT_YET
104 & umask_wet, vmask_wet, &
105#endif
106 & Hz, tl_Hz, &
107 & pmon_u, pnom_v, pm, pn, &
108#ifdef DIFF_3DCOEF
109# ifdef TS_U3ADV_SPLIT_NOT_YET
110 & diff3d_u, diff3d_v, &
111# else
112 & diff3d_r, &
113# endif
114#else
115 & diff4, &
116#endif
117#ifdef TS_MIX_CLIMA
118 & tclm, &
119#endif
120#ifdef DIAGNOSTICS_TS
121!! & DiaTwrk, &
122#endif
123 & t, tl_t)
124!***********************************************************************
125!
126 USE mod_param
127 USE mod_ncparam
128 USE mod_scalars
129!
130! Imported variable declarations.
131!
132 integer, intent(in) :: ng, tile
133 integer, intent(in) :: LBi, UBi, LBj, UBj
134 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
135 integer, intent(in) :: nrhs, nstp, nnew
136
137#ifdef ASSUMED_SHAPE
138# ifdef MASKING
139 real(r8), intent(in) :: umask(LBi:,LBj:)
140 real(r8), intent(in) :: vmask(LBi:,LBj:)
141# endif
142# ifdef WET_DRY_NOT_YET
143 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
144 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
145# endif
146# ifdef DIFF_3DCOEF
147# ifdef TS_U3ADV_SPLIT_NOT_YET
148 real(r8), intent(in) :: diff3d_u(LBi:,LBj:,:)
149 real(r8), intent(in) :: diff3d_v(LBi:,LBj:,:)
150# else
151 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
152# endif
153# else
154 real(r8), intent(in) :: diff4(LBi:,LBj:,:)
155# endif
156 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
157 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
158 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
159 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
160 real(r8), intent(in) :: pm(LBi:,LBj:)
161 real(r8), intent(in) :: pn(LBi:,LBj:)
162 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
163# ifdef TS_MIX_CLIMA
164 real(r8), intent(in) :: tclm(LBi:,LBj:,:,:)
165# endif
166# ifdef DIAGNOSTICS_TS
167!! real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
168# endif
169 real(r8), intent(inout) :: tl_t(LBi:,LBj:,:,:,:)
170#else
171# ifdef MASKING
172 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
173 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
174# endif
175# ifdef WET_DRY_NOT_YET
176 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
177 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
178# endif
179# ifdef DIFF_3DCOEF
180# ifdef TS_U3ADV_SPLIT_NOT_YET
181 real(r8), intent(in) :: diff3d_u(LBi:UBi,LBj:UBj,N(ng))
182 real(r8), intent(in) :: diff3d_v(LBi:UBi,LBj:UBj,N(ng))
183# else
184 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
185# endif
186# else
187 real(r8), intent(in) :: diff4(LBi:UBi,LBj:UBj,NT(ng))
188# endif
189 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
190 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
191 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
192 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
193 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
194 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
195 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
196# ifdef TS_MIX_CLIMA
197 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
198# endif
199# ifdef DIAGNOSTICS_TS
200!! real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
201!! & NDT)
202# endif
203 real(r8), intent(inout) :: tl_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
204#endif
205!
206! Local variable declarations.
207!
208 integer :: Imin, Imax, Jmin, Jmax
209 integer :: i, itrc, j, k
210
211 real(r8) :: cff, cff1, tl_cff, tl_cff1
212
213 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE
214 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX
215 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: LapT
216
217 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FE
218 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FX
219 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_LapT
220
221#include "set_bounds.h"
222!
223!-----------------------------------------------------------------------
224! Compute horizontal biharmonic diffusion along constant S-surfaces.
225! The biharmonic operator is computed by applying the harmonic
226! operator twice.
227#ifdef TS_MIX_STABILITY
228! In order to increase stability, the biharmonic operator is applied
229! as: 3/4 t(:,:,:,nrhs,:) + 1/4 t(:,:,:,nstp,:).
230#endif
231!-----------------------------------------------------------------------
232!
233! Set local I- and J-ranges.
234!
235 IF (ewperiodic(ng)) THEN
236 imin=istr-1
237 imax=iend+1
238 ELSE
239 imin=max(istr-1,1)
240 imax=min(iend+1,lm(ng))
241 END IF
242 IF (nsperiodic(ng)) THEN
243 jmin=jstr-1
244 jmax=jend+1
245 ELSE
246 jmin=max(jstr-1,1)
247 jmax=min(jend+1,mm(ng))
248 END IF
249!
250! Compute horizontal tracer flux in the XI- and ETA-directions.
251!
252 DO itrc=1,nt(ng)
253 DO k=1,n(ng)
254 DO j=jmin,jmax
255 DO i=imin,imax+1
256#ifdef DIFF_3DCOEF
257# ifdef TS_U3ADV_SPLIT_NOT_YET
258 cff=0.5_r8*diff3d_u(i,j,k)*pmon_u(i,j)
259# else
260 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
261 & pmon_u(i,j)
262# endif
263#else
264 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i-1,j,itrc))* &
265 & pmon_u(i,j)
266#endif
267#ifdef MASKING
268 cff=cff*umask(i,j)
269#endif
270#ifdef WET_DRY_NOT_YET
271 cff=cff*umask_wet(i,j)
272#endif
273#if defined TS_MIX_STABILITY
274 fx(i,j)=cff* &
275 & (hz(i,j,k)+hz(i-1,j,k))* &
276 & (0.75_r8*(t(i ,j,k,nrhs,itrc)- &
277 & t(i-1,j,k,nrhs,itrc))+ &
278 & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
279 & t(i-1,j,k,nstp,itrc)))
280 tl_fx(i,j)=cff* &
281 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
282 & (0.75_r8*(t(i ,j,k,nrhs,itrc)- &
283 & t(i-1,j,k,nrhs,itrc))+ &
284 & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
285 & t(i-1,j,k,nstp,itrc)))+ &
286 & (hz(i,j,k)+hz(i-1,j,k))* &
287 & (0.75_r8*(tl_t(i ,j,k,nrhs,itrc)- &
288 & tl_t(i-1,j,k,nrhs,itrc))+ &
289 & 0.25_r8*(tl_t(i ,j,k,nstp,itrc)- &
290 & tl_t(i-1,j,k,nstp,itrc))))- &
291# ifdef TL_IOMS
292 & fx(i,j)
293# endif
294#elif defined TS_MIX_CLIMA
295 IF (ltracerclm(itrc,ng)) THEN
296 fx(i,j)=cff* &
297 & (hz(i,j,k)+hz(i-1,j,k))* &
298 & ((t(i ,j,k,nrhs,itrc)-tclm(i ,j,k,itrc))- &
299 & (t(i-1,j,k,nrhs,itrc)-tclm(i-1,j,k,itrc)))
300 tl_fx(i,j)=cff* &
301 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
302 & ((t(i ,j,k,nrhs,itrc)- &
303 & tclm(i ,j,k,itrc))- &
304 & (t(i-1,j,k,nrhs,itrc)- &
305 & tclm(i-1,j,k,itrc)))+ &
306 & (hz(i,j,k)+hz(i-1,j,k))* &
307 & (tl_t(i ,j,k,nrhs,itrc)- &
308 & tl_t(i-1,j,k,nrhs,itrc)))- &
309# ifdef TL_IOMS
310 & fx(i,j)
311# endif
312 ELSE
313 fx(i,j)=cff* &
314 & (hz(i,j,k)+hz(i-1,j,k))* &
315 & (t(i ,j,k,nrhs,itrc)- &
316 & t(i-1,j,k,nrhs,itrc))
317 tl_fx(i,j)=cff* &
318 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
319 & (t(i ,j,k,nrhs,itrc)- &
320 & t(i-1,j,k,nrhs,itrc))+ &
321 & (hz(i,j,k)+hz(i-1,j,k))* &
322 & (tl_t(i ,j,k,nrhs,itrc)- &
323 & tl_t(i-1,j,k,nrhs,itrc)))- &
324# ifdef TL_IOMS
325 & fx(i,j)
326# endif
327 END IF
328#else
329 fx(i,j)=cff* &
330 & (hz(i,j,k)+hz(i-1,j,k))* &
331 & (t(i ,j,k,nrhs,itrc)- &
332 & t(i-1,j,k,nrhs,itrc))
333 tl_fx(i,j)=cff* &
334 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
335 & (t(i ,j,k,nrhs,itrc)- &
336 & t(i-1,j,k,nrhs,itrc))+ &
337 & (hz(i,j,k)+hz(i-1,j,k))* &
338 & (tl_t(i ,j,k,nrhs,itrc)- &
339 & tl_t(i-1,j,k,nrhs,itrc)))- &
340# ifdef TL_IOMS
341 & fx(i,j)
342# endif
343#endif
344 END DO
345 END DO
346 DO j=jmin,jmax+1
347 DO i=imin,imax
348#ifdef DIFF_3DCOEF
349# ifdef TS_U3ADV_SPLIT_NOT_YET
350 cff=0.5_r8*diff3d_v(i,j,k)*pnom_v(i,j)
351# else
352 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
353 & pnom_v(i,j)
354# endif
355#else
356 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i,j-1,itrc))* &
357 & pnom_v(i,j)
358#endif
359#ifdef MASKING
360 cff=cff*vmask(i,j)
361#endif
362#ifdef WET_DRY_NOT_YET
363 cff=cff*vmask_wet(i,j)
364#endif
365#if defined TS_MIX_STABILITY
366 fe(i,j)=cff* &
367 & (hz(i,j,k)+hz(i,j-1,k))* &
368 & (0.75_r8*(t(i,j ,k,nrhs,itrc)- &
369 & t(i,j-1,k,nrhs,itrc))+ &
370 & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
371 & t(i,j-1,k,nstp,itrc)))
372 tl_fe(i,j)=cff* &
373 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
374 & (0.75_r8*(t(i,j ,k,nrhs,itrc)- &
375 & t(i,j-1,k,nrhs,itrc))+ &
376 & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
377 & t(i,j-1,k,nstp,itrc)))+ &
378 & (hz(i,j,k)+hz(i,j-1,k))* &
379 & (0.75_r8*(tl_t(i,j ,k,nrhs,itrc)- &
380 & tl_t(i,j-1,k,nrhs,itrc))+ &
381 & 0.25_r8*(tl_t(i,j ,k,nstp,itrc)- &
382 & tl_t(i,j-1,k,nstp,itrc))))- &
383# ifdef TL_IOMS
384 & fe(i,j)
385# endif
386#elif defined TS_MIX_CLIMA
387 IF (ltracerclm(itrc,ng)) THEN
388 fe(i,j)=cff* &
389 & (hz(i,j,k)+hz(i,j-1,k))* &
390 & ((t(i,j ,k,nrhs,itrc)-tclm(i,j ,k,itrc))- &
391 & (t(i,j-1,k,nrhs,itrc)-tclm(i,j-1,k,itrc)))
392 tl_fe(i,j)=cff* &
393 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
394 & ((t(i,j ,k,nrhs,itrc)- &
395 & tclm(i,j ,k,itrc))- &
396 & (t(i,j-1,k,nrhs,itrc)- &
397 & tclm(i,j-1,k,itrc)))+ &
398 & (hz(i,j,k)+hz(i,j-1,k))* &
399 & (tl_t(i,j ,k,nrhs,itrc)- &
400 & tl_t(i,j-1,k,nrhs,itrc)))- &
401# ifdef TL_IOMS
402 & fe(i,j)
403# endif
404 ELSE
405 fe(i,j)=cff* &
406 & (hz(i,j,k)+hz(i,j-1,k))* &
407 & (t(i,j ,k,nrhs,itrc)- &
408 & t(i,j-1,k,nrhs,itrc))
409 tl_fe(i,j)=cff* &
410 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
411 & (t(i,j ,k,nrhs,itrc)- &
412 & t(i,j-1,k,nrhs,itrc))+ &
413 & (hz(i,j,k)+hz(i,j-1,k))* &
414 & (tl_t(i,j ,k,nrhs,itrc)- &
415 & tl_t(i,j-1,k,nrhs,itrc)))- &
416# ifdef TL_IOMS
417 & fe(i,j)
418# endif
419 END IF
420#else
421 fe(i,j)=cff* &
422 & (hz(i,j,k)+hz(i,j-1,k))* &
423 & (t(i,j ,k,nrhs,itrc)- &
424 & t(i,j-1,k,nrhs,itrc))
425 tl_fe(i,j)=cff* &
426 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
427 & (t(i,j ,k,nrhs,itrc)- &
428 & t(i,j-1,k,nrhs,itrc))+ &
429 & (hz(i,j,k)+hz(i,j-1,k))* &
430 & (tl_t(i,j ,k,nrhs,itrc)- &
431 & tl_t(i,j-1,k,nrhs,itrc)))- &
432# ifdef TL_IOMS
433 & fe(i,j)
434# endif
435#endif
436 END DO
437 END DO
438!
439! Compute first harmonic operator and multiply by the metrics of the
440! second harmonic operator.
441!
442 DO j=jmin,jmax
443 DO i=imin,imax
444 cff=1.0_r8/hz(i,j,k)
445 tl_cff=-cff*cff*tl_hz(i,j,k)+ &
446#ifdef TL_IOMS
447 & 2.0_r8*cff
448#endif
449 lapt(i,j)=pm(i,j)*pn(i,j)*cff* &
450 & (fx(i+1,j)-fx(i,j)+ &
451 & fe(i,j+1)-fe(i,j))
452 tl_lapt(i,j)=pm(i,j)*pn(i,j)* &
453 & (tl_cff* &
454 & (fx(i+1,j)-fx(i,j)+ &
455 & fe(i,j+1)-fe(i,j))+ &
456 & cff* &
457 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
458 & tl_fe(i,j+1)-tl_fe(i,j)))- &
459#ifdef TL_IOMS
460 & lapt(i,j)
461#endif
462 END DO
463 END DO
464!
465! Apply boundary conditions (except periodic; closed or gradient)
466! to the first harmonic operator.
467!
468 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
469 IF (domain(ng)%Western_Edge(tile)) THEN
470 IF (tl_lbc(iwest,istvar(itrc),ng)%closed) THEN
471 DO j=jmin,jmax
472 lapt(istr-1,j)=0.0_r8
473 tl_lapt(istr-1,j)=0.0_r8
474 END DO
475 ELSE
476 DO j=jmin,jmax
477 lapt(istr-1,j)=lapt(istr,j)
478 tl_lapt(istr-1,j)=tl_lapt(istr,j)
479 END DO
480 END IF
481 END IF
482 END IF
483!
484 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
485 IF (domain(ng)%Eastern_Edge(tile)) THEN
486 IF (tl_lbc(ieast,istvar(itrc),ng)%closed) THEN
487 DO j=jmin,jmax
488 lapt(iend+1,j)=0.0_r8
489 tl_lapt(iend+1,j)=0.0_r8
490 END DO
491 ELSE
492 DO j=jmin,jmax
493 lapt(iend+1,j)=lapt(iend,j)
494 tl_lapt(iend+1,j)=tl_lapt(iend,j)
495 END DO
496 END IF
497 END IF
498 END IF
499!
500 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
501 IF (domain(ng)%Southern_Edge(tile)) THEN
502 IF (tl_lbc(isouth,istvar(itrc),ng)%closed) THEN
503 DO i=imin,imax
504 lapt(i,jstr-1)=0.0_r8
505 tl_lapt(i,jstr-1)=0.0_r8
506 END DO
507 ELSE
508 DO i=imin,imax
509 lapt(i,jstr-1)=lapt(i,jstr)
510 tl_lapt(i,jstr-1)=tl_lapt(i,jstr)
511 END DO
512 END IF
513 END IF
514 END IF
515!
516 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
517 IF (domain(ng)%Northern_Edge(tile)) THEN
518 IF (tl_lbc(inorth,istvar(itrc),ng)%closed) THEN
519 DO i=imin,imax
520 lapt(i,jend+1)=0.0_r8
521 tl_lapt(i,jend+1)=0.0_r8
522 END DO
523 ELSE
524 DO i=imin,imax
525 lapt(i,jend+1)=lapt(i,jend)
526 tl_lapt(i,jend+1)=tl_lapt(i,jend)
527 END DO
528 END IF
529 END IF
530 END IF
531!
532! Compute FX=d(LapT)/d(xi) and FE=d(LapT)/d(eta) terms.
533!
534 DO j=jstr,jend
535 DO i=istr,iend+1
536#ifdef DIFF_3DCOEF
537# ifdef TS_U3ADV_SPLIT_NOT_YET
538 cff=0.5_r8*diff3d_u(i,j,k)*pmon_u(i,j)
539# else
540 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
541 & pmon_u(i,j)
542# endif
543#else
544 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i-1,j,itrc))* &
545 & pmon_u(i,j)
546#endif
547!^ FX(i,j)=cff* &
548!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
549!^ & (LapT(i,j)-LapT(i-1,j))
550!^
551 tl_fx(i,j)=cff* &
552 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
553 & (lapt(i,j)-lapt(i-1,j))+ &
554 & (hz(i,j,k)+hz(i-1,j,k))* &
555 & (tl_lapt(i,j)-tl_lapt(i-1,j)))- &
556# ifdef TL_IOMS
557 & cff* &
558 & (hz(i,j,k)+hz(i-1,j,k))* &
559 & (lapt(i,j)-lapt(i-1,j))
560# endif
561#ifdef MASKING
562!^ FX(i,j)=FX(i,j)*umask(i,j)
563!^
564 tl_fx(i,j)=tl_fx(i,j)*umask(i,j)
565#endif
566#ifdef WET_DRY_NOT_YET
567 fx(i,j)=fx(i,j)*umask_wet(i,j)
568#endif
569 END DO
570 END DO
571 DO j=jstr,jend+1
572 DO i=istr,iend
573#ifdef DIFF_3DCOEF
574# ifdef TS_U3ADV_SPLIT_NOT_YET
575 cff=0.5_r8*diff3d_v(i,j,k)*pnom_v(i,j)
576# else
577 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
578 & pnom_v(i,j)
579# endif
580#else
581 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i,j-1,itrc))* &
582 & pnom_v(i,j)
583#endif
584!^ FE(i,j)=cff* &
585!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
586!^ & (LapT(i,j)-LapT(i,j-1))
587!^
588 tl_fe(i,j)=cff* &
589 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
590 & (lapt(i,j)-lapt(i,j-1))+ &
591 & (hz(i,j,k)+hz(i,j-1,k))* &
592 & (tl_lapt(i,j)-tl_lapt(i,j-1)))- &
593#ifdef TL_IOMS
594 & cff* &
595 & (hz(i,j,k)+hz(i,j-1,k))* &
596 & (lapt(i,j)-lapt(i,j-1))
597#endif
598#ifdef MASKING
599!^ FE(i,j)=FE(i,j)*vmask(i,j)
600!^
601 tl_fe(i,j)=tl_fe(i,j)*vmask(i,j)
602#endif
603#ifdef WET_DRY_NOT_YET
604 fe(i,j)=fe(i,j)*vmask_wet(i,j)
605#endif
606 END DO
607 END DO
608!
609! Time-step biharmonic, S-surfaces diffusion term (m Tunits).
610!
611 DO j=jstr,jend
612 DO i=istr,iend
613!^ cff=dt(ng)*pm(i,j)*pn(i,j)* &
614!^ & (FX(i+1,j)-FX(i,j)+ &
615!^ & FE(i,j+1)-FE(i,j))
616!^
617 tl_cff=dt(ng)*pm(i,j)*pn(i,j)* &
618 & (tl_fx(i+1,j)-tl_fx(i,j)+ &
619 & tl_fe(i,j+1)-tl_fe(i,j))
620!^ t(i,j,k,nnew,itrc)=t(i,j,k,nnew,itrc)-cff
621!^
622 tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)-tl_cff
623#ifdef DIAGNOSTICS_TS
624!! DiaTwrk(i,j,k,itrc,iThdif)=-cff
625#endif
626 END DO
627 END DO
628 END DO
629 END DO
630!
631 RETURN
632 END SUBROUTINE rp_t3dmix4_s_tile
633
634 END MODULE rp_t3dmix4_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
integer, dimension(:), allocatable istvar
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter irpm
Definition mod_param.F:664
type(t_lbc), dimension(:,:,:), allocatable tl_lbc
Definition mod_param.F:379
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable lm
Definition mod_param.F:455
integer, dimension(:), allocatable mm
Definition mod_param.F:456
real(dp), dimension(:), allocatable dt
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
logical, dimension(:,:), allocatable compositegrid
integer, parameter isouth
integer, parameter ieast
logical, dimension(:,:), allocatable ltracerclm
integer, parameter inorth
integer, dimension(:), allocatable nrhs
integer, dimension(:), allocatable nnew
integer, dimension(:), allocatable nstp
subroutine rp_t3dmix4_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_u, diff3d_v, diff3d_r, diff4, tclm, t, tl_t)
subroutine, public rp_t3dmix4(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