ROMS
Loading...
Searching...
No Matches
ad_t3dmix4_s.h
Go to the documentation of this file.
1 MODULE ad_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 adjoint horizontal biharmonic mixing !
11! of tracers along S-coordinate levels surfaces. !
12! !
13!=======================================================================
14!
15 implicit none
16!
17 PRIVATE
18 PUBLIC ad_t3dmix4
19!
20 CONTAINS
21!
22!***********************************************************************
23 SUBROUTINE ad_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, iadm, 27, __line__, myfile)
51#endif
52 CALL ad_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) % ad_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_Y ET
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) % ad_t)
88#ifdef PROFILE
89 CALL wclock_off (ng, iadm, 27, __line__, myfile)
90#endif
91!
92 RETURN
93 END SUBROUTINE ad_t3dmix4
94!
95!***********************************************************************
96 SUBROUTINE ad_t3dmix4_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, ad_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, ad_t)
124!***********************************************************************
125!
126 USE mod_param
127 USE mod_scalars
128!
129! Imported variable declarations.
130!
131 integer, intent(in) :: ng, tile
132 integer, intent(in) :: LBi, UBi, LBj, UBj
133 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
134 integer, intent(in) :: nrhs, nstp, nnew
135
136#ifdef ASSUMED_SHAPE
137# ifdef MASKING
138 real(r8), intent(in) :: umask(LBi:,LBj:)
139 real(r8), intent(in) :: vmask(LBi:,LBj:)
140# endif
141# ifdef WET_DRY_NOT_YET
142 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
143 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
144# endif
145# ifdef DIFF_3DCOEF
146# ifdef TS_U3ADV_SPLIT_NOT_YET
147 real(r8), intent(in) :: diff3d_u(LBi:,LBj:,:)
148 real(r8), intent(in) :: diff3d_v(LBi:,LBj:,:)
149# else
150 real(r8), intent(in) :: diff3d_r(LBi:,LBj:,:)
151# endif
152# else
153 real(r8), intent(in) :: diff4(LBi:,LBj:,:)
154# endif
155 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
156 real(r8), intent(in) :: pmon_u(LBi:,LBj:)
157 real(r8), intent(in) :: pnom_v(LBi:,LBj:)
158 real(r8), intent(in) :: pm(LBi:,LBj:)
159 real(r8), intent(in) :: pn(LBi:,LBj:)
160 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
161# ifdef TS_MIX_CLIMA
162 real(r8), intent(in) :: tclm(LBi:,LBj:,:,:)
163# endif
164# ifdef DIAGNOSTICS_TS
165!! real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
166# endif
167 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
168 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
169#else
170# ifdef MASKING
171 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
172 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
173# endif
174# ifdef WET_DRY_NOT_YET
175 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
176 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
177# endif
178# ifdef DIFF_3DCOEF
179# ifdef TS_U3ADV_SPLIT_NOT_YET
180 real(r8), intent(in) :: diff3d_u(LBi:UBi,LBj:UBj,N(ng))
181 real(r8), intent(in) :: diff3d_v(LBi:UBi,LBj:UBj,N(ng))
182# else
183 real(r8), intent(in) :: diff3d_r(LBi:UBi,LBj:UBj,N(ng))
184# endif
185# else
186 real(r8), intent(in) :: diff4(LBi:UBi,LBj:UBj,NT(ng))
187# endif
188 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
189 real(r8), intent(in) :: pmon_u(LBi:UBi,LBj:UBj)
190 real(r8), intent(in) :: pnom_v(LBi:UBi,LBj:UBj)
191 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
192 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
193 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
194# ifdef TS_MIX_CLIMA
195 real(r8), intent(in) :: tclm(LBi:UBi,LBj:UBj,N(ng),NT(ng))
196# endif
197# ifdef DIAGNOSTICS_TS
198!! real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
199!! & NDT)
200# endif
201 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
202 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
203#endif
204!
205! Local variable declarations.
206!
207 integer :: Imin, Imax, Jmin, Jmax
208 integer :: i, itrc, j, k
209
210 real(r8) :: cff, cff1
211 real(r8) :: adfac, adfac1, adfac2, adfac3, adfac4, ad_cff, ad_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) :: ad_FE
218 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FX
219 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_LapT
220
221#include "set_bounds.h"
222!
223!-----------------------------------------------------------------------
224! Initialize adjoint private variables.
225!-----------------------------------------------------------------------
226!
227 ad_cff=0.0_r8
228 ad_cff1=0.0_r8
229
230 ad_fe(imins:imaxs,jmins:jmaxs)=0.0_r8
231 ad_fx(imins:imaxs,jmins:jmaxs)=0.0_r8
232
233 ad_lapt(imins:imaxs,jmins:jmaxs)=0.0_r8
234!
235!----------------------------------------------------------------------
236! Compute adjoint horizontal biharmonic diffusion along constant
237! S-surfaces. The biharmonic operator is computed by applying the
238! harmonic operator twice.
239!----------------------------------------------------------------------
240!
241! Set local I- and J-ranges.
242!
243 IF (ewperiodic(ng)) THEN
244 imin=istr-1
245 imax=iend+1
246 ELSE
247 imin=max(istr-1,1)
248 imax=min(iend+1,lm(ng))
249 END IF
250 IF (nsperiodic(ng)) THEN
251 jmin=jstr-1
252 jmax=jend+1
253 ELSE
254 jmin=max(jstr-1,1)
255 jmax=min(jend+1,mm(ng))
256 END IF
257!
258! Compute FX=d(LapT)/d(xi) and FE=d(LapT)/d(eta) BASIC STATE terms.
259!
260 DO itrc=1,nt(ng)
261 DO k=1,n(ng)
262 DO j=jmin,jmax
263 DO i=imin,imax+1
264#ifdef DIFF_3DCOEF
265# ifdef TS_U3ADV_SPLIT_NOT_YET
266 cff=0.5_r8*diff3d_u(i,j,k)*pmon_u(i,j)
267# else
268 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
269 & pmon_u(i,j)
270# endif
271#else
272 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i-1,j,itrc))* &
273 & pmon_u(i,j)
274#endif
275#ifdef MASKING
276 cff=cff*umask(i,j)
277#endif
278#ifdef WET_DRY
279 cff=cff*umask_wet(i,j)
280#endif
281#if defined TS_MIX_STABILITY
282 fx(i,j)=cff* &
283 & (hz(i,j,k)+hz(i-1,j,k))* &
284 & (0.75_r8*(t(i ,j,k,nrhs,itrc)- &
285 & t(i-1,j,k,nrhs,itrc))+ &
286 & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
287 & t(i-1,j,k,nstp,itrc)))
288#elif defined TS_MIX_CLIMA
289 IF (ltracerclm(itrc,ng)) THEN
290 fx(i,j)=cff* &
291 & (hz(i,j,k)+hz(i-1,j,k))* &
292 & ((t(i ,j,k,nrhs,itrc)-tclm(i ,j,k,itrc))- &
293 & (t(i-1,j,k,nrhs,itrc)-tclm(i-1,j,k,itrc)))
294 ELSE
295 fx(i,j)=cff* &
296 & (hz(i,j,k)+hz(i-1,j,k))* &
297 & (t(i ,j,k,nrhs,itrc)- &
298 & t(i-1,j,k,nrhs,itrc))
299 END IF
300#else
301 fx(i,j)=cff* &
302 & (hz(i,j,k)+hz(i-1,j,k))* &
303 & (t(i ,j,k,nrhs,itrc)- &
304 & t(i-1,j,k,nrhs,itrc))
305#endif
306 END DO
307 END DO
308 DO j=jmin,jmax+1
309 DO i=imin,imax
310#ifdef DIFF_3DCOEF
311# ifdef TS_U3ADV_SPLIT_NOT_YET
312 cff=0.5_r8*diff3d_v(i,j,k)*pnom_v(i,j)
313# else
314 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
315 & pnom_v(i,j)
316# endif
317#else
318 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i,j-1,itrc))* &
319 & pnom_v(i,j)
320#endif
321#ifdef MASKING
322 cff=cff*vmask(i,j)
323#endif
324#ifdef WET_DRY
325 cff=cff*vmask_wet(i,j)
326#endif
327#if defined TS_MIX_STABILITY
328 fe(i,j)=cff* &
329 & (hz(i,j,k)+hz(i,j-1,k))* &
330 & (0.75_r8*(t(i,j ,k,nrhs,itrc)- &
331 & t(i,j-1,k,nrhs,itrc))+ &
332 & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
333 & t(i,j-1,k,nstp,itrc)))
334#elif defined TS_MIX_CLIMA
335 IF (ltracerclm(itrc,ng)) THEN
336 fe(i,j)=cff* &
337 & (hz(i,j,k)+hz(i,j-1,k))* &
338 & ((t(i,j ,k,nrhs,itrc)-tclm(i,j ,k,itrc))- &
339 & (t(i,j-1,k,nrhs,itrc)-tclm(i,j-1,k,itrc)))
340 ELSE
341 fe(i,j)=cff* &
342 & (hz(i,j,k)+hz(i,j-1,k))* &
343 & (t(i,j ,k,nrhs,itrc)- &
344 & t(i,j-1,k,nrhs,itrc))
345 END IF
346#else
347 fe(i,j)=cff* &
348 & (hz(i,j,k)+hz(i,j-1,k))* &
349 & (t(i,j ,k,nrhs,itrc)- &
350 & t(i,j-1,k,nrhs,itrc))
351#endif
352 END DO
353 END DO
354!
355! Compute first BASIC STATE harmonic operator and multiply by the
356! metrics of the second harmonic operator.
357!
358 DO j=jmin,jmax
359 DO i=imin,imax
360 cff=1.0_r8/hz(i,j,k)
361 lapt(i,j)=pm(i,j)*pn(i,j)*cff* &
362 & (fx(i+1,j)-fx(i,j)+ &
363 & fe(i,j+1)-fe(i,j))
364 END DO
365 END DO
366!
367! Apply boundary conditions (except periodic; closed or gradient)
368! to the first harmonic operator.
369!
370 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
371 IF (domain(ng)%Western_Edge(tile)) THEN
372 IF (ad_lbc(iwest,istvar(itrc),ng)%closed) THEN
373 DO j=jmin,jmax
374 lapt(istr-1,j)=0.0_r8
375 END DO
376 ELSE
377 DO j=jmin,jmax
378 lapt(istr-1,j)=lapt(istr,j)
379 END DO
380 END IF
381 END IF
382 END IF
383!
384 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
385 IF (domain(ng)%Eastern_Edge(tile)) THEN
386 IF (ad_lbc(ieast,istvar(itrc),ng)%closed) THEN
387 DO j=jmin,jmax
388 lapt(iend+1,j)=0.0_r8
389 END DO
390 ELSE
391 DO j=jmin,jmax
392 lapt(iend+1,j)=lapt(iend,j)
393 END DO
394 END IF
395 END IF
396 END IF
397!
398 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
399 IF (domain(ng)%Southern_Edge(tile)) THEN
400 IF (ad_lbc(isouth,istvar(itrc),ng)%closed) THEN
401 DO i=imin,imax
402 lapt(i,jstr-1)=0.0_r8
403 END DO
404 ELSE
405 DO i=imin,imax
406 lapt(i,jstr-1)=lapt(i,jstr)
407 END DO
408 END IF
409 END IF
410 END IF
411!
412 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
413 IF (domain(ng)%Northern_Edge(tile)) THEN
414 IF (ad_lbc(inorth,istvar(itrc),ng)%closed) THEN
415 DO i=imin,imax
416 lapt(i,jend+1)=0.0_r8
417 END DO
418 ELSE
419 DO i=imin,imax
420 lapt(i,jend+1)=lapt(i,jend)
421 END DO
422 END IF
423 END IF
424 END IF
425!
426! Time-step biharmonic, S-surfaces diffusion term (m Tunits).
427!
428 DO j=jstr,jend
429 DO i=istr,iend
430#ifdef DIAGNOSTICS_TS
431!! DiaTwrk(i,j,k,itrc,iThdif)=-cff
432#endif
433!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)-tl_cff
434!^
435 ad_cff=ad_cff-ad_t(i,j,k,nnew,itrc)
436!^ tl_cff=dt(ng)*pm(i,j)*pn(i,j)* &
437!^ & (tl_FX(i+1,j)-tl_FX(i,j)+ &
438!^ & tl_FE(i,j+1)-tl_FE(i,j))
439!^
440 adfac=dt(ng)*pm(i,j)*pn(i,j)*ad_cff
441 ad_fx(i ,j)=ad_fx(i ,j)-adfac
442 ad_fx(i+1,j)=ad_fx(i+1,j)+adfac
443 ad_fe(i,j )=ad_fe(i,j )-adfac
444 ad_fe(i,j+1)=ad_fe(i,j+1)+adfac
445 ad_cff=0.0_r8
446 END DO
447 END DO
448!
449! Compute ad_FE=d(ad_LapT)/d(eta) and ad_FX=d(ad_LapT)/d(xi) terms.
450!
451 DO j=jstr,jend+1
452 DO i=istr,iend
453#ifdef DIFF_3DCOEF
454# ifdef TS_U3ADV_SPLIT_NOT_YET
455 cff=0.5_r8*diff3d_v(i,j,k)*pnom_v(i,j)
456# else
457 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
458 & pnom_v(i,j)
459# endif
460#else
461 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i,j-1,itrc))* &
462 & pnom_v(i,j)
463#endif
464#ifdef WET_DRY_NOT_YET
465 fe(i,j)=fe(i,j)*vmask_wet(i,j)
466#endif
467#ifdef MASKING
468!^ tl_FE(i,j)=tl_FE(i,j)*vmask(i,j)
469!^
470 ad_fe(i,j)=ad_fe(i,j)*vmask(i,j)
471#endif
472!^ tl_FE(i,j)=cff* &
473!^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* &
474!^ & (LapT(i,j)-LapT(i,j-1))+ &
475!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
476!^ & (tl_LapT(i,j)-tl_LapT(i,j-1)))
477!^
478 adfac=cff*ad_fe(i,j)
479 adfac1=adfac*(lapt(i,j)-lapt(i,j-1))
480 adfac2=adfac*(hz(i,j,k)+hz(i,j-1,k))
481 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac1
482 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac1
483 ad_lapt(i,j-1)=ad_lapt(i,j-1)-adfac2
484 ad_lapt(i,j )=ad_lapt(i,j )+adfac2
485 ad_fe(i,j)=0.0_r8
486 END DO
487 END DO
488 DO j=jstr,jend
489 DO i=istr,iend+1
490#ifdef DIFF_3DCOEF
491# ifdef TS_U3ADV_SPLIT_NOT_YET
492 cff=0.5_r8*diff3d_u(i,j,k)*pmon_u(i,j)
493# else
494 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
495 & pmon_u(i,j)
496# endif
497#else
498 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i-1,j,itrc))* &
499 & pmon_u(i,j)
500#endif
501#ifdef WET_DRY_NOT_YET
502 fx(i,j)=fx(i,j)*umask_wet(i,j)
503#endif
504#ifdef MASKING
505!^ tl_FX(i,j)=tl_FX(i,j)*umask(i,j)
506!^
507 ad_fx(i,j)=ad_fx(i,j)*umask(i,j)
508#endif
509!^ tl_FX(i,j)=cff* &
510!^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* &
511!^ & (LapT(i,j)-LapT(i-1,j))+ &
512!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
513!^ & (tl_LapT(i,j)-tl_LapT(i-1,j)))
514!^
515 adfac=cff*ad_fx(i,j)
516 adfac1=adfac*(lapt(i,j)-lapt(i-1,j))
517 adfac2=adfac*(hz(i,j,k)+hz(i-1,j,k))
518 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac1
519 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac1
520 ad_lapt(i-1,j)=ad_lapt(i-1,j)-adfac2
521 ad_lapt(i ,j)=ad_lapt(i ,j)+adfac2
522 ad_fx(i,j)=0.0_r8
523 END DO
524 END DO
525!
526! Apply boundary conditions (except periodic; closed or gradient)
527! to the first harmonic operator.
528!
529 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
530 IF (domain(ng)%Northern_Edge(tile)) THEN
531 IF (ad_lbc(inorth,istvar(itrc),ng)%closed) THEN
532 DO i=imin,imax
533!^ tl_LapT(i,Jend+1)=0.0_r8
534!^
535 ad_lapt(i,jend+1)=0.0_r8
536 END DO
537 ELSE
538 DO i=imin,imax
539!^ tl_LapT(i,Jend+1)=tl_LapT(i,Jend)
540!^
541 ad_lapt(i,jend)=ad_lapt(i,jend)+ad_lapt(i,jend+1)
542 ad_lapt(i,jend+1)=0.0_r8
543 END DO
544 END IF
545 END IF
546 END IF
547!
548 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
549 IF (domain(ng)%Southern_Edge(tile)) THEN
550 IF (ad_lbc(isouth,istvar(itrc),ng)%closed) THEN
551 DO i=imin,imax
552!^ tl_LapT(i,Jstr-1)=0.0_r8
553!^
554 ad_lapt(i,jstr-1)=0.0_r8
555 END DO
556 ELSE
557 DO i=imin,imax
558!^ tl_LapT(i,Jstr-1)=tl_LapT(i,Jstr)
559!^
560 ad_lapt(i,jstr)=ad_lapt(i,jstr)+ad_lapt(i,jstr-1)
561 ad_lapt(i,jstr-1)=0.0_r8
562 END DO
563 END IF
564 END IF
565 END IF
566!
567 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
568 IF (domain(ng)%Eastern_Edge(tile)) THEN
569 IF (ad_lbc(ieast,istvar(itrc),ng)%closed) THEN
570 DO j=jmin,jmax
571!^ tl_LapT(Iend+1,j)=0.0_r8
572!^
573 ad_lapt(iend+1,j)=0.0_r8
574 END DO
575 ELSE
576 DO j=jmin,jmax
577!^ tl_LapT(Iend+1,j)=tl_LapT(Iend,j)
578!^
579 ad_lapt(iend,j)=ad_lapt(iend,j)+ad_lapt(iend+1,j)
580 ad_lapt(iend+1,j)=0.0_r8
581 END DO
582 END IF
583 END IF
584 END IF
585!
586 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
587 IF (domain(ng)%Western_Edge(tile)) THEN
588 IF (ad_lbc(iwest,istvar(itrc),ng)%closed) THEN
589 DO j=jmin,jmax
590!^ tl_LapT(Istr-1,j)=0.0_r8
591!^
592 ad_lapt(istr-1,j)=0.0_r8
593 END DO
594 ELSE
595 DO j=jmin,jmax
596!^ tl_LapT(Istr-1,j)=tl_LapT(Istr,j)
597!^
598 ad_lapt(istr,j)=ad_lapt(istr,j)+ad_lapt(istr-1,j)
599 ad_lapt(istr-1,j)=0.0_r8
600 END DO
601 END IF
602 END IF
603 END IF
604!
605! Compute first harmonic operator and multiply by the metrics of the
606! second harmonic operator.
607!
608 DO j=jmin,jmax
609 DO i=imin,imax
610 cff=1.0_r8/hz(i,j,k)
611!^ tl_LapT(i,j)=pm(i,j)*pn(i,j)* &
612!^ & (tl_cff* &
613!^ & (FX(i+1,j)-FX(i,j)+ &
614!^ & FE(i,j+1)-FE(i,j))+ &
615!^ & cff* &
616!^ & (tl_FX(i+1,j)-tl_FX(i,j)+ &
617!^ & tl_FE(i,j+1)-tl_FE(i,j)))
618!^
619 adfac=pm(i,j)*pn(i,j)*ad_lapt(i,j)
620 adfac1=cff*adfac
621 ad_cff=ad_cff+ &
622 & adfac*(fx(i+1,j)-fx(i,j)+ &
623 & fe(i,j+1)-fe(i,j))
624 ad_fx(i ,j)=ad_fx(i ,j)-adfac1
625 ad_fx(i+1,j)=ad_fx(i+1,j)+adfac1
626 ad_fe(i,j )=ad_fe(i,j )-adfac1
627 ad_fe(i,j+1)=ad_fe(i,j+1)+adfac1
628 ad_lapt(i,j)=0.0_r8
629!^ tl_cff=-cff*cff*tl_Hz(i,j,k)
630!^
631 ad_hz(i,j,k)=ad_hz(i,j,k)-cff*cff*ad_cff
632 ad_cff=0.0_r8
633 END DO
634 END DO
635!
636! Compute horizontal tracer flux in the ETA- and XI-directions.
637!
638 DO j=jmin,jmax+1
639 DO i=imin,imax
640#ifdef DIFF_3DCOEF
641# ifdef TS_U3ADV_SPLIT_NOT_YET
642 cff=0.5_r8*diff3d_v(i,j,k)*pnom_v(i,j)
643# else
644 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i,j-1,k))* &
645 & pnom_v(i,j)
646# endif
647#else
648 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i,j-1,itrc))* &
649 & pnom_v(i,j)
650#endif
651#ifdef WET_DRY_NOT_YET
652 cff=cff*vmask_wet(i,j)
653#endif
654#ifdef MASKING
655 cff=cff*vmask(i,j)
656#endif
657#if defined TS_MIX_STABILITY
658!^ tl_FE(i,j)=cff* &
659!^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* &
660!^ & (0.75_r8*(t(i,j ,k,nrhs,itrc)- &
661!^ & t(i,j-1,k,nrhs,itrc))+ &
662!^ & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
663!^ & t(i,j-1,k,nstp,itrc)))+ &
664!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
665!^ & (0.75_r8*(tl_t(i,j ,k,nrhs,itrc)- &
666!^ & tl_t(i,j-1,k,nrhs,itrc))+ &
667!^ & 0.25_r8*(tl_t(i,j ,k,nstp,itrc)- &
668!^ & tl_t(i,j-1,k,nstp,itrc))))
669!^
670 adfac=cff*ad_fe(i,j)
671 adfac1=adfac*(0.75_r8*(t(i,j ,k,nrhs,itrc)- &
672 & t(i,j-1,k,nrhs,itrc))+ &
673 & 0.25_r8*(t(i,j ,k,nstp,itrc)- &
674 & t(i,j-1,k,nstp,itrc)))
675 adfac2=adfac*(hz(i,j,k)+hz(i,j-1,k))
676 adfac3=adfac2*0.75_r8
677 adfac4=adfac2*0.25_r8
678 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac1
679 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac1
680 ad_t(i,j-1,k,nrhs,itrc)=ad_t(i,j-1,k,nrhs,itrc)-adfac3
681 ad_t(i,j ,k,nrhs,itrc)=ad_t(i,j ,k,nrhs,itrc)+adfac3
682 ad_t(i,j-1,k,nstp,itrc)=ad_t(i,j-1,k,nstp,itrc)-adfac4
683 ad_t(i,j ,k,nstp,itrc)=ad_t(i,j ,k,nstp,itrc)+adfac4
684 ad_fe(i,j)=0.0_r8
685#elif defined TS_MIX_CLIMA
686 IF (ltracerclm(itrc,ng)) THEN
687!^ tl_FE(i,j)=cff* &
688!^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* &
689!^ & ((t(i,j ,k,nrhs,itrc)- &
690!^ & tclm(i,j ,k,itrc))- &
691!^ & (t(i,j-1,k,nrhs,itrc)- &
692!^ & tclm(i,j-1,k,itrc)))+ &
693!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
694!^ & (tl_t(i,j ,k,nrhs,itrc)- &
695!^ & tl_t(i,j-1,k,nrhs,itrc)))
696!^
697 adfac=cff*ad_fe(i,j)
698 adfac1=adfac*((t(i,j ,k,nrhs,itrc)- &
699 & tclm(i,j ,k,itrc))- &
700 & (t(i,j-1,k,nrhs,itrc)- &
701 & tclm(i,j-1,k,itrc)))
702 adfac2=adfac*(hz(i,j,k)+hz(i,j-1,k))
703 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac1
704 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac1
705 ad_t(i,j-1,k,nrhs,itrc)=ad_t(i,j-1,k,nrhs,itrc)-adfac2
706 ad_t(i,j ,k,nrhs,itrc)=ad_t(i,j ,k,nrhs,itrc)+adfac2
707 ad_fe(i,j)=0.0_r8
708 ELSE
709!^ tl_FE(i,j)=cff* &
710!^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* &
711!^ & (t(i,j ,k,nrhs,itrc)- &
712!^ & t(i,j-1,k,nrhs,itrc))+ &
713!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
714!^ & (tl_t(i,j ,k,nrhs,itrc)- &
715!^ & tl_t(i,j-1,k,nrhs,itrc)))
716!^
717 adfac=cff*ad_fe(i,j)
718 adfac1=adfac*(t(i,j ,k,nrhs,itrc)- &
719 & t(i,j-1,k,nrhs,itrc))
720 adfac2=adfac*(hz(i,j,k)+hz(i,j-1,k))
721 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac1
722 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac1
723 ad_t(i,j-1,k,nrhs,itrc)=ad_t(i,j-1,k,nrhs,itrc)-adfac2
724 ad_t(i,j ,k,nrhs,itrc)=ad_t(i,j ,k,nrhs,itrc)+adfac2
725 ad_fe(i,j)=0.0_r8
726 END IF
727#else
728!^ tl_FE(i,j)=cff* &
729!^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* &
730!^ & (t(i,j ,k,nrhs,itrc)- &
731!^ & t(i,j-1,k,nrhs,itrc))+ &
732!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
733!^ & (tl_t(i,j ,k,nrhs,itrc)- &
734!^ & tl_t(i,j-1,k,nrhs,itrc)))
735!^
736 adfac=cff*ad_fe(i,j)
737 adfac1=adfac*(t(i,j ,k,nrhs,itrc)- &
738 & t(i,j-1,k,nrhs,itrc))
739 adfac2=adfac*(hz(i,j,k)+hz(i,j-1,k))
740 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac1
741 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac1
742 ad_t(i,j-1,k,nrhs,itrc)=ad_t(i,j-1,k,nrhs,itrc)-adfac2
743 ad_t(i,j ,k,nrhs,itrc)=ad_t(i,j ,k,nrhs,itrc)+adfac2
744 ad_fe(i,j)=0.0_r8
745#endif
746 END DO
747 END DO
748
749 DO j=jmin,jmax
750 DO i=imin,imax+1
751#ifdef DIFF_3DCOEF
752# ifdef TS_U3ADV_SPLIT_NOT_YET
753 cff=0.5_r8*diff3d_u(i,j,k)*pmon_u(i,j)
754# else
755 cff=0.25_r8*(diff3d_r(i,j,k)+diff3d_r(i-1,j,k))* &
756 & pmon_u(i,j)
757# endif
758#else
759 cff=0.25_r8*(diff4(i,j,itrc)+diff4(i-1,j,itrc))* &
760 & pmon_u(i,j)
761#endif
762#ifdef WET_DRY_NOT_YET
763 cff=cff*umask_wet(i,j)
764#endif
765#ifdef MASKING
766 cff=cff*umask(i,j)
767#endif
768#if defined TS_MIX_STABILITY
769!^ tl_FX(i,j)=cff* &
770!^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* &
771!^ & (0.75_r8*(t(i ,j,k,nrhs,itrc)- &
772!^ & t(i-1,j,k,nrhs,itrc))+ &
773!^ & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
774!^ & t(i-1,j,k,nstp,itrc)))+ &
775!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
776!^ & (0.75_r8*(tl_t(i ,j,k,nrhs,itrc)- &
777!^ & tl_t(i-1,j,k,nrhs,itrc))+ &
778!^ & 0.25_r8*(tl_t(i ,j,k,nstp,itrc)- &
779!^ & tl_t(i-1,j,k,nstp,itrc))))
780!^
781 adfac=cff*ad_fx(i,j)
782 adfac1=adfac*(0.75_r8*(t(i ,j,k,nrhs,itrc)- &
783 & t(i-1,j,k,nrhs,itrc))+ &
784 & 0.25_r8*(t(i ,j,k,nstp,itrc)- &
785 & t(i-1,j,k,nstp,itrc)))
786 adfac2=adfac*(hz(i,j,k)+hz(i-1,j,k))
787 adfac3=adfac2*0.75_r8
788 adfac4=adfac2*0.25_r8
789 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac1
790 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac1
791 ad_t(i-1,j,k,nrhs,itrc)=ad_t(i-1,j,k,nrhs,itrc)-adfac3
792 ad_t(i ,j,k,nrhs,itrc)=ad_t(i ,j,k,nrhs,itrc)+adfac3
793 ad_t(i-1,j,k,nstp,itrc)=ad_t(i-1,j,k,nstp,itrc)-adfac4
794 ad_t(i ,j,k,nstp,itrc)=ad_t(i ,j,k,nstp,itrc)+adfac4
795 ad_fx(i,j)=0.0_r8
796#elif defined TS_MIX_CLIMA
797 IF (ltracerclm(itrc,ng)) THEN
798!^ tl_FX(i,j)=cff* &
799!^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* &
800!^ & ((t(i ,j,k,nrhs,itrc)- &
801!^ & tclm(i ,j,k,itrc))- &
802!^ & (t(i-1,j,k,nrhs,itrc)- &
803!^ & tclm(i-1,j,k,itrc)))+ &
804!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
805!^ & (tl_t(i ,j,k,nrhs,itrc)- &
806!^ & tl_t(i-1,j,k,nrhs,itrc)))
807!^
808 adfac=cff*ad_fx(i,j)
809 adfac1=adfac*((t(i ,j,k,nrhs,itrc)- &
810 & tclm(i ,j,k,itrc))- &
811 & (t(i-1,j,k,nrhs,itrc)- &
812 & tclm(i-1,j,k,itrc)))
813 adfac2=adfac*(hz(i,j,k)+hz(i-1,j,k))
814 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac1
815 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac1
816 ad_t(i-1,j,k,nrhs,itrc)=ad_t(i-1,j,k,nrhs,itrc)-adfac2
817 ad_t(i ,j,k,nrhs,itrc)=ad_t(i ,j,k,nrhs,itrc)+adfac2
818 ad_fx(i,j)=0.0_r8
819 ELSE
820!^ tl_FX(i,j)=cff* &
821!^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* &
822!^ & (t(i ,j,k,nrhs,itrc)- &
823!^ & t(i-1,j,k,nrhs,itrc))+ &
824!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
825!^ & (tl_t(i ,j,k,nrhs,itrc)- &
826!^ & tl_t(i-1,j,k,nrhs,itrc)))
827!^
828 adfac=cff*ad_fx(i,j)
829 adfac1=adfac*(t(i ,j,k,nrhs,itrc)- &
830 & t(i-1,j,k,nrhs,itrc))
831 adfac2=adfac*(hz(i,j,k)+hz(i-1,j,k))
832 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac1
833 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac1
834 ad_t(i-1,j,k,nrhs,itrc)=ad_t(i-1,j,k,nrhs,itrc)-adfac2
835 ad_t(i ,j,k,nrhs,itrc)=ad_t(i ,j,k,nrhs,itrc)+adfac2
836 ad_fx(i,j)=0.0_r8
837 END IF
838#else
839!^ tl_FX(i,j)=cff* &
840!^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* &
841!^ & (t(i ,j,k,nrhs,itrc)- &
842!^ & t(i-1,j,k,nrhs,itrc))+ &
843!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
844!^ & (tl_t(i ,j,k,nrhs,itrc)- &
845!^ & tl_t(i-1,j,k,nrhs,itrc)))
846!^
847 adfac=cff*ad_fx(i,j)
848 adfac1=adfac*(t(i ,j,k,nrhs,itrc)- &
849 & t(i-1,j,k,nrhs,itrc))
850 adfac2=adfac*(hz(i,j,k)+hz(i-1,j,k))
851 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac1
852 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac1
853 ad_t(i-1,j,k,nrhs,itrc)=ad_t(i-1,j,k,nrhs,itrc)-adfac2
854 ad_t(i ,j,k,nrhs,itrc)=ad_t(i ,j,k,nrhs,itrc)+adfac2
855 ad_fx(i,j)=0.0_r8
856#endif
857 END DO
858 END DO
859 END DO
860 END DO
861!
862 RETURN
863 END SUBROUTINE ad_t3dmix4_s_tile
864
865 END MODULE ad_t3dmix4_mod
subroutine, public ad_t3dmix4(ng, tile)
subroutine ad_t3dmix4_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nstp, nnew, umask, vmask, umask_wet, vmask_wet, hz, ad_hz, pmon_u, pnom_v, pm, pn, diff3d_u, diff3d_v, diff3d_r, diff4, tclm, t, ad_t)
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
type(t_lbc), dimension(:,:,:), allocatable ad_lbc
Definition mod_param.F:378
integer, parameter iadm
Definition mod_param.F:665
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
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