ROMS
Loading...
Searching...
No Matches
ad_step3d_t_mod Module Reference

Functions/Subroutines

subroutine, public ad_step3d_t (ng, tile)
 
subroutine ad_step3d_t_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nstp, nnew, rmask, umask, vmask, rmask_wet, umask_wet, vmask_wet, omn, om_u, om_v, on_u, on_v, pm, pn, hz, ad_hz, huon, ad_huon, hvom, ad_hvom, z_r, ad_z_r, akt, ad_akt, w, ad_w, t, ad_t)
 

Function/Subroutine Documentation

◆ ad_step3d_t()

subroutine, public ad_step3d_t_mod::ad_step3d_t ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 46 of file ad_step3d_t.F.

47!***********************************************************************
48!
49 USE mod_param
50# ifdef DIAGNOSTICS_TS
51!! USE mod_diags
52# endif
53 USE mod_grid
54 USE mod_mixing
55 USE mod_ocean
56 USE mod_stepping
57!
58! Imported variable declarations.
59!
60 integer, intent(in) :: ng, tile
61!
62! Local variable declarations.
63!
64 character (len=*), parameter :: MyFile = &
65 & __FILE__
66!
67# include "tile.h"
68!
69# ifdef PROFILE
70 CALL wclock_on (ng, iadm, 35, __line__, myfile)
71# endif
72 CALL ad_step3d_t_tile (ng, tile, &
73 & lbi, ubi, lbj, ubj, &
74 & imins, imaxs, jmins, jmaxs, &
75 & nrhs(ng), nstp(ng), nnew(ng), &
76# ifdef MASKING
77 & grid(ng) % rmask, &
78 & grid(ng) % umask, &
79 & grid(ng) % vmask, &
80# endif
81# ifdef WET_DRY
82 & grid(ng) % rmask_wet, &
83 & grid(ng) % umask_wet, &
84 & grid(ng) % vmask_wet, &
85# endif
86 & grid(ng) % omn, &
87 & grid(ng) % om_u, &
88 & grid(ng) % om_v, &
89 & grid(ng) % on_u, &
90 & grid(ng) % on_v, &
91 & grid(ng) % pm, &
92 & grid(ng) % pn, &
93 & grid(ng) % Hz, &
94 & grid(ng) % ad_Hz, &
95 & grid(ng) % Huon, &
96 & grid(ng) % ad_Huon, &
97 & grid(ng) % Hvom, &
98 & grid(ng) % ad_Hvom, &
99 & grid(ng) % z_r, &
100 & grid(ng) % ad_z_r, &
101 & mixing(ng) % Akt, &
102 & mixing(ng) % ad_Akt, &
103 & ocean(ng) % W, &
104 & ocean(ng) % ad_W, &
105# if defined FLOATS_NOT_YET && defined FLOAT_VWALK
106 & mixing(ng) % dAktdz, &
107# endif
108# ifdef DIAGNOSTICS_TS
109!! & DIAGS(ng) % DiaTwrk, &
110# endif
111 & ocean(ng) % t, &
112 & ocean(ng) % ad_t)
113# ifdef PROFILE
114 CALL wclock_off (ng, iadm, 35, __line__, myfile)
115# endif
116!
117 RETURN
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 iadm
Definition mod_param.F:665
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

References ad_step3d_t_tile(), mod_grid::grid, mod_param::iadm, mod_mixing::mixing, mod_stepping::nnew, mod_stepping::nrhs, mod_stepping::nstp, mod_ocean::ocean, wclock_off(), and wclock_on().

Referenced by ad_main3d().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ad_step3d_t_tile()

subroutine ad_step3d_t_mod::ad_step3d_t_tile ( integer, intent(in) ng,
integer, intent(in) tile,
integer, intent(in) lbi,
integer, intent(in) ubi,
integer, intent(in) lbj,
integer, intent(in) ubj,
integer, intent(in) imins,
integer, intent(in) imaxs,
integer, intent(in) jmins,
integer, intent(in) jmaxs,
integer, intent(in) nrhs,
integer, intent(in) nstp,
integer, intent(in) nnew,
real(r8), dimension(lbi:,lbj:), intent(in) rmask,
real(r8), dimension(lbi:,lbj:), intent(in) umask,
real(r8), dimension(lbi:,lbj:), intent(in) vmask,
real(r8), dimension(lbi:,lbj:), intent(in) rmask_wet,
real(r8), dimension(lbi:,lbj:), intent(in) umask_wet,
real(r8), dimension(lbi:,lbj:), intent(in) vmask_wet,
real(r8), dimension(lbi:,lbj:), intent(in) omn,
real(r8), dimension(lbi:,lbj:), intent(in) om_u,
real(r8), dimension(lbi:,lbj:), intent(in) om_v,
real(r8), dimension(lbi:,lbj:), intent(in) on_u,
real(r8), dimension(lbi:,lbj:), intent(in) on_v,
real(r8), dimension(lbi:,lbj:), intent(in) pm,
real(r8), dimension(lbi:,lbj:), intent(in) pn,
real(r8), dimension(lbi:,lbj:,:), intent(in) hz,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_hz,
real(r8), dimension(lbi:,lbj:,:), intent(in) huon,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_huon,
real(r8), dimension(lbi:,lbj:,:), intent(in) hvom,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_hvom,
real(r8), dimension(lbi:,lbj:,:), intent(in) z_r,
real(r8), dimension(lbi:,lbj:,:), intent(inout) ad_z_r,
real(r8), dimension(lbi:,lbj:,0:,:), intent(in) akt,
real(r8), dimension(lbi:,lbj:,0:,:), intent(inout) ad_akt,
real(r8), dimension(lbi:,lbj:,0:), intent(in) w,
real(r8), dimension(lbi:,lbj:,0:), intent(inout) ad_w,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(in) t,
real(r8), dimension(lbi:,lbj:,:,:,:), intent(inout) ad_t )
private

Definition at line 121 of file ad_step3d_t.F.

147!***********************************************************************
148!
149 USE mod_param
150 USE mod_clima
151 USE mod_ncparam
152 USE mod_scalars
153 USE mod_sources
154!
156# ifdef DISTRIBUTE
159# endif
160 USE ad_t3dbc_mod, ONLY : ad_t3dbc_tile
161!
162! Imported variable declarations.
163!
164 integer, intent(in) :: ng, tile
165 integer, intent(in) :: LBi, UBi, LBj, UBj
166 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
167 integer, intent(in) :: nrhs, nstp, nnew
168!
169# ifdef ASSUMED_SHAPE
170# ifdef MASKING
171 real(r8), intent(in) :: rmask(LBi:,LBj:)
172 real(r8), intent(in) :: umask(LBi:,LBj:)
173 real(r8), intent(in) :: vmask(LBi:,LBj:)
174# endif
175# ifdef WET_DRY
176 real(r8), intent(in) :: rmask_wet(LBi:,LBj:)
177 real(r8), intent(in) :: umask_wet(LBi:,LBj:)
178 real(r8), intent(in) :: vmask_wet(LBi:,LBj:)
179# endif
180 real(r8), intent(in) :: omn(LBi:,LBj:)
181 real(r8), intent(in) :: om_u(LBi:,LBj:)
182 real(r8), intent(in) :: om_v(LBi:,LBj:)
183 real(r8), intent(in) :: on_u(LBi:,LBj:)
184 real(r8), intent(in) :: on_v(LBi:,LBj:)
185 real(r8), intent(in) :: pm(LBi:,LBj:)
186 real(r8), intent(in) :: pn(LBi:,LBj:)
187 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
188 real(r8), intent(in) :: Huon(LBi:,LBj:,:)
189 real(r8), intent(in) :: Hvom(LBi:,LBj:,:)
190 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
191# ifdef SUN
192 real(r8), intent(in) :: Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
193 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
194# else
195 real(r8), intent(in) :: Akt(LBi:,LBj:,0:,:)
196 real(r8), intent(in) :: t(LBi:,LBj:,:,:,:)
197# endif
198 real(r8), intent(in) :: W(LBi:,LBj:,0:)
199
200# ifdef DIAGNOSTICS_TS
201!! real(r8), intent(inout) :: DiaTwrk(LBi:,LBj:,:,:,:)
202# endif
203 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
204 real(r8), intent(inout) :: ad_Huon(LBi:,LBj:,:)
205 real(r8), intent(inout) :: ad_Hvom(LBi:,LBj:,:)
206 real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:)
207# ifdef SUN
208 real(r8), intent(inout) :: ad_Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
209# else
210 real(r8), intent(inout) :: ad_Akt(LBi:,LBj:,0:,:)
211# endif
212 real(r8), intent(inout) :: ad_W(LBi:,LBj:,0:)
213# ifdef SUN
214 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
215# else
216 real(r8), intent(inout) :: ad_t(LBi:,LBj:,:,:,:)
217# endif
218# if defined FLOATS_NOT_YET && defined FLOAT_VWALK
219 real(r8), intent(out) :: dAktdz(LBi:,LBj:,:)
220# endif
221
222# else
223
224# ifdef MASKING
225 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
226 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
227 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
228# endif
229# ifdef WET_DRY
230 real(r8), intent(in) :: rmask_wet(LBi:UBi,LBj:UBj)
231 real(r8), intent(in) :: umask_wet(LBi:UBi,LBj:UBj)
232 real(r8), intent(in) :: vmask_wet(LBi:UBi,LBj:UBj)
233# endif
234 real(r8), intent(in) :: omn(LBi:UBi,LBj:UBj)
235 real(r8), intent(in) :: om_u(LBi:UBi,LBj:UBj)
236 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
237 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
238 real(r8), intent(in) :: on_v(LBi:UBi,LBj:UBj)
239 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
240 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
241 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
242 real(r8), intent(in) :: Huon(LBi:UBi,LBj:UBj,N(ng))
243 real(r8), intent(in) :: Hvom(LBi:UBi,LBj:UBj,N(ng))
244 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
245 real(r8), intent(in) :: Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
246 real(r8), intent(in) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
247 real(r8), intent(in) :: W(LBi:UBi,LBj:UBj,0:N(ng))
248
249# ifdef DIAGNOSTICS_TS
250!! real(r8), intent(inout) :: DiaTwrk(LBi:UBi,LBj:UBj,N(ng),NT(ng), &
251!! & NDT)
252# endif
253 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
254 real(r8), intent(inout) :: ad_Huon(LBi:UBi,LBj:UBj,N(ng))
255 real(r8), intent(inout) :: ad_Hvom(LBi:UBi,LBj:UBj,N(ng))
256 real(r8), intent(inout) :: ad_z_r(LBi:UBi,LBj:UBj,N(ng))
257 real(r8), intent(inout) :: ad_Akt(LBi:UBi,LBj:UBj,0:N(ng),NAT)
258 real(r8), intent(inout) :: ad_W(LBi:UBi,LBj:UBj,0:N(ng))
259 real(r8), intent(inout) :: ad_t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
260# if defined FLOATS_NOT_YET && defined FLOAT_VWALK
261 real(r8), intent(out) :: dAktdz(LBi:UBi,LBj:UBj,N(ng))
262# endif
263# endif
264!
265! Local variable declarations.
266!
267 logical :: LapplySrc, Lhsimt
268!
269 integer :: JminT, JmaxT
270 integer :: Isrc, Jsrc
271 integer :: i, ic, is, itrc, j, k, ltrc
272# if defined AGE_MEAN && defined T_PASSIVE
273 integer :: iage
274# endif
275# ifdef DIAGNOSTICS_TS
276 integer :: idiag
277# endif
278
279 real(r8), parameter :: eps = 1.0e-16_r8
280
281 real(r8) :: cff, cff1, cff2, cff3
282 real(r8) :: ad_cff, ad_cff1, ad_cff2, ad_cff3
283 real(r8) :: adfac, adfac1, adfac2
284
285 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: CF
286 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: BC
287 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: DC
288# ifdef SPLINES_VDIFF
289 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: DC1
290# endif
291 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: FC
292
293 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_CF
294 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_BC
295 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_DC
296 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_FC
297
298 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FE
299 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FX
300 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: curv
301 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: grad
302
303 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FE
304 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FX
305 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_curv
306 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_grad
307
308 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: oHz
309 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: ad_oHz
310
311# include "set_bounds.h"
312!
313!-----------------------------------------------------------------------
314! Initialize adjoint private variables.
315!-----------------------------------------------------------------------
316!
317 lhsimt =any(ad_hadvection(:,ng)%HSIMT).and. &
318 & any(ad_vadvection(:,ng)%HSIMT)
319!
320! Initialize.
321!
322 ad_cff=0.0_r8
323 ad_cff1=0.0_r8
324 ad_cff2=0.0_r8
325 ad_cff3=0.0_r8
326 DO j=jmins,jmaxs
327 DO i=imins,imaxs
328 ad_fe(i,j)=0.0_r8
329 ad_fx(i,j)=0.0_r8
330 ad_curv(i,j)=0.0_r8
331 ad_grad(i,j)=0.0_r8
332 END DO
333 DO k=1,n(ng)
334 DO i=imins,imaxs
335 ad_ohz(i,j,k)=0.0_r8
336 END DO
337 END DO
338 END DO
339 DO k=0,n(ng)
340 DO i=imins,imaxs
341 ad_cf(i,k)=0.0_r8
342 ad_bc(i,k)=0.0_r8
343 ad_dc(i,k)=0.0_r8
344 ad_fc(i,k)=0.0_r8
345 END DO
346 END DO
347
348# if defined FLOATS_NOT_YET && defined FLOAT_VWALK
349!
350!-----------------------------------------------------------------------
351! Compute vertical gradient in vertical T-diffusion coefficient for
352! floats random walk.
353!-----------------------------------------------------------------------
354!
355! Exchange boundary data.
356!
357# ifdef DISTRIBUTE
358 CALL mp_exchange3d (ng, tile, inlm, 1, &
359 & lbi, ubi, lbj, ubj, 1, n(ng), &
360 & nghostpoints, &
361 & ewperiodic(ng), nsperiodic(ng), &
362 & daktdz)
363# endif
364 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
365 CALL exchange_r3d_tile (ng, tile, &
366 & lbi, ubi, lbj, ubj, 1, n(ng), &
367 & daktdz)
368 END IF
369!
370 DO j=jstrr,jendr
371 DO i=istrr,iendr
372 DO k=1,n(ng)
373 daktdz(i,j,k)=(akt(i,j,k,1)-akt(i,j,k-1,1))/hz(i,j,k)
374 END DO
375 END DO
376 END DO
377# endif
378!
379!-----------------------------------------------------------------------
380! Apply adjoint lateral boundary conditions and, if appropriate, nudge
381! to tracer data and apply Land/Sea mask.
382!-----------------------------------------------------------------------
383!
384# ifdef DISTRIBUTE
385! Adjoint of exchange boundary data.
386!
387!^ CALL mp_exchange4d (ng, tile, iTLM, 1, &
388!^ & LBi, UBi, LBj, UBj, 1, N(ng), 1, NT(ng), &
389!^ & NghostPoints, &
390!^ & EWperiodic(ng), NSperiodic(ng), &
391!^ & tl_t(:,:,:,nnew,:))
392!^
393 CALL ad_mp_exchange4d (ng, tile, iadm, 1, &
394 & lbi, ubi, lbj, ubj, 1, n(ng), 1, nt(ng), &
395 & nghostpoints, &
396 & ewperiodic(ng), nsperiodic(ng), &
397 & ad_t(:,:,:,nnew,:))
398!
399# endif
400!
401! Initialize tracer counter index. The "tclm" array is only allocated
402! to the NTCLM fields that need to be processed. This is done to
403! reduce memory.
404!
405 ic=0
406!
407 DO itrc=1,nt(ng)
408!
409! Set compact reduced memory tracer index for nudging coefficients and
410! climatology arrays.
411!
412 IF (ltracerclm(itrc,ng).and.lnudgetclm(itrc,ng)) THEN
413 ic=ic+1
414 END IF
415!
416! Apply adjoint periodic boundary conditions.
417!
418 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
419!^ CALL exchange_r3d_tile (ng, tile, &
420!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
421!^ & tl_t(:,:,:,nnew,itrc))
422!^
423 CALL ad_exchange_r3d_tile (ng, tile, &
424 & lbi, ubi, lbj, ubj, 1, n(ng), &
425 & ad_t(:,:,:,nnew,itrc))
426 END IF
427
428# ifdef DIAGNOSTICS_TS
429!!
430!! Compute time-rate-of-change diagnostic term.
431!!
432!! DO k=1,N(ng)
433!! DO j=JstrR,JendR
434!! DO i=IstrR,IendR
435!! DiaTwrk(i,j,k,itrc,iTrate)=t(i,j,k,nnew,itrc)- &
436!! & t(i,j,k,nstp,itrc)
437!! DiaTwrk(i,j,k,itrc,iTrate)=t(i,j,k,nnew,itrc)- &
438!! & DiaTwrk(i,j,k,itrc,iTrate)
439!! END DO
440!! END DO
441!! END DO
442# endif
443# ifdef MASKING
444!
445! Apply Land/Sea mask.
446!
447 DO k=1,n(ng)
448 DO j=jstrr,jendr
449 DO i=istrr,iendr
450!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)*rmask(i,j)
451!^
452 ad_t(i,j,k,nnew,itrc)=ad_t(i,j,k,nnew,itrc)*rmask(i,j)
453 END DO
454 END DO
455 END DO
456# endif
457!
458! Adjoint of nudge towards tracer climatology.
459!
460 IF (ltracerclm(itrc,ng).and.lnudgetclm(itrc,ng)) THEN
461 DO k=1,n(ng)
462 DO j=jstrr,jendr
463 DO i=istrr,iendr
464!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)- &
465!^ & dt(ng)* &
466!^ & CLIMA(ng)%Tnudgcof(i,j,k,ic)* &
467!^ & tl_t(i,j,k,nnew,itrc)
468!^
469 ad_t(i,j,k,nnew,itrc)=ad_t(i,j,k,nnew,itrc)- &
470 & dt(ng)* &
471 & clima(ng)%Tnudgcof(i,j,k,ic)* &
472 & ad_t(i,j,k,nnew,itrc)
473 END DO
474 END DO
475 END DO
476 END IF
477!
478! Set adjoint lateral boundary conditions.
479!
480!^ CALL tl_t3dbc_tile (ng, tile, itrc, ic, &
481!^ & LBi, UBi, LBj, UBj, N(ng), NT(ng), &
482!^ & IminS, ImaxS, JminS, JmaxS, &
483!^ & nstp, nnew, &
484!^ & tl_t)
485!^
486 CALL ad_t3dbc_tile (ng, tile, itrc, ic, &
487 & lbi, ubi, lbj, ubj, n(ng), nt(ng), &
488 & imins, imaxs, jmins, jmaxs, &
489 & nstp, nnew, &
490 & ad_t)
491 END DO
492
493# if defined AGE_MEAN && defined T_PASSIVE
494!
495!-----------------------------------------------------------------------
496! If inert passive tracer and Mean Age, compute age concentration (even
497! inert index) forced by the right-hand-side term that is concentration
498! of an associated conservative passive tracer (odd inert index). Mean
499! Age is age concentration divided by conservative passive tracer
500! concentration. Code implements NPT/2 mean age tracer pairs.
501!
502! Implemented and tested by W.G. Zhang and J. Wilkin. See following
503! reference for details.
504!
505! Zhang et al. (2010): Simulation of water age and residence time in
506! the New York Bight, JPO, 40,965-982, doi:10.1175/2009JPO4249.1
507!-----------------------------------------------------------------------
508!
509 DO itrc=1,npt,2
510 iage=inert(itrc+1) ! even inert tracer index
511 DO k=1,n(ng)
512 DO j=jstr,jend
513 DO i=istr,iend
514 IF ((ad_hadvection(inert(itrc),ng)%MPDATA).and. &
515 & (ad_vadvection(inert(itrc),ng)%MPDATA)) THEN
516 CONTINUE ! not supported
517 ELSE
518!^ tl_t(i,j,k,nnew,iage)=tl_t(i,j,k,nnew,iage)+ &
519!^ & dt(ng)* &
520!^ & tl_t(i,j,k,3,inert(itrc))
521!^
522 ad_t(i,j,k,3,inert(itrc))=ad_t(i,j,k,3,inert(itrc))+ &
523 & dt(ng)*ad_t(i,j,k,nnew,iage)
524 END IF
525 END DO
526 END DO
527 END DO
528 END DO
529# endif
530!
531!-----------------------------------------------------------------------
532! Time-step adjoint vertical diffusion term.
533!-----------------------------------------------------------------------
534!
535! Compute BASIC STATE reciprocal thickness, 1/Hz.
536!
537 IF (lhsimt) THEN
538 DO k=1,n(ng)
539 DO j=jstrm2,jendp2
540 DO i=istrm2,iendp2
541 ohz(i,j,k)=1.0_r8/hz(i,j,k)
542 END DO
543 END DO
544 END DO
545 ELSE
546 DO k=1,n(ng)
547 DO j=jstr,jend
548 DO i=istr,iend
549 ohz(i,j,k)=1.0_r8/hz(i,j,k)
550 END DO
551 END DO
552 END DO
553 END IF
554!
555! Compute adjoint vertical diffusion term.
556!
557 j_loop2 : DO j=jstr,jend ! start pipelined J-loop
558 DO itrc=1,nt(ng)
559 ltrc=min(nat,itrc)
560
561# ifdef SPLINES_VDIFF
562 IF (.not.((ad_hadvection(itrc,ng)%MPDATA).and. &
563 & (ad_vadvection(itrc,ng)%MPDATA))) THEN
564!
565! Use conservative, parabolic spline reconstruction of BASIC STATE
566! vertical diffusion derivatives. Solve BASIC STATE tridiagonal
567! system.
568!
569 cff1=1.0_r8/6.0_r8
570 DO k=1,n(ng)-1
571 DO i=istr,iend
572 fc(i,k)=cff1*hz(i,j,k )- &
573 & dt(ng)*akt(i,j,k-1,ltrc)*ohz(i,j,k )
574 cf(i,k)=cff1*hz(i,j,k+1)- &
575 & dt(ng)*akt(i,j,k+1,ltrc)*ohz(i,j,k+1)
576 END DO
577 END DO
578 DO i=istr,iend
579 cf(i,0)=0.0_r8
580 dc(i,0)=0.0_r8
581 END DO
582!
583! LU decomposition and forward substitution.
584!
585 cff1=1.0_r8/3.0_r8
586 DO k=1,n(ng)-1
587 DO i=istr,iend
588 bc(i,k)=cff1*(hz(i,j,k)+hz(i,j,k+1))+ &
589 & dt(ng)*akt(i,j,k,ltrc)*(ohz(i,j,k)+ohz(i,j,k+1))
590 cff=1.0_r8/(bc(i,k)-fc(i,k)*cf(i,k-1))
591 cf(i,k)=cff*cf(i,k)
592 dc(i,k)=cff*(t(i,j,k+1,nnew,itrc)-t(i,j,k,nnew,itrc)- &
593 & fc(i,k)*dc(i,k-1))
594 END DO
595 END DO
596!
597! Backward substitution. Save DC for the adjoint code below.
598!
599 DO i=istr,iend
600 dc(i,n(ng))=0.0_r8
601 END DO
602 DO k=n(ng)-1,1,-1
603 DO i=istr,iend
604 dc(i,k)=dc(i,k)-cf(i,k)*dc(i,k+1)
605 END DO
606 END DO
607!
608! Multiply DC by Akt and save it on DC1.
609!
610 DO k=0,n(ng)
611 DO i=istr,iend
612 dc1(i,k)=dc(i,k)*akt(i,j,k,ltrc)
613 END DO
614 END DO
615!
616! Time-step adjoint diffusion term.
617!
618!^ DO k=1,N(ng)
619!^
620 DO k=n(ng),1,-1
621 DO i=istr,iend
622# ifdef DIAGNOSTICS_TS
623!! DiaTwrk(i,j,k,itrc,iTvdif)=DiaTwrk(i,j,k,itrc,iTvdif)+ &
624!! & cff1
625# endif
626!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)+tl_cff1
627!^
628 ad_cff1=ad_cff1+ad_t(i,j,k,nnew,itrc)
629!^ tl_cff1=dt(ng)*(tl_oHz(i,j,k)*(DC(i,k)-DC(i,k-1))+ &
630!^ & oHz(i,j,k)*(tl_DC(i,k)-tl_DC(i,k-1)))
631!^ use DC1 instead
632 adfac=dt(ng)*ad_cff1
633 adfac1=adfac*ohz(i,j,k)
634 ad_dc(i,k-1)=ad_dc(i,k-1)-adfac1
635 ad_dc(i,k )=ad_dc(i,k )+adfac1
636 ad_ohz(i,j,k)=ad_ohz(i,j,k)+(dc1(i,k)-dc1(i,k-1))*adfac
637 ad_cff1=0.0_r8
638!^ tl_DC(i,k)=tl_DC(i,k)*Akt(i,j,k,ltrc)+ &
639!^ & DC(i,k)*tl_Akt(i,j,k,ltrc)
640!^ use DC here
641 ad_dc(i,k)=ad_dc(i,k)*akt(i,j,k,ltrc)
642 ad_akt(i,j,k,ltrc)=ad_akt(i,j,k,ltrc)+dc(i,k)*ad_dc(i,k)
643 END DO
644 END DO
645!
646! Adjoint back substitution
647!
648 DO k=1,n(ng)-1
649 DO i=istr,iend
650!^ tl_DC(i,k)=tl_DC(i,k)-CF(i,k)*tl_DC(i,k+1)
651!^
652 ad_dc(i,k+1)=ad_dc(i,k+1)-cf(i,k)*ad_dc(i,k)
653 END DO
654 END DO
655 DO i=istr,iend
656!^ tl_DC(i,N(ng))=0.0_r8
657!^
658 ad_dc(i,n(ng))=0.0_r8
659 END DO
660!
661! Adjoint LU decomposition and forward substitution.
662!
663 cff1=1.0_r8/3.0_r8
664 DO k=n(ng)-1,1,-1
665 DO i=istr,iend
666 bc(i,k)=cff1*(hz(i,j,k)+hz(i,j,k+1))+ &
667 & dt(ng)*akt(i,j,k,ltrc)*(ohz(i,j,k)+ohz(i,j,k+1))
668 cff=1.0_r8/(bc(i,k)-fc(i,k)*cf(i,k-1))
669!^ tl_DC(i,k)=cff*(tl_t(i,j,k+1,nnew,itrc)- &
670!^ & tl_t(i,j,k ,nnew,itrc)- &
671!^ & (tl_FC(i,k)*DC(i,k-1)+ &
672!^ & tl_BC(i,k)*DC(i,k )+ &
673!^ & tl_CF(i,k)*DC(i,k+1))- &
674!^ & FC(i,k)*tl_DC(i,k-1))
675!^
676 adfac=cff*ad_dc(i,k)
677 ad_dc(i,k-1)=ad_dc(i,k-1)-fc(i,k)*adfac
678 ad_cf(i,k)=ad_cf(i,k)-dc(i,k+1)*adfac
679 ad_bc(i,k)=ad_bc(i,k)-dc(i,k )*adfac
680 ad_fc(i,k)=ad_fc(i,k)-dc(i,k-1)*adfac
681 ad_t(i,j,k ,nnew,itrc)=ad_t(i,j,k ,nnew,itrc)-adfac
682 ad_t(i,j,k+1,nnew,itrc)=ad_t(i,j,k+1,nnew,itrc)+adfac
683 ad_dc(i,k)=0.0_r8
684!^ tl_BC(i,k)=cff1*(tl_Hz(i,j,k)+tl_Hz(i,j,k+1))+ &
685!^ & dt(ng)*(tl_Akt(i,j,k,ltrc)* &
686!^ & (oHz(i,j,k)+oHz(i,j,k+1))+ &
687!^ & Akt(i,j,k,ltrc)* &
688!^ & (tl_oHz(i,j,k)+tl_oHz(i,j,k+1)))
689!^
690 adfac=cff1*ad_bc(i,k)
691 adfac1=dt(ng)*ad_bc(i,k)
692 adfac2=adfac1*akt(i,j,k,ltrc)
693 ad_ohz(i,j,k )=ad_ohz(i,j,k )+adfac2
694 ad_ohz(i,j,k+1)=ad_ohz(i,j,k+1)+adfac2
695 ad_akt(i,j,k,ltrc)=ad_akt(i,j,k,ltrc)+ &
696 & (ohz(i,j,k)+ohz(i,j,k+1))*adfac1
697 ad_hz(i,j,k )=ad_hz(i,j,k )+adfac
698 ad_hz(i,j,k+1)=ad_hz(i,j,k+1)+adfac
699 ad_bc(i,k)=0.0_r8
700 END DO
701 END DO
702!
703! Use conservative, parabolic spline reconstruction of tangent linear
704! vertical diffusion derivatives. Then, time step vertical diffusion
705! term implicitly.
706!
707! Note that the BASIC STATE "t" must in Tunits when used in the
708! tangent spline routine below, which it does in the present code.
709!
710 DO i=istr,iend
711!^ tl_CF(i,0)=0.0_r8
712!^
713 ad_cf(i,0)=0.0_r8
714!^ tl_DC(i,0)=0.0_r8
715!^
716 ad_dc(i,0)=0.0_r8
717 END DO
718 cff1=1.0_r8/6.0_r8
719 DO k=1,n(ng)-1
720 DO i=istr,iend
721!^ tl_CF(i,k)=cff1*tl_Hz(i,j,k+1)- &
722!^ & dt(ng)*(tl_Akt(i,j,k+1,ltrc)*oHz(i,j,k+1)+ &
723!^ & Akt(i,j,k+1,ltrc)*tl_oHz(i,j,k+1))
724!^
725 adfac=dt(ng)*ad_cf(i,k)
726 ad_ohz(i,j,k+1)=ad_ohz(i,j,k+1)- &
727 & akt(i,j,k+1,ltrc)*adfac
728 ad_akt(i,j,k+1,ltrc)=ad_akt(i,j,k+1,ltrc)- &
729 & ohz(i,j,k+1)*adfac
730 ad_hz(i,j,k+1)=ad_hz(i,j,k+1)+cff1*ad_cf(i,k)
731 ad_cf(i,k)=0.0_r8
732!^ tl_FC(i,k)=cff1*tl_Hz(i,j,k )- &
733!^ & dt(ng)*(tl_Akt(i,j,k-1,ltrc)*oHz(i,j,k )+ &
734!^ & Akt(i,j,k-1,ltrc)*tl_oHz(i,j,k ))
735!^
736 adfac=dt(ng)*ad_fc(i,k)
737 ad_ohz(i,j,k )=ad_ohz(i,j,k )- &
738 & akt(i,j,k-1,ltrc)*adfac
739 ad_akt(i,j,k-1,ltrc)=ad_akt(i,j,k-1,ltrc)- &
740 & ohz(i,j,k )*adfac
741 ad_hz(i,j,k )=ad_hz(i,j,k )+cff1*ad_fc(i,k)
742 ad_fc(i,k)=0.0_r8
743 END DO
744 END DO
745 ELSE
746# endif
747!
748! Compute off-diagonal BASIC STATE coefficients FC [lambda*dt*Akt/Hz]
749! for the implicit vertical diffusion terms at future time step,
750! located at horizontal RHO-points and vertical W-points.
751! Also set FC at the top and bottom levels.
752!
753 cff=-dt(ng)*lambda
754 DO k=1,n(ng)-1
755 DO i=istr,iend
756 cff1=1.0_r8/(z_r(i,j,k+1)-z_r(i,j,k))
757 fc(i,k)=cff*cff1*akt(i,j,k,ltrc)
758 END DO
759 END DO
760 DO i=istr,iend
761 fc(i,0)=0.0_r8
762 fc(i,n(ng))=0.0_r8
763 END DO
764!
765! Compute diagonal matrix coefficients BC.
766!
767 DO k=1,n(ng)
768 DO i=istr,iend
769 bc(i,k)=hz(i,j,k)-fc(i,k)-fc(i,k-1)
770 END DO
771 END DO
772!
773! Compute new solution by back substitution.
774! (DC is a tangent linear variable here).
775!
776 DO i=istr,iend
777 cff=1.0_r8/bc(i,1)
778 cf(i,1)=cff*fc(i,1)
779 END DO
780 DO k=2,n(ng)-1
781 DO i=istr,iend
782 cff=1.0_r8/(bc(i,k)-fc(i,k-1)*cf(i,k-1))
783 cf(i,k)=cff*fc(i,k)
784 END DO
785 END DO
786!^ DO k=N(ng)-1,1,-1
787!^
788 DO k=1,n(ng)-1
789 DO i=istr,iend
790# ifdef DIAGNOSTICS_TS
791!! DiaTwrk(i,j,k,itrc,iTvdif)=DiaTwrk(i,j,k,itrc,iTvdif)+ &
792!! & t(i,j,k,nnew,itrc)-cff1
793# endif
794!^ tl_t(i,j,k,nnew,itrc)=DC(i,k)
795!^
796 ad_dc(i,k)=ad_dc(i,k)+ad_t(i,j,k,nnew,itrc)
797 ad_t(i,j,k,nnew,itrc)=0.0_r8
798!^ DC(i,k)=DC(i,k)-CF(i,k)*DC(i,k+1)
799!^
800 ad_dc(i,k+1)=-cf(i,k)*ad_dc(i,k)
801# ifdef DIAGNOSTICS_TS
802!! cff1=t(i,j,k,nnew,itrc)*oHz(i,j,k)
803# endif
804 END DO
805 END DO
806 DO i=istr,iend
807# ifdef DIAGNOSTICS_TS
808!! DiaTwrk(i,j,N(ng),itrc,iTvdif)= &
809!! & DiaTwrk(i,j,N(ng),itrc,iTvdif)+ &
810!! & t(i,j,N(ng),nnew,itrc)-cff1
811# endif
812!^ tl_t(i,j,N(ng),nnew,itrc)=DC(i,N(ng))
813!^
814 ad_dc(i,n(ng))=ad_dc(i,n(ng))+ad_t(i,j,n(ng),nnew,itrc)
815 ad_t(i,j,n(ng),nnew,itrc)=0.0_r8
816!^ DC(i,N(ng))=(DC(i,N(ng))-FC(i,N(ng)-1)*DC(i,N(ng)-1))/ &
817!^ & (BC(i,N(ng))-FC(i,N(ng)-1)*CF(i,N(ng)-1))
818!^
819 adfac=ad_dc(i,n(ng))/ &
820 & (bc(i,n(ng))-fc(i,n(ng)-1)*cf(i,n(ng)-1))
821 ad_dc(i,n(ng)-1)=ad_dc(i,n(ng)-1)-fc(i,n(ng)-1)*adfac
822 ad_dc(i,n(ng) )=adfac
823 END DO
824!
825! Solve the adjoint tridiagonal system.
826! (DC is a tangent linear variable here).
827!
828 DO k=n(ng)-1,2,-1
829 DO i=istr,iend
830 cff=1.0_r8/(bc(i,k)-fc(i,k-1)*cf(i,k-1))
831!^ DC(i,k)=cff*(DC(i,k)-FC(i,k-1)*DC(i,k-1))
832!^
833 adfac=cff*ad_dc(i,k)
834 ad_dc(i,k-1)=ad_dc(i,k-1)-fc(i,k-1)*adfac
835 ad_dc(i,k )=adfac
836 END DO
837 END DO
838 DO i=istr,iend
839 cff=1.0_r8/bc(i,1)
840!^ DC(i,1)=cff*DC(i,1)
841!^
842 ad_dc(i,1)=cff*ad_dc(i,1)
843 END DO
844!
845 DO i=istr,iend
846!^ DC(i,N(ng))=tl_t(i,j,N(ng),nnew,itrc)- &
847!^ & (tl_FC(i,N(ng)-1)*t(i,j,N(ng)-1,nnew,itrc)+ &
848!^ & tl_BC(i,N(ng) )*t(i,j,N(ng) ,nnew,itrc))
849!^
850 ad_bc(i,n(ng) )=-t(i,j,n(ng) ,nnew,itrc)*ad_dc(i,n(ng))
851 ad_fc(i,n(ng)-1)=-t(i,j,n(ng)-1,nnew,itrc)*ad_dc(i,n(ng))
852 ad_t(i,j,n(ng),nnew,itrc)=ad_dc(i,n(ng))
853 ad_dc(i,n(ng))=0.0_r8
854!^ DC(i,1)=tl_t(i,j,1,nnew,itrc)- &
855!^ & (tl_BC(i,1)*t(i,j,1,nnew,itrc)+ &
856!^ & tl_FC(i,1)*t(i,j,2,nnew,itrc))
857!^
858 ad_fc(i,1)=-t(i,j,2,nnew,itrc)*ad_dc(i,1)
859 ad_bc(i,1)=-t(i,j,1,nnew,itrc)*ad_dc(i,1)
860 ad_t(i,j,1,nnew,itrc)=ad_dc(i,1)
861 ad_dc(i,1)=0.0_r8
862 END DO
863 DO k=2,n(ng)-1
864 DO i=istr,iend
865!^ DC(i,k)=tl_t(i,j,k,nnew,itrc)- &
866!^ & (tl_FC(i,k-1)*t(i,j,k-1,nnew,itrc)+ &
867!^ & tl_BC(i,k )*t(i,j,k ,nnew,itrc)+ &
868!^ & tl_FC(i,k )*t(i,j,k+1,nnew,itrc))
869!^
870 ad_fc(i,k-1)=ad_fc(i,k-1)- &
871 & t(i,j,k-1,nnew,itrc)*ad_dc(i,k)
872 ad_fc(i,k )=ad_fc(i,k )- &
873 & t(i,j,k+1,nnew,itrc)*ad_dc(i,k)
874 ad_bc(i,k)=ad_bc(i,k)- &
875 & t(i,j,k ,nnew,itrc)*ad_dc(i,k)
876 ad_t(i,j,k,nnew,itrc)=ad_t(i,j,k,nnew,itrc)+ad_dc(i,k)
877 ad_dc(i,k)=0.0_r8
878 END DO
879 END DO
880!
881! Compute diagonal matrix coefficients BC.
882!
883 DO k=1,n(ng)
884 DO i=istr,iend
885!^ tl_BC(i,k)=tl_Hz(i,j,k)-tl_FC(i,k)-tl_FC(i,k-1)
886!^
887 ad_fc(i,k-1)=ad_fc(i,k-1)-ad_bc(i,k)
888 ad_fc(i,k )=ad_fc(i,k )-ad_bc(i,k)
889 ad_hz(i,j,k)=ad_hz(i,j,k)+ad_bc(i,k)
890 ad_bc(i,k)=0.0_r8
891 END DO
892 END DO
893!
894! Compute off-diagonal coefficients FC [lambda*dt*Akt/Hz] for the
895! implicit vertical diffusion terms at future time step, located
896! at horizontal RHO-points and vertical W-points.
897! Also set FC at the top and bottom levels.
898!
899! NOTE: The original code solves the tridiagonal system A*t=r where
900! A is a matrix and t and r are vectors. We need to solve the
901! tangent linear C of this system which is A*tl_t+tl_A*t=tl_r.
902! Here, tl_A*t and tl_r are known, so we must solve for tl_t
903! by inverting A*tl_t=tl_r-tl_A*t.
904!
905 DO i=istr,iend
906!^ tl_FC(i,N(ng))=0.0_r8
907!^
908 ad_fc(i,n(ng))=0.0_r8
909!^ tl_FC(i,0)=0.0_r8
910!^
911 ad_fc(i,0)=0.0_r8
912 END DO
913 cff=-dt(ng)*lambda
914 DO k=1,n(ng)-1
915 DO i=istr,iend
916 cff1=1.0_r8/(z_r(i,j,k+1)-z_r(i,j,k))
917!^ tl_FC(i,k)=cff*(tl_cff1*Akt(i,j,k,ltrc)+ &
918!^ & cff1*tl_Akt(i,j,k,ltrc))
919!^
920 adfac=cff*ad_fc(i,k)
921 ad_akt(i,j,k,ltrc)=ad_akt(i,j,k,ltrc)+cff1*adfac
922 ad_cff1=ad_cff1+akt(i,j,k,ltrc)*adfac
923 ad_fc(i,k)=0.0_r8
924!^ tl_cff1=-cff1*cff1*(tl_z_r(i,j,k+1)-tl_z_r(i,j,k))
925!^
926 adfac=-cff1*cff1*ad_cff1
927 ad_z_r(i,j,k )=ad_z_r(i,j,k )-adfac
928 ad_z_r(i,j,k+1)=ad_z_r(i,j,k+1)+adfac
929 ad_cff1=0.0_r8
930 END DO
931 END DO
932# ifdef SPLINES_VDIFF
933 END IF
934# endif
935 END DO
936 END DO j_loop2
937!
938!-----------------------------------------------------------------------
939! Adjoint of add tracer divergence due to cell-centered (LwSrc) point
940! sources.
941!-----------------------------------------------------------------------
942!
943! When LTracerSrc is .true. the inflowing concentration is Tsrc.
944! When LtracerSrc is .false. we add tracer mass to compensate for the
945! added volume to keep the receiving cell concentration unchanged.
946! J. Levin (Jupiter Intelligence Inc.) and J. Wilkin
947!
948! Dsrc(is) = 2, flow across grid cell w-face (positive or negative)
949!
950 IF (lwsrc(ng)) THEN
951 DO itrc=1,nt(ng)
952 IF (.not.((ad_hadvection(itrc,ng)%MPDATA).and. &
953 & (ad_vadvection(itrc,ng)%MPDATA))) THEN
954 DO is=1,nsrc(ng)
955 IF (int(sources(ng)%Dsrc(is)).eq.2) THEN
956 isrc=sources(ng)%Isrc(is)
957 jsrc=sources(ng)%Jsrc(is)
958 IF (((istr.le.isrc).and.(isrc.le.iend+1)).and. &
959 & ((jstr.le.jsrc).and.(jsrc.le.jend+1))) THEN
960 DO k=1,n(ng)
961 cff=dt(ng)*pm(isrc,jsrc)*pn(isrc,jsrc)
962# ifdef SPLINES_VDIFF
963 cff=cff*ohz(isrc,jsrc,k)
964# endif
965 IF (ltracersrc(itrc,ng)) THEN
966 cff3=sources(ng)%Tsrc(is,k,itrc)
967 ELSE
968 cff3=t(isrc,jsrc,k,3,itrc)
969 END IF
970# ifdef SPLINES_VDIFF
971!^ tl_t(Isrc,Jsrc,k,nnew,itrc)= &
972!^ & tl_t(Isrc,Jsrc,k,nnew,itrc)+ &
973!^ & cff*(SOURCES(ng)%tl_Qsrc(is,k)* &
974!^ & cff3+ &
975!^ & SOURCES(ng)%Qsrc(is,k)* &
976!^ & tl_cff3)+ &
977!^ & tl_cff*SOURCES(ng)%Qsrc(is,k)* &
978!^ & cff3
979!^
980 adfac=cff*ad_t(isrc,jsrc,k,nnew,itrc)
981 sources(ng)%ad_Qsrc(is,k)=sources(ng)%ad_Qsrc(is,k)+&
982 & cff3*adfac
983 ad_cff3=ad_cff3+ &
984 & sources(ng)%Qsrc(is,k)*adfac
985 ad_cff=ad_cff+ &
986 & sources(ng)%Qsrc(is,k)*cff3* &
987 & ad_t(isrc,jsrc,k,nnew,itrc)
988# else
989!^ tl_t(Isrc,Jsrc,k,nnew,itrc)= &
990!^ & tl_t(Isrc,Jsrc,k,nnew,itrc)+ &
991!^ & cff*(SOURCES(ng)%tl_Qsrc(is,k)* &
992!^ & cff3+ &
993!^ & SOURCES(ng)%Qsrc(is,k)* &
994!^ & tl_cff3)
995!^
996 adfac=cff*ad_t(isrc,jsrc,k,nnew,itrc)
997 sources(ng)%ad_Qsrc(is,k)=sources(ng)%ad_Qsrc(is,k)+&
998 & cff3*adfac
999 ad_cff3=ad_cff3+ &
1000 & sources(ng)%Qsrc(is,k)*adfac
1001# endif
1002 IF (ltracersrc(itrc,ng)) THEN
1003!^ tl_cff3=SOURCES(ng)%tl_Tsrc(is,k,itrc)
1004!^
1005 sources(ng)%ad_Tsrc(is,k,itrc)= &
1006 & sources(ng)%ad_Tsrc(is,k,itrc)+ &
1007 ad_cff3
1008 ELSE
1009!^ tl_cff3=tl_t(Isrc,Jsrc,k,3,itrc)
1010!^
1011 ad_t(isrc,jsrc,k,3,itrc)=ad_t(isrc,jsrc,k,3,itrc)+&
1012 & ad_cff3
1013 END IF
1014 ad_cff3=0.0_r8
1015# ifdef SPLINES_VDIFF
1016!^ tl_cff=cff*tl_oHz(Isrc,Jsrc,k)
1017!^
1018 ad_ohz(isrc,jsrc,k)=ad_ohz(isrc,jsrc,k)+ &
1019 & cff*ad_cff
1020 ad_cff=0.0_r8
1021# endif
1022 END DO
1023 END IF
1024 END IF
1025 END DO
1026 END IF
1027 END DO
1028 END IF
1029!
1030!-----------------------------------------------------------------------
1031! Time-step adjoint vertical advection term.
1032!-----------------------------------------------------------------------
1033!
1034 t_loop2 : DO itrc=1,nt(ng)
1035 IF (ad_hadvection(itrc,ng)%MPDATA) THEN
1036 jmint=jstrvm2
1037 jmaxt=jendp2i
1038 ELSE
1039 jmint=jstr
1040 jmaxt=jend
1041 END IF
1042!
1043 j_loop1 : DO j=jmint,jmaxt ! start pipelined J-loop
1044!
1045! Time-step adjoint vertical advection term. Compute first BASIC
1046! STATE vertical advection flux, FC.
1047!
1048 vadv_flux_basic : IF (ad_vadvection(itrc,ng)%SPLINES) THEN
1049!
1050! Build conservative parabolic splines for the BASIC STATE vertical
1051! derivatives "FC" of the tracer. Then, the interfacial "FC" values
1052! are converted to vertical advective flux.
1053!
1054 DO i=istr,iend
1055# ifdef NEUMANN
1056 fc(i,0)=1.5_r8*t(i,j,1,3,itrc)
1057 cf(i,1)=0.5_r8
1058# else
1059 fc(i,0)=2.0_r8*t(i,j,1,3,itrc)
1060 cf(i,1)=1.0_r8
1061# endif
1062 END DO
1063 DO k=1,n(ng)-1
1064 DO i=istr,iend
1065 cff=1.0_r8/(2.0_r8*hz(i,j,k)+ &
1066 & hz(i,j,k+1)*(2.0_r8-cf(i,k)))
1067 cf(i,k+1)=cff*hz(i,j,k)
1068 fc(i,k)=cff*(3.0_r8*(hz(i,j,k )*t(i,j,k+1,3,itrc)+ &
1069 & hz(i,j,k+1)*t(i,j,k ,3,itrc))- &
1070 & hz(i,j,k+1)*fc(i,k-1))
1071 END DO
1072 END DO
1073 DO i=istr,iend
1074# ifdef NEUMANN
1075 fc(i,n(ng))=(3.0_r8*t(i,j,n(ng),3,itrc)-fc(i,n(ng)-1))/ &
1076 & (2.0_r8-cf(i,n(ng)))
1077# else
1078 fc(i,n(ng))=(2.0_r8*t(i,j,n(ng),3,itrc)-fc(i,n(ng)-1))/ &
1079 & (1.0_r8-cf(i,n(ng)))
1080# endif
1081 END DO
1082 DO k=n(ng)-1,0,-1
1083 DO i=istr,iend
1084 fc(i,k)=fc(i,k)-cf(i,k+1)*fc(i,k+1)
1085 fc(i,k+1)=w(i,j,k+1)*fc(i,k+1)
1086 END DO
1087 END DO
1088 DO i=istr,iend
1089 fc(i,n(ng))=0.0_r8
1090 fc(i,0)=0.0_r8
1091 END DO
1092!
1093 ELSE IF (ad_vadvection(itrc,ng)%AKIMA4) THEN
1094!
1095! Fourth-order, BASIC STATE Akima vertical advective flux.
1096!
1097 DO k=1,n(ng)-1
1098 DO i=istr,iend
1099 fc(i,k)=t(i,j,k+1,3,itrc)- &
1100 & t(i,j,k ,3,itrc)
1101 END DO
1102 END DO
1103 DO i=istr,iend
1104 fc(i,0)=fc(i,1)
1105 fc(i,n(ng))=fc(i,n(ng)-1)
1106 END DO
1107 DO k=1,n(ng)
1108 DO i=istr,iend
1109 cff=2.0_r8*fc(i,k)*fc(i,k-1)
1110 IF (cff.gt.eps) THEN
1111 cf(i,k)=cff/(fc(i,k)+fc(i,k-1))
1112 ELSE
1113 cf(i,k)=0.0_r8
1114 END IF
1115 END DO
1116 END DO
1117 cff1=1.0_r8/3.0_r8
1118 DO k=1,n(ng)-1
1119 DO i=istr,iend
1120 fc(i,k)=w(i,j,k)* &
1121 & 0.5_r8*(t(i,j,k ,3,itrc)+ &
1122 & t(i,j,k+1,3,itrc)- &
1123 & cff1*(cf(i,k+1)-cf(i,k)))
1124 END DO
1125 END DO
1126 DO i=istr,iend
1127# ifdef SED_MORPH
1128 fc(i,0)=w(i,j,0)*t(i,j,1,3,itrc)
1129# else
1130 fc(i,0)=0.0_r8
1131# endif
1132 fc(i,n(ng))=0.0_r8
1133 END DO
1134!
1135 ELSE IF (ad_vadvection(itrc,ng)%CENTERED2) THEN
1136!
1137! Second-order, BASIC STATE central differences vertical advective
1138! flux.
1139!
1140 DO k=1,n(ng)-1
1141 DO i=istr,iend
1142 fc(i,k)=w(i,j,k)* &
1143 & 0.5_r8*(t(i,j,k ,3,itrc)+ &
1144 & t(i,j,k+1,3,itrc))
1145 END DO
1146 END DO
1147 DO i=istr,iend
1148# ifdef SED_MORPH
1149 fc(i,0)=w(i,j,0)*t(i,j,1,3,itrc)
1150# else
1151 fc(i,0)=0.0_r8
1152# endif
1153 fc(i,n(ng))=0.0_r8
1154 END DO
1155!
1156 ELSE IF (ad_vadvection(itrc,ng)%MPDATA) THEN
1157!
1158! First_order, BASIC STATE upstream differences vertical advective
1159! flux.
1160!
1161 CONTINUE ! not supported
1162!
1163 ELSE IF (ad_vadvection(itrc,ng)%HSIMT) THEN
1164!
1165! Third High-order Spatial Interpolation at the Middle Temporal level
1166! (HSIMT; Wu and Zhu, 2010) with a Total Variation Diminishing (TVD)
1167! limiter vertical advection flux (Tunits m3/s).
1168!
1169 CONTINUE ! not supported
1170!
1171 ELSE IF ((ad_vadvection(itrc,ng)%CENTERED4).or. &
1172 & (ad_vadvection(itrc,ng)%SPLIT_U3)) THEN
1173!
1174! Fourth-order, BASIC STATE central differences vertical advective
1175! flux. (Not really needed, HGA).
1176!
1177 cff1=0.5_r8
1178 cff2=7.0_r8/12.0_r8
1179 cff3=1.0_r8/12.0_r8
1180 DO k=2,n(ng)-2
1181 DO i=istr,iend
1182 fc(i,k)=w(i,j,k)* &
1183 & (cff2*(t(i,j,k ,3,itrc)+ &
1184 & t(i,j,k+1,3,itrc))- &
1185 & cff3*(t(i,j,k-1,3,itrc)+ &
1186 & t(i,j,k+2,3,itrc)))
1187 END DO
1188 END DO
1189 DO i=istr,iend
1190# ifdef SED_MORPH
1191 fc(i,0)=w(i,j,0)*2.0_r8* &
1192 & (cff2*t(i,j,1,3,itrc)- &
1193 & cff3*t(i,j,2,3,itrc))
1194# else
1195 fc(i,0)=0.0_r8
1196# endif
1197 fc(i,1)=w(i,j,1)* &
1198 & (cff1*t(i,j,1,3,itrc)+ &
1199 & cff2*t(i,j,2,3,itrc)- &
1200 & cff3*t(i,j,3,3,itrc))
1201 fc(i,n(ng)-1)=w(i,j,n(ng)-1)* &
1202 & (cff1*t(i,j,n(ng) ,3,itrc)+ &
1203 & cff2*t(i,j,n(ng)-1,3,itrc)- &
1204 & cff3*t(i,j,n(ng)-2,3,itrc))
1205 fc(i,n(ng))=0.0_r8
1206 END DO
1207 END IF vadv_flux_basic
1208!
1209! Time-step vertical advection term.
1210# ifdef SPLINES_VDIFF
1211! The BASIC STATE "t" used below must be in transport units, but "t"
1212! retrived is in Tunits so we multiply by "Hz".
1213# endif
1214!
1215 vadv_stepping : IF (ad_vadvection(itrc,ng)%MPDATA) THEN
1216 CONTINUE ! not supported
1217 ELSE
1218 DO i=istr,iend
1219 cf(i,0)=dt(ng)*pm(i,j)*pn(i,j)
1220 END DO
1221 DO k=1,n(ng)
1222 DO i=istr,iend
1223 cff1=cf(i,0)*(fc(i,k)-fc(i,k-1))
1224# ifdef DIAGNOSTICS_TS
1225!! DiaTwrk(i,j,k,itrc,iTvadv)=-cff1
1226!! DO idiag=1,NDT
1227!! DiaTwrk(i,j,k,itrc,idiag)=DiaTwrk(i,j,k,itrc,idiag)* &
1228!! & oHz(i,j,k)
1229!! END DO
1230# endif
1231# ifdef SPLINES_VDIFF
1232!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)* &
1233!^ & oHz(i,j,k)+ &
1234!^ & (t(i,j,k,nnew,itrc)*Hz(i,j,k))* &
1235!^ & tl_oHz(i,j,k)
1236!^
1237 ad_ohz(i,j,k)=ad_ohz(i,j,k)+ &
1238 & (t(i,j,k,nnew,itrc)*hz(i,j,k))* &
1239 & ad_t(i,j,k,nnew,itrc)
1240 ad_t(i,j,k,nnew,itrc)=ad_t(i,j,k,nnew,itrc)*ohz(i,j,k)
1241# endif
1242!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)-tl_cff1
1243!^
1244 ad_cff1=ad_cff1-ad_t(i,j,k,nnew,itrc)
1245!^ tl_cff1=CF(i,0)*(tl_FC(i,k)-tl_FC(i,k-1))
1246!^
1247 adfac=cf(i,0)*ad_cff1
1248 ad_fc(i,k-1)=ad_fc(i,k-1)-adfac
1249 ad_fc(i,k )=ad_fc(i,k )+adfac
1250 ad_cff1=0.0_r8
1251 END DO
1252 END DO
1253 END IF vadv_stepping
1254!
1255! Compute adjoint of vertical advection fluxes.
1256!
1257 vadv_flux : IF (ad_vadvection(itrc,ng)%SPLINES) THEN
1258!
1259! Build conservative parabolic splines for the vertical derivatives
1260! "FC" of the tracer. Then, the interfacial "FC" values are
1261! converted to vertical advective flux.
1262!
1263 DO i=istr,iend
1264# ifdef NEUMANN
1265 fc(i,0)=1.5_r8*t(i,j,1,3,itrc)
1266 cf(i,1)=0.5_r8
1267# else
1268 fc(i,0)=2.0_r8*t(i,j,1,3,itrc)
1269 cf(i,1)=1.0_r8
1270# endif
1271 END DO
1272 DO k=1,n(ng)-1
1273 DO i=istr,iend
1274 cff=1.0_r8/(2.0_r8*hz(i,j,k)+ &
1275 & hz(i,j,k+1)*(2.0_r8-cf(i,k)))
1276 cf(i,k+1)=cff*hz(i,j,k)
1277 fc(i,k)=cff*(3.0_r8*(hz(i,j,k )*t(i,j,k+1,3,itrc)+ &
1278 & hz(i,j,k+1)*t(i,j,k ,3,itrc))- &
1279 & hz(i,j,k+1)*fc(i,k-1))
1280 END DO
1281 END DO
1282 DO i=istr,iend
1283# ifdef NEUMANN
1284 fc(i,n(ng))=(3.0_r8*t(i,j,n(ng),3,itrc)-fc(i,n(ng)-1))/ &
1285 & (2.0_r8-cf(i,n(ng)))
1286# else
1287 fc(i,n(ng))=(2.0_r8*t(i,j,n(ng),3,itrc)-fc(i,n(ng)-1))/ &
1288 & (1.0_r8-cf(i,n(ng)))
1289# endif
1290 END DO
1291 DO k=n(ng)-1,0,-1
1292 DO i=istr,iend
1293 fc(i,k)=fc(i,k)-cf(i,k+1)*fc(i,k+1)
1294 END DO
1295 END DO
1296!
1297! Now the adjoint splines code.
1298!
1299 DO i=istr,iend
1300!^ tl_FC(i,N(ng))=0.0_r8
1301!^
1302 ad_fc(i,n(ng))=0.0_r8
1303!^ tl_FC(i,0)=0.0_r8
1304!^
1305 ad_fc(i,0)=0.0_r8
1306 END DO
1307!
1308! Adjoint back substitution.
1309!
1310 DO k=0,n(ng)-1
1311 DO i=istr,iend
1312!^ tl_FC(i,k+1)=tl_W(i,j,k+1)*FC(i,k+1)+ &
1313!^ & W(i,j,k+1)*tl_FC(i,k+1)
1314!^
1315 ad_w(i,j,k+1)=ad_w(i,j,k+1)+fc(i,k+1)*ad_fc(i,k+1)
1316 ad_fc(i,k+1)=w(i,j,k+1)*ad_fc(i,k+1)
1317!^ tl_FC(i,k)=tl_FC(i,k)-CF(i,k+1)*tl_FC(i,k+1)
1318!^
1319 ad_fc(i,k+1)=ad_fc(i,k+1)-cf(i,k+1)*ad_fc(i,k)
1320 END DO
1321 END DO
1322 DO i=istr,iend
1323# ifdef NEUMANN
1324!^ tl_FC(i,N(ng))=(3.0_r8*tl_t(i,j,N(ng),3,itrc)- &
1325!^ & tl_FC(i,N(ng)-1))/ &
1326!^ & (2.0_r8-CF(i,N(ng)))
1327!^
1328 adfac=ad_fc(i,n(ng))/(2.0_r8-cf(i,n(ng)))
1329 ad_t(i,j,n(ng),3,itrc)=ad_t(i,j,n(ng),3,itrc)+3.0_r8*adfac
1330 ad_fc(i,n(ng)-1)=ad_fc(i,n(ng)-1)-adfac
1331 ad_fc(i,n(ng))=0.0_r8
1332# else
1333!^ tl_FC(i,N(ng))=(2.0_r8*tl_t(i,j,N(ng),3,itrc)- &
1334!^ & tl_FC(i,N(ng)-1))/ &
1335!^ & (1.0_r8-CF(i,N(ng)))
1336!^
1337 adfac=ad_fc(i,n(ng))/(1.0_r8-cf(i,n(ng)))
1338 ad_t(i,j,n(ng),3,itrc)=ad_t(i,j,n(ng),3,itrc)+2.0_r8*adfac
1339 ad_fc(i,n(ng)-1)=ad_fc(i,n(ng)-1)-adfac
1340 ad_fc(i,n(ng))=0.0
1341# endif
1342 END DO
1343 DO k=n(ng)-1,1,-1
1344 DO i=istr,iend
1345 cff=1.0_r8/(2.0_r8*hz(i,j,k)+ &
1346 & hz(i,j,k+1)*(2.0_r8-cf(i,k)))
1347!^ tl_FC(i,k)=cff* &
1348!^ & (3.0_r8*(Hz(i,j,k )*tl_t(i,j,k+1,3,itrc)+ &
1349!^ & Hz(i,j,k+1)*tl_t(i,j,k ,3,itrc)+ &
1350!^ & tl_Hz(i,j,k )*t(i,j,k+1,3,itrc)+ &
1351!^ & tl_Hz(i,j,k+1)*t(i,j,k ,3,itrc))- &
1352!^ & (tl_Hz(i,j,k+1)*FC(i,k-1)+ &
1353!^ & 2.0_r8*(tl_Hz(i,j,k )+ &
1354!^ & tl_Hz(i,j,k+1))*FC(i,k)+ &
1355!^ & tl_Hz(i,j,k )*FC(i,k+1))- &
1356!^ & Hz(i,j,k+1)*tl_FC(i,k-1))
1357!^
1358 adfac=cff*ad_fc(i,k)
1359 adfac1=3.0_r8*adfac
1360 adfac2=2.0_r8*adfac
1361 ad_t(i,j,k ,3,itrc)=ad_t(i,j,k ,3,itrc)+ &
1362 & hz(i,j,k+1)*adfac1
1363 ad_t(i,j,k+1,3,itrc)=ad_t(i,j,k+1,3,itrc)+ &
1364 & hz(i,j,k )*adfac1
1365 ad_hz(i,j,k )=ad_hz(i,j,k )+ &
1366 & t(i,j,k+1,3,itrc)*adfac1- &
1367 & fc(i,k )*adfac2- &
1368 & fc(i,k+1)*adfac
1369 ad_hz(i,j,k+1)=ad_hz(i,j,k+1)+ &
1370 & t(i,j,k ,3,itrc)*adfac1- &
1371 & fc(i,k-1)*adfac- &
1372 & fc(i,k )*adfac2
1373 ad_fc(i,k-1)=ad_fc(i,k-1)-hz(i,j,k+1)*adfac
1374 ad_fc(i,k)=0.0_r8
1375 END DO
1376 END DO
1377 DO i=istr,iend
1378# ifdef NEUMANN
1379!^ tl_FC(i,0)=1.5_r8*tl_t(i,j,1,3,itrc)
1380!^
1381 ad_t(i,j,1,3,itrc)=ad_t(i,j,1,3,itrc)+1.5_r8*ad_fc(i,0)
1382 ad_fc(i,0)=0.0_r8
1383# else
1384!^ tl_FC(i,0)=2.0_r8*tl_t(i,j,1,3,itrc)
1385!^
1386 ad_t(i,j,1,3,itrc)=ad_t(i,j,1,3,itrc)+2.0_r8*ad_fc(i,0)
1387 ad_fc(i,0)=0.0_r8
1388# endif
1389 END DO
1390!
1391 ELSE IF (ad_vadvection(itrc,ng)%AKIMA4) THEN
1392!
1393! Fourth-order, Akima adjoint vertical advective flux.
1394!
1395 DO k=1,n(ng)-1
1396 DO i=istr,iend
1397 fc(i,k)=t(i,j,k+1,3,itrc)- &
1398 & t(i,j,k ,3,itrc)
1399 END DO
1400 END DO
1401 DO i=istr,iend
1402 fc(i,0)=fc(i,1)
1403 fc(i,n(ng))=fc(i,n(ng)-1)
1404 END DO
1405 DO k=1,n(ng)
1406 DO i=istr,iend
1407 cff=2.0_r8*fc(i,k)*fc(i,k-1)
1408 IF (cff.gt.eps) THEN
1409 cf(i,k)=cff/(fc(i,k)+fc(i,k-1))
1410 ELSE
1411 cf(i,k)=0.0_r8
1412 END IF
1413 END DO
1414 END DO
1415 DO i=istr,iend
1416!^ tl_FC(i,N(ng))=0.0_r8
1417!^
1418 ad_fc(i,n(ng))=0.0_r8
1419# ifdef SED_MORPH
1420!^ tl_FC(i,0)=tl_W(i,j,0)*t(i,j,1,3,itrc)+ &
1421!^ & W(i,j,0)*tl_t(i,j,1,3,itrc)
1422!^
1423 ad_t(i,j,1,3,itrc)=ad_t(i,j,1,3,itrc)+w(i,j,0)*ad_fc(i,0)
1424 ad_w(i,j,0)=ad_w(i,j,0)+t(i,j,1,3,itrc)*ad_fc(i,0)
1425 ad_fc(i,0)=0.0_r8
1426# else
1427!^ tl_FC(i,0)=0.0_r8
1428!^
1429 ad_fc(i,0)=0.0_r8
1430# endif
1431 END DO
1432 cff1=1.0_r8/3.0_r8
1433 DO k=1,n(ng)-1
1434 DO i=istr,iend
1435!^ tl_FC(i,k)=0.5_r8* &
1436!^ & (tl_W(i,j,k)* &
1437!^ & (t(i,j,k ,3,itrc)+ &
1438!^ & t(i,j,k+1,3,itrc)- &
1439!^ & cff1*(CF(i,k+1)-CF(i,k)))+ &
1440!^ & W(i,j,k)* &
1441!^ & (tl_t(i,j,k ,3,itrc)+ &
1442!^ & tl_t(i,j,k+1,3,itrc)- &
1443!^ & cff1*(tl_CF(i,k+1)-tl_CF(i,k))))
1444!^
1445 adfac=0.5_r8*ad_fc(i,k)
1446 adfac1=adfac*w(i,j,k)
1447 adfac2=adfac1*cff1
1448 ad_cf(i,k )=ad_cf(i,k )+adfac2
1449 ad_cf(i,k+1)=ad_cf(i,k+1)-adfac2
1450 ad_t(i,j,k ,3,itrc)=ad_t(i,j,k ,3,itrc)+adfac1
1451 ad_t(i,j,k+1,3,itrc)=ad_t(i,j,k+1,3,itrc)+adfac1
1452 ad_w(i,j,k)=ad_w(i,j,k)+ &
1453 & (t(i,j,k ,3,itrc)+ &
1454 & t(i,j,k+1,3,itrc)- &
1455 & cff1*(cf(i,k+1)-cf(i,k)))*adfac
1456 ad_fc(i,k)=0.0_r8
1457 END DO
1458 END DO
1459 DO k=1,n(ng)
1460 DO i=istr,iend
1461 cff=2.0_r8*fc(i,k)*fc(i,k-1)
1462 IF (cff.gt.eps) THEN
1463!^ tl_CF(i,k)=((FC(i,k)+FC(i,k-1))*tl_cff- &
1464!^ & cff*(tl_FC(i,k)+tl_FC(i,k-1)))/ &
1465!^ & ((FC(i,k)+FC(i,k-1))*(FC(i,k)+FC(i,k-1)))
1466!^
1467 adfac=ad_cf(i,k)/ &
1468 & ((fc(i,k)+fc(i,k-1))*(fc(i,k)+fc(i,k-1)))
1469 adfac1=adfac*cff
1470 ad_fc(i,k-1)=ad_fc(i,k-1)-adfac1
1471 ad_fc(i,k )=ad_fc(i,k )-adfac1
1472 ad_cff=ad_cff+(fc(i,k)+fc(i,k-1))*adfac
1473 ad_cf(i,k)=0.0_r8
1474 ELSE
1475!^ tl_CF(i,k)=0.0_r8
1476!^
1477 ad_cf(i,k)=0.0_r8
1478 END IF
1479!^ tl_cff=2.0_r8*(tl_FC(i,k)*FC(i,k-1)+ &
1480!^ & FC(i,k)*tl_FC(i,k-1))
1481!^
1482 adfac=2.0_r8*ad_cff
1483 ad_fc(i,k-1)=ad_fc(i,k-1)+fc(i,k )*adfac
1484 ad_fc(i,k )=ad_fc(i,k )+fc(i,k-1)*adfac
1485 ad_cff=0.0_r8
1486 END DO
1487 END DO
1488 DO i=istr,iend
1489!^ tl_FC(i,N(ng))=tl_FC(i,N(ng)-1)
1490!^
1491 ad_fc(i,n(ng)-1)=ad_fc(i,n(ng)-1)+ad_fc(i,n(ng))
1492 ad_fc(i,n(ng))=0.0_r8
1493!^ tl_FC(i,0)=tl_FC(i,1)
1494!^
1495 ad_fc(i,1)=ad_fc(i,1)+ad_fc(i,0)
1496 ad_fc(i,0)=0.0_r8
1497 END DO
1498 DO k=1,n(ng)-1
1499 DO i=istr,iend
1500!^ tl_FC(i,k)=tl_t(i,j,k+1,3,itrc)- &
1501!^ & tl_t(i,j,k ,3,itrc)
1502!^
1503 ad_t(i,j,k ,3,itrc)=ad_t(i,j,k ,3,itrc)-ad_fc(i,k)
1504 ad_t(i,j,k+1,3,itrc)=ad_t(i,j,k+1,3,itrc)+ad_fc(i,k)
1505 ad_fc(i,k)=0.0_r8
1506 END DO
1507 END DO
1508!
1509 ELSE IF (ad_vadvection(itrc,ng)%CENTERED2) THEN
1510!
1511! Second-order, central differences adjoint vertical advective flux.
1512!
1513 DO i=istr,iend
1514!^ tl_FC(i,N(ng))=0.0_r8
1515!^
1516 ad_fc(i,n(ng))=0.0_r8
1517# ifdef SED_MORPH
1518!^ tl_FC(i,0)=tl_W(i,j,0)*t(i,j,1,3,itrc)+ &
1519!^ & W(i,j,0)*tl_t(i,j,1,3,itrc)
1520!^
1521 ad_t(i,j,1,3,itrc)=ad_t(i,j,1,3,itrc)+w(i,j,0)*ad_fc(i,0)
1522 ad_w(i,j,0)=ad_w(i,j,0)+t(i,j,1,3,itrc)*ad_fc(i,0)
1523 ad_fc(i,0)=0.0_r8
1524# else
1525!^ tl_FC(i,0)=0.0_r8
1526!^
1527 ad_fc(i,0)=0.0_r8
1528# endif
1529 END DO
1530 DO k=1,n(ng)-1
1531 DO i=istr,iend
1532!^ tl_FC(i,k)=0.5_r8* &
1533!^ & (tl_W(i,j,k)* &
1534!^ & (t(i,j,k ,3,itrc)+ &
1535!^ & t(i,j,k+1,3,itrc))+ &
1536!^ & W(i,j,k)* &
1537!^ & (tl_t(i,j,k ,3,itrc)+ &
1538!^ & tl_t(i,j,k+1,3,itrc)))
1539!^
1540 adfac=0.5_r8*ad_fc(i,k)
1541 adfac1=adfac*w(i,j,k)
1542 ad_t(i,j,k ,3,itrc)=ad_t(i,j,k ,3,itrc)+adfac1
1543 ad_t(i,j,k+1,3,itrc)=ad_t(i,j,k+1,3,itrc)+adfac1
1544 ad_w(i,j,k)=ad_w(i,j,k)+ &
1545 & (t(i,j,k ,3,itrc)+ &
1546 & t(i,j,k+1,3,itrc))*adfac
1547 ad_fc(i,k)=0.0_r8
1548 END DO
1549 END DO
1550!
1551 ELSE IF (ad_vadvection(itrc,ng)%MPDATA) THEN
1552!
1553! First_order, upstream differences vertical advective flux.
1554!
1555 CONTINUE ! not supported
1556!
1557 ELSE IF (ad_vadvection(itrc,ng)%HSIMT) THEN
1558!
1559! Third High-order Spatial Interpolation at the Middle Temporal level
1560! (HSIMT; Wu and Zhu, 2010) with a Total Variation Diminishing (TVD)
1561! limiter vertical advection flux (Tunits m3/s).
1562!
1563 CONTINUE ! not supported
1564!
1565 ELSE IF ((ad_vadvection(itrc,ng)%CENTERED4).or. &
1566 & (ad_vadvection(itrc,ng)%SPLIT_U3)) THEN
1567!
1568! Fourth-order, central differences adjoint vertical advective flux.
1569!
1570 cff1=0.5_r8
1571 cff2=7.0_r8/12.0_r8
1572 cff3=1.0_r8/12.0_r8
1573 DO i=istr,iend
1574!^ tl_FC(i,N(ng))=0.0_r8
1575!^
1576 ad_fc(i,n(ng))=0.0_r8
1577!^ tl_FC(i,N(ng)-1)=tl_W(i,j,N(ng)-1)* &
1578!^ & (cff1*t(i,j,N(ng) ,3,itrc)+ &
1579!^ & cff2*t(i,j,N(ng)-1,3,itrc)- &
1580!^ & cff3*t(i,j,N(ng)-2,3,itrc))+ &
1581!^ & W(i,j,N(ng)-1)* &
1582!^ & (cff1*tl_t(i,j,N(ng) ,3,itrc)+ &
1583!^ & cff2*tl_t(i,j,N(ng)-1,3,itrc)- &
1584!^ & cff3*tl_t(i,j,N(ng)-2,3,itrc))
1585!^
1586 adfac=w(i,j,n(ng)-1)*ad_fc(i,n(ng)-1)
1587 ad_w(i,j,n(ng)-1)=ad_w(i,j,n(ng)-1)+ &
1588 & (cff1*t(i,j,n(ng) ,3,itrc)+ &
1589 & cff2*t(i,j,n(ng)-1,3,itrc)- &
1590 & cff3*t(i,j,n(ng)-2,3,itrc))* &
1591 & ad_fc(i,n(ng)-1)
1592 ad_t(i,j,n(ng)-2,3,itrc)=ad_t(i,j,n(ng)-2,3,itrc)- &
1593 & cff3*adfac
1594 ad_t(i,j,n(ng)-1,3,itrc)=ad_t(i,j,n(ng)-1,3,itrc)+ &
1595 & cff2*adfac
1596 ad_t(i,j,n(ng) ,3,itrc)=ad_t(i,j,n(ng) ,3,itrc)+ &
1597 & cff1*adfac
1598 ad_fc(i,n(ng)-1)=0.0_r8
1599!^ tl_FC(i,1)=tl_W(i,j,1)* &
1600!^ & (cff1*t(i,j,1,3,itrc)+ &
1601!^ & cff2*t(i,j,2,3,itrc)- &
1602!^ & cff3*t(i,j,3,3,itrc))+ &
1603!^ & W(i,j,1)* &
1604!^ & (cff1*tl_t(i,j,1,3,itrc)+ &
1605!^ & cff2*tl_t(i,j,2,3,itrc)- &
1606!^ & cff3*tl_t(i,j,3,3,itrc))
1607!^
1608 adfac=w(i,j,1)*ad_fc(i,1)
1609 ad_w(i,j,1)=ad_w(i,j,1)+ &
1610 & (cff1*t(i,j,1,3,itrc)+ &
1611 & cff2*t(i,j,2,3,itrc)- &
1612 & cff3*t(i,j,3,3,itrc))*ad_fc(i,1)
1613 ad_t(i,j,1,3,itrc)=ad_t(i,j,1,3,itrc)+cff1*adfac
1614 ad_t(i,j,2,3,itrc)=ad_t(i,j,2,3,itrc)+cff2*adfac
1615 ad_t(i,j,3,3,itrc)=ad_t(i,j,3,3,itrc)-cff3*adfac
1616 ad_fc(i,1)=0.0_r8
1617# ifdef SED_MORPH
1618!^ tl_FC(i,0)=2.0_r8* &
1619!^ & (tl_W(i,j,0)* &
1620!^ & (cff2*t(i,j,1,3,itrc)- &
1621!^ & cff3*t(i,j,2,3,itrc))+ &
1622!^ & W(i,j,0)* &
1623!^ & (cff2*tl_t(i,j,1,3,itrc)- &
1624!^ & cff3*tl_t(i,j,2,3,itrc)))
1625!^
1626 adfac=2.0_r8*ad_fc(i,0)
1627 adfac1=adfac*w(i,j,0)
1628 ad_t(i,j,1,3,itrc)=ad_t(i,j,1,3,itrc)+cff2*adfac1
1629 ad_t(i,j,2,3,itrc)=ad_t(i,j,2,3,itrc)-cff3*adfac1
1630 ad_w(i,j,0)=ad_w(i,j,0)+ &
1631 & (cff2*t(i,j,1,3,itrc)- &
1632 & cff3*t(i,j,2,3,itrc))*adfac
1633 ad_fc(i,0)=0.0_r8
1634# else
1635!^ tl_FC(i,0)=0.0_r8
1636!^
1637 ad_fc(i,0)=0.0_r8
1638# endif
1639 END DO
1640 DO k=2,n(ng)-2
1641 DO i=istr,iend
1642!^ tl_FC(i,k)=tl_W(i,j,k)* &
1643!^ & (cff2*(t(i,j,k ,3,itrc)+ &
1644!^ & t(i,j,k+1,3,itrc))- &
1645!^ & cff3*(t(i,j,k-1,3,itrc)+ &
1646!^ & t(i,j,k+2,3,itrc)))+ &
1647!^ & W(i,j,k)* &
1648!^ & (cff2*(tl_t(i,j,k ,3,itrc)+ &
1649!^ & tl_t(i,j,k+1,3,itrc))- &
1650!^ & cff3*(tl_t(i,j,k-1,3,itrc)+ &
1651!^ & tl_t(i,j,k+2,3,itrc)))
1652!^
1653 adfac=w(i,j,k)*ad_fc(i,k)
1654 adfac1=adfac*cff2
1655 adfac2=adfac*cff3
1656 ad_w(i,j,k)=ad_w(i,j,k)+ &
1657 & (cff2*(t(i,j,k ,3,itrc)+ &
1658 & t(i,j,k+1,3,itrc))- &
1659 & cff3*(t(i,j,k-1,3,itrc)+ &
1660 & t(i,j,k+2,3,itrc)))*ad_fc(i,k)
1661 ad_t(i,j,k-1,3,itrc)=ad_t(i,j,k-1,3,itrc)-adfac2
1662 ad_t(i,j,k ,3,itrc)=ad_t(i,j,k ,3,itrc)+adfac1
1663 ad_t(i,j,k+1,3,itrc)=ad_t(i,j,k+1,3,itrc)+adfac1
1664 ad_t(i,j,k+2,3,itrc)=ad_t(i,j,k+2,3,itrc)-adfac2
1665 ad_fc(i,k)=0.0_r8
1666 END DO
1667 END DO
1668 END IF vadv_flux
1669 END DO j_loop1
1670 END DO t_loop2
1671!
1672!-----------------------------------------------------------------------
1673! Time-step adjoint horizontal advection term.
1674!-----------------------------------------------------------------------
1675!
1676 t_loop1 : DO itrc=1,nt(ng)
1677!
1678 k_loop : DO k=1,n(ng)
1679!
1680! Time-step adjoint horizontal advection term.
1681!
1682 hadv_stepping : IF (ad_hadvection(itrc,ng)%MPDATA) THEN
1683 CONTINUE ! not supported
1684 ELSE
1685 DO j=jstr,jend
1686 DO i=istr,iend
1687 cff=dt(ng)*pm(i,j)*pn(i,j)
1688# ifdef DIAGNOSTICS_TS
1689!! DiaTwrk(i,j,k,itrc,iThadv)=-cff3
1690!! DiaTwrk(i,j,k,itrc,iTyadv)=-cff2
1691!! DiaTwrk(i,j,k,itrc,iTxadv)=-cff1
1692# endif
1693!^ tl_t(i,j,k,nnew,itrc)=tl_t(i,j,k,nnew,itrc)-tl_cff3
1694!^
1695 ad_cff3=ad_cff3-ad_t(i,j,k,nnew,itrc)
1696!^ tl_cff3=tl_cff1+tl_cff2
1697!^
1698 ad_cff1=ad_cff1+ad_cff3
1699 ad_cff2=ad_cff2+ad_cff3
1700 ad_cff3=0.0_r8
1701!^ tl_cff2=cff*(tl_FE(i,j+1)-tl_FE(i,j))
1702!^
1703 adfac=cff*ad_cff2
1704 ad_fe(i,j )=ad_fe(i,j )-adfac
1705 ad_fe(i,j+1)=ad_fe(i,j+1)+adfac
1706 ad_cff2=0.0_r8
1707!^ tl_cff1=cff*(tl_FX(i+1,j)-tl_FX(i,j))
1708!^
1709 adfac=cff*ad_cff1
1710 ad_fx(i ,j)=ad_fx(i ,j)-adfac
1711 ad_fx(i+1,j)=ad_fx(i+1,j)+adfac
1712 ad_cff1=0.0_r8
1713 END DO
1714 END DO
1715 END IF hadv_stepping
1716!
1717! Apply adjoint tracers point sources to the horizontal advection
1718! terms, if any.
1719!
1720! Dsrc(is) = 0, flow across grid cell u-face (positive or negative)
1721! Dsrc(is) = 1, flow across grid cell v-face (positive or negative)
1722!
1723 IF (luvsrc(ng)) THEN
1724 DO is=1,nsrc(ng)
1725 isrc=sources(ng)%Isrc(is)
1726 jsrc=sources(ng)%Jsrc(is)
1727 IF (int(sources(ng)%Dsrc(is)).eq.0) THEN
1728 IF ((ad_hadvection(itrc,ng)%MPDATA).or. &
1729 & (ad_hadvection(itrc,ng)%HSIMT)) THEN
1730 lapplysrc=(istrum2.le.isrc).and. &
1731 & (isrc.le.iendp3).and. &
1732 & (jstrvm2.le.jsrc).and. &
1733 & (jsrc.le.jendp2i)
1734 ELSE
1735 lapplysrc=(istr.le.isrc).and. &
1736 & (isrc.le.iend+1).and. &
1737 & (jstr.le.jsrc).and. &
1738 & (jsrc.le.jend)
1739 END IF
1740 IF (lapplysrc) THEN
1741 IF (ltracersrc(itrc,ng)) THEN
1742!^ tl_FX(Isrc,Jsrc)=tl_Huon(Isrc,Jsrc,k)* &
1743!^ & SOURCES(ng)%Tsrc(is,k,itrc)+ &
1744!^ & Huon(Isrc,Jsrc,k)* &
1745!^ & SOURCES(ng)%tl_Tsrc(is,k,itrc)
1746!^
1747 ad_huon(isrc,jsrc,k)=ad_huon(isrc,jsrc,k)+ &
1748 & sources(ng)%Tsrc(is,k,itrc)* &
1749 & ad_fx(isrc,jsrc)
1750 sources(ng)%ad_Tsrc(is,k,itrc)= &
1751 & sources(ng)%ad_Tsrc(is,k,itrc)+ &
1752 & huon(isrc,jsrc,k)* &
1753 & ad_fx(isrc,jsrc)
1754 ad_fx(isrc,jsrc)=0.0_r8
1755# ifdef MASKING
1756 ELSE
1757 IF ((rmask(isrc ,jsrc).eq.0.0_r8).and. &
1758 & (rmask(isrc-1,jsrc).eq.1.0_r8)) THEN
1759!^ tl_FX(Isrc,Jsrc)=tl_Huon(Isrc,Jsrc,k)* &
1760!^ & t(Isrc-1,Jsrc,k,3,itrc)+ &
1761!^ & Huon(Isrc,Jsrc,k)* &
1762!^ & tl_t(Isrc-1,Jsrc,k,3,itrc)
1763!^
1764 ad_t(isrc-1,jsrc,k,3,itrc)= &
1765 & ad_t(isrc-1,jsrc,k,3,itrc)+ &
1766 & huon(isrc,jsrc,k)* &
1767 & ad_fx(isrc,jsrc)
1768 ad_huon(isrc,jsrc,k)=ad_huon(isrc,jsrc,k)+ &
1769 & t(isrc-1,jsrc,k,3,itrc)* &
1770 & ad_fx(isrc,jsrc)
1771 ad_fx(isrc,jsrc)=0.0_r8
1772 ELSE IF ((rmask(isrc ,jsrc).eq.1.0_r8).and. &
1773 & (rmask(isrc-1,jsrc).eq.0.0_r8)) THEN
1774!^ tl_FX(Isrc,Jsrc)=tl_Huon(Isrc,Jsrc,k)* &
1775!^ & t(Isrc ,Jsrc,k,3,itrc)+ &
1776!^ & Huon(Isrc,Jsrc,k)*
1777!^ & tl_t(Isrc ,Jsrc,k,3,itrc)
1778!^
1779 ad_t(isrc ,jsrc,k,3,itrc)= &
1780 & ad_t(isrc ,jsrc,k,3,itrc)+ &
1781 & huon(isrc,jsrc,k)* &
1782 & ad_fx(isrc,jsrc)
1783 ad_huon(isrc,jsrc,k)=ad_huon(isrc,jsrc,k)+ &
1784 & t(isrc ,jsrc,k,3,itrc)* &
1785 & ad_fx(isrc,jsrc)
1786 ad_fx(isrc,jsrc)=0.0_r8
1787 END IF
1788# endif
1789 END IF
1790 END IF
1791 ELSE IF (int(sources(ng)%Dsrc(is)).eq.1) THEN
1792 IF ((ad_hadvection(itrc,ng)%MPDATA).or. &
1793 & (ad_hadvection(itrc,ng)%HSIMT)) THEN
1794 lapplysrc=(istrum2.le.isrc).and. &
1795 & (isrc.le.iendp2i).and. &
1796 & (jstrvm2.le.jsrc).and. &
1797 & (jsrc.le.jendp3)
1798 ELSE
1799 lapplysrc=(istr.le.isrc).and. &
1800 & (isrc.le.iend).and. &
1801 & (jstr.le.jsrc).and. &
1802 & (jsrc.le.jend+1)
1803 END IF
1804 IF (lapplysrc) THEN
1805 IF (ltracersrc(itrc,ng)) THEN
1806!^ tl_FE(Isrc,Jsrc)=tl_Hvom(Isrc,Jsrc,k)* &
1807!^ & SOURCES(ng)%Tsrc(is,k,itrc)
1808!^ & Hvom(Isrc,Jsrc,k)* &
1809!^ & SOURCES(ng)%tl_Tsrc(is,k,itrc)
1810!^
1811 ad_hvom(isrc,jsrc,k)=ad_hvom(isrc,jsrc,k)+ &
1812 & sources(ng)%Tsrc(is,k,itrc)* &
1813 & ad_fe(isrc,jsrc)
1814 sources(ng)%ad_Tsrc(is,k,itrc)= &
1815 & sources(ng)%ad_Tsrc(is,k,itrc)+ &
1816 & hvom(isrc,jsrc,k)* &
1817 & ad_fe(isrc,jsrc)
1818 ad_fe(isrc,jsrc)=0.0_r8
1819# ifdef MASKING
1820 ELSE
1821 IF ((rmask(isrc,jsrc ).eq.0.0_r8).and. &
1822 & (rmask(isrc,jsrc-1).eq.1.0_r8)) THEN
1823!^ tl_FE(Isrc,Jsrc)=tl_Hvom(Isrc,Jsrc,k)* &
1824!^ & t(Isrc,Jsrc-1,k,3,itrc)+ &
1825!^ & Hvom(Isrc,Jsrc,k)* &
1826!^ & tl_t(Isrc,Jsrc-1,k,3,itrc)
1827!^
1828 ad_t(isrc,jsrc-1,k,3,itrc)= &
1829 & ad_t(isrc,jsrc-1,k,3,itrc)+ &
1830 & hvom(isrc,jsrc,k)* &
1831 & ad_fe(isrc,jsrc)
1832 ad_hvom(isrc,jsrc,k)=ad_hvom(isrc,jsrc,k)+ &
1833 & t(isrc,jsrc-1,k,3,itrc)* &
1834 & ad_fe(isrc,jsrc)
1835 ad_fe(isrc,jsrc)=0.0_r8
1836 ELSE IF ((rmask(isrc,jsrc ).eq.1.0_r8).and. &
1837 & (rmask(isrc,jsrc-1).eq.0.0_r8)) THEN
1838!^ tl_FE(Isrc,Jsrc)=tl_Hvom(Isrc,Jsrc,k)* &
1839!^ & t(Isrc,Jsrc ,k,3,itrc)+ &
1840!^ & Hvom(Isrc,Jsrc,k)*
1841!^ & tl_t(Isrc,Jsrc ,k,3,itrc)
1842!^
1843 ad_t(isrc,jsrc ,k,3,itrc)= &
1844 & ad_t(isrc,jsrc ,k,3,itrc)+ &
1845 & hvom(isrc,jsrc,k)* &
1846 & ad_fe(isrc,jsrc)
1847 ad_hvom(isrc,jsrc,k)=ad_hvom(isrc,jsrc,k)+ &
1848 & t(isrc,jsrc ,k,3,itrc)* &
1849 & ad_fe(isrc,jsrc)
1850 ad_fe(isrc,jsrc)=0.0_r8
1851 END IF
1852# endif
1853 END IF
1854 END IF
1855 END IF
1856 END DO
1857 END IF
1858!
1859! Compute adjoint of tracer horizontal advection fluxes.
1860!
1861 hadv_flux : IF (ad_hadvection(itrc,ng)%CENTERED2) THEN
1862!
1863! Second-order, centered differences adjoint horizontal advective fluxes.
1864!
1865 DO j=jstr,jend+1
1866 DO i=istr,iend
1867!^ tl_FE(i,j)=0.5_r8* &
1868!^ & (tl_Hvom(i,j,k)*(t(i,j-1,k,3,itrc)+ &
1869!^ & t(i,j ,k,3,itrc))+ &
1870!^ & Hvom(i,j,k)*(tl_t(i,j-1,k,3,itrc)+ &
1871!^ & tl_t(i,j ,k,3,itrc)))
1872!^
1873 adfac=0.5_r8*ad_fe(i,j)
1874 adfac1=adfac*hvom(i,j,k)
1875 ad_t(i,j-1,k,3,itrc)=ad_t(i,j-1,k,3,itrc)+adfac1
1876 ad_t(i,j ,k,3,itrc)=ad_t(i,j ,k,3,itrc)+adfac1
1877 ad_hvom(i,j,k)=ad_hvom(i,j,k)+ &
1878 & (t(i,j-1,k,3,itrc)+ &
1879 & t(i,j ,k,3,itrc))*adfac
1880 ad_fe(i,j)=0.0_r8
1881 END DO
1882 END DO
1883 DO j=jstr,jend
1884 DO i=istr,iend+1
1885!^ tl_FX(i,j)=0.5_r8* &
1886!^ & (tl_Huon(i,j,k)*(t(i-1,j,k,3,itrc)+ &
1887!^ & t(i ,j,k,3,itrc))+ &
1888!^ & Huon(i,j,k)*(tl_t(i-1,j,k,3,itrc)+ &
1889!^ & tl_t(i ,j,k,3,itrc)))
1890!^
1891 adfac=0.5_r8*ad_fx(i,j)
1892 adfac1=adfac*huon(i,j,k)
1893 ad_t(i-1,j,k,3,itrc)=ad_t(i-1,j,k,3,itrc)+adfac1
1894 ad_t(i ,j,k,3,itrc)=ad_t(i ,j,k,3,itrc)+adfac1
1895 ad_huon(i,j,k)=ad_huon(i,j,k)+ &
1896 & (t(i-1,j,k,3,itrc)+ &
1897 & t(i ,j,k,3,itrc))*adfac
1898 ad_fx(i,j)=0.0_r8
1899 END DO
1900 END DO
1901!
1902 ELSE IF (ad_hadvection(itrc,ng)%MPDATA) THEN
1903!
1904! First-order, upstream differences adjoint horizontal advective
1905! fluxes.
1906!
1907 CONTINUE ! not supported
1908!
1909 ELSE IF (ad_hadvection(itrc,ng)%HSIMT) THEN
1910!
1911! Third High-order Spatial Interpolation at the Middle Temporal level
1912! (HSIMT; Wu and Zhu, 2010) with a Total Variation Diminishing (TVD)
1913! limiter horizontal advection fluxes.
1914!
1915 CONTINUE ! not supported
1916!
1917 ELSE IF ((ad_hadvection(itrc,ng)%AKIMA4).or. &
1918 & (ad_hadvection(itrc,ng)%CENTERED4).or. &
1919 & (ad_hadvection(itrc,ng)%SPLIT_U3).or. &
1920 & (ad_hadvection(itrc,ng)%UPSTREAM3)) THEN
1921!
1922! Fourth-order Akima, fourth-order centered differences, or third-order
1923! upstream-biased horizontal advective fluxes.
1924!
1925 DO j=jstrm1,jendp2
1926 DO i=istr,iend
1927 fe(i,j)=t(i,j ,k,3,itrc)- &
1928 & t(i,j-1,k,3,itrc)
1929# ifdef MASKING
1930 fe(i,j)=fe(i,j)*vmask(i,j)
1931# endif
1932 END DO
1933 END DO
1934 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
1935 IF (domain(ng)%Southern_Edge(tile)) THEN
1936 DO i=istr,iend
1937 fe(i,jstr-1)=fe(i,jstr)
1938 END DO
1939 END IF
1940 END IF
1941 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
1942 IF (domain(ng)%Northern_Edge(tile)) THEN
1943 DO i=istr,iend
1944 fe(i,jend+2)=fe(i,jend+1)
1945 END DO
1946 END IF
1947 END IF
1948!
1949 DO j=jstr-1,jend+1
1950 DO i=istr,iend
1951 IF (ad_hadvection(itrc,ng)%UPSTREAM3) THEN
1952 curv(i,j)=fe(i,j+1)-fe(i,j)
1953 ELSE IF (ad_hadvection(itrc,ng)%AKIMA4) THEN
1954 cff=2.0_r8*fe(i,j+1)*fe(i,j)
1955 IF (cff.gt.eps) THEN
1956 grad(i,j)=cff/(fe(i,j+1)+fe(i,j))
1957 ELSE
1958 grad(i,j)=0.0_r8
1959 END IF
1960 ELSE IF ((ad_hadvection(itrc,ng)%CENTERED4).or. &
1961 & (ad_hadvection(itrc,ng)%SPLIT_U3)) THEN
1962 grad(i,j)=0.5_r8*(fe(i,j+1)+fe(i,j))
1963 END IF
1964 END DO
1965 END DO
1966!
1967 cff1=1.0_r8/6.0_r8
1968 cff2=1.0_r8/3.0_r8
1969 DO j=jstr,jend+1
1970 DO i=istr,iend
1971 IF (ad_hadvection(itrc,ng)%UPSTREAM3) THEN
1972!^ tl_FE(i,j)=0.5_r8* &
1973!^ & (tl_Hvom(i,j,k)* &
1974!^ & (t(i,j-1,k,3,itrc)+ &
1975!^ & t(i,j ,k,3,itrc))+ &
1976!^ & Hvom(i,j,k)* &
1977!^ & (tl_t(i,j-1,k,3,itrc)+ &
1978!^ & tl_t(i,j ,k,3,itrc)))- &
1979!^ & cff1* &
1980!^ & (tl_curv(i,j-1)*MAX(Hvom(i,j,k),0.0_r8)+ &
1981!^ & curv(i,j-1)* &
1982!^ & (0.5_r8+SIGN(0.5_r8, Hvom(i,j,k)))* &
1983!^ & tl_Hvom(i,j,k)+ &
1984!^ & tl_curv(i,j )*MIN(Hvom(i,j,k),0.0_r8)+ &
1985!^ & curv(i,j )* &
1986!^ & (0.5_r8+SIGN(0.5_r8,-Hvom(i,j,k)))* &
1987!^ & tl_Hvom(i,j,k))
1988!^
1989 adfac=0.5_r8*ad_fe(i,j)
1990 adfac1=adfac*hvom(i,j,k)
1991 adfac2=cff1*ad_fe(i,j)
1992 ad_hvom(i,j,k)=ad_hvom(i,j,k)+ &
1993 & (t(i,j-1,k,3,itrc)+ &
1994 & t(i,j ,k,3,itrc))*adfac- &
1995 & (curv(i,j-1)* &
1996 & (0.5_r8+sign(0.5_r8, hvom(i,j,k)))+ &
1997 & curv(i,j )* &
1998 & (0.5_r8+sign(0.5_r8,-hvom(i,j,k))))* &
1999 & adfac2
2000 ad_t(i,j-1,k,3,itrc)=ad_t(i,j-1,k,3,itrc)+adfac1
2001 ad_t(i,j ,k,3,itrc)=ad_t(i,j ,k,3,itrc)+adfac1
2002 ad_curv(i,j-1)=ad_curv(i,j-1)- &
2003 & max(hvom(i,j,k),0.0_r8)*adfac2
2004 ad_curv(i,j )=ad_curv(i,j )- &
2005 & min(hvom(i,j,k),0.0_r8)*adfac2
2006 ad_fe(i,j)=0.0_r8
2007 ELSE IF ((ad_hadvection(itrc,ng)%AKIMA4).or. &
2008 & (ad_hadvection(itrc,ng)%CENTERED4).or. &
2009 & (ad_hadvection(itrc,ng)%SPLIT_U3)) THEN
2010!^ tl_FE(i,j)=0.5_r8* &
2011!^ & (tl_Hvom(i,j,k)* &
2012!^ & (t(i,j-1,k,3,itrc)+ &
2013!^ & t(i,j ,k,3,itrc)- &
2014!^ & cff2*(grad(i,j )- &
2015!^ & grad(i,j-1)))+ &
2016!^ & Hvom(i,j,k)* &
2017!^ & (tl_t(i,j-1,k,3,itrc)+ &
2018!^ & tl_t(i,j ,k,3,itrc)- &
2019!^ & cff2*(tl_grad(i,j )- &
2020!^ & tl_grad(i,j-1))))
2021!^
2022 adfac=0.5_r8*ad_fe(i,j)
2023 adfac1=adfac*hvom(i,j,k)
2024 adfac2=adfac1*cff2
2025 ad_hvom(i,j,k)=ad_hvom(i,j,k)+ &
2026 & adfac*(t(i,j-1,k,3,itrc)+ &
2027 & t(i,j ,k,3,itrc)- &
2028 & cff2*(grad(i,j )- &
2029 & grad(i,j-1)))
2030 ad_t(i,j-1,k,3,itrc)=ad_t(i,j-1,k,3,itrc)+adfac1
2031 ad_t(i,j ,k,3,itrc)=ad_t(i,j ,k,3,itrc)+adfac1
2032 ad_grad(i,j-1)=ad_grad(i,j-1)+adfac2
2033 ad_grad(i,j )=ad_grad(i,j )-adfac2
2034 ad_fe(i,j)=0.0_r8
2035 END IF
2036 END DO
2037 END DO
2038!
2039 DO j=jstr-1,jend+1
2040 DO i=istr,iend
2041 IF (ad_hadvection(itrc,ng)%UPSTREAM3) THEN
2042!^ tl_curv(i,j)=tl_FE(i,j+1)-tl_FE(i,j)
2043!^
2044 ad_fe(i,j )=ad_fe(i,j )-ad_curv(i,j)
2045 ad_fe(i,j+1)=ad_fe(i,j+1)+ad_curv(i,j)
2046 ad_curv(i,j)=0.0_r8
2047 ELSE IF (ad_hadvection(itrc,ng)%AKIMA4) THEN
2048 cff=2.0_r8*fe(i,j+1)*fe(i,j)
2049 IF (cff.gt.eps) THEN
2050!^ tl_grad(i,j)=((FE(i,j+1)+FE(i,j))*tl_cff- &
2051!^ & cff*(tl_FE(i,j+1)+tl_FE(i,j)))/ &
2052!^ & ((FE(i,j+1)+FE(i,j))* &
2053!^ & (FE(i,j+1)+FE(i,j)))
2054!^
2055 adfac=ad_grad(i,j)/ &
2056 & ((fe(i,j+1)+fe(i,j))*(fe(i,j+1)+fe(i,j)))
2057 adfac1=adfac*cff
2058 ad_fe(i,j )=ad_fe(i,j)-adfac1
2059 ad_fe(i,j+1)=ad_fe(i,j+1)-adfac1
2060 ad_cff=ad_cff+(fe(i,j+1)+fe(i,j))*adfac
2061 ad_grad(i,j)=0.0_r8
2062 ELSE
2063!^ tl_grad(i,j)=0.0_r8
2064!^
2065 ad_grad(i,j)=0.0_r8
2066 END IF
2067!^ tl_cff=2.0_r8*(tl_FE(i,j+1)*FE(i,j)+ &
2068!^ & FE(i,j+1)*tl_FE(i,j))
2069!^
2070 adfac=2.0_r8*ad_cff
2071 ad_fe(i,j )=ad_fe(i,j )+fe(i,j+1)*adfac
2072 ad_fe(i,j+1)=ad_fe(i,j+1)+fe(i,j )*adfac
2073 ad_cff=0.0_r8
2074 ELSE IF ((ad_hadvection(itrc,ng)%CENTERED4).or. &
2075 & (ad_hadvection(itrc,ng)%SPLIT_U3)) THEN
2076!^ tl_grad(i,j)=0.5_r8*(tl_FE(i,j+1)+tl_FE(i,j))
2077!^
2078 adfac=0.5_r8*ad_grad(i,j)
2079 ad_fe(i,j )=ad_fe(i,j )+adfac
2080 ad_fe(i,j+1)=ad_fe(i,j+1)+adfac
2081 ad_grad(i,j)=0.0_r8
2082 END IF
2083 END DO
2084 END DO
2085 IF (.not.(compositegrid(inorth,ng).or.nsperiodic(ng))) THEN
2086 IF (domain(ng)%Northern_Edge(tile)) THEN
2087 DO i=istr,iend
2088!^ tl_FE(i,Jend+2)=tl_FE(i,Jend+1)
2089!^
2090 ad_fe(i,jend+1)=ad_fe(i,jend+1)+ad_fe(i,jend+2)
2091 ad_fe(i,jend+2)=0.0_r8
2092 END DO
2093 END IF
2094 END IF
2095 IF (.not.(compositegrid(isouth,ng).or.nsperiodic(ng))) THEN
2096 IF (domain(ng)%Southern_Edge(tile)) THEN
2097 DO i=istr,iend
2098!^ tl_FE(i,Jstr-1)=tl_FE(i,Jstr)
2099!^
2100 ad_fe(i,jstr)=ad_fe(i,jstr)+ad_fe(i,jstr-1)
2101 ad_fe(i,jstr-1)=0.0_r8
2102 END DO
2103 END IF
2104 END IF
2105!
2106 DO j=jstrm1,jendp2
2107 DO i=istr,iend
2108# ifdef MASKING
2109!^ tl_FE(i,j)=tl_FE(i,j)*vmask(i,j)
2110!^
2111 ad_fe(i,j)=ad_fe(i,j)*vmask(i,j)
2112# endif
2113!^ tl_FE(i,j)=tl_t(i,j ,k,3,itrc)- &
2114!^ & tl_t(i,j-1,k,3,itrc)
2115!^
2116 ad_t(i,j-1,k,3,itrc)=ad_t(i,j-1,k,3,itrc)-ad_fe(i,j)
2117 ad_t(i,j ,k,3,itrc)=ad_t(i,j ,k,3,itrc)+ad_fe(i,j)
2118 ad_fe(i,j)=0.0_r8
2119 END DO
2120 END DO
2121!
2122 DO j=jstr,jend
2123 DO i=istrm1,iendp2
2124 fx(i,j)=t(i ,j,k,3,itrc)- &
2125 & t(i-1,j,k,3,itrc)
2126# ifdef MASKING
2127 fx(i,j)=fx(i,j)*umask(i,j)
2128# endif
2129 END DO
2130 END DO
2131 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
2132 IF (domain(ng)%Western_Edge(tile)) THEN
2133 DO j=jstr,jend
2134 fx(istr-1,j)=fx(istr,j)
2135 END DO
2136 END IF
2137 END IF
2138 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
2139 IF (domain(ng)%Eastern_Edge(tile)) THEN
2140 DO j=jstr,jend
2141 fx(iend+2,j)=fx(iend+1,j)
2142 END DO
2143 END IF
2144 END IF
2145!
2146 DO j=jstr,jend
2147 DO i=istr-1,iend+1
2148 IF (ad_hadvection(itrc,ng)%UPSTREAM3) THEN
2149 curv(i,j)=fx(i+1,j)-fx(i,j)
2150 ELSE IF (ad_hadvection(itrc,ng)%AKIMA4) THEN
2151 cff=2.0_r8*fx(i+1,j)*fx(i,j)
2152 IF (cff.gt.eps) THEN
2153 grad(i,j)=cff/(fx(i+1,j)+fx(i,j))
2154 ELSE
2155 grad(i,j)=0.0_r8
2156 END IF
2157 ELSE IF ((ad_hadvection(itrc,ng)%CENTERED4).or. &
2158 & (ad_hadvection(itrc,ng)%SPLIT_U3)) THEN
2159 grad(i,j)=0.5_r8*(fx(i+1,j)+fx(i,j))
2160 END IF
2161 END DO
2162 END DO
2163!
2164 cff1=1.0_r8/6.0_r8
2165 cff2=1.0_r8/3.0_r8
2166 DO j=jstr,jend
2167 DO i=istr,iend+1
2168 IF (ad_hadvection(itrc,ng)%UPSTREAM3) THEN
2169!^ tl_FX(i,j)=0.5_r8* &
2170!^ & (tl_Huon(i,j,k)* &
2171!^ & (t(i-1,j,k,3,itrc)+ &
2172!^ & t(i ,j,k,3,itrc))+ &
2173!^ & Huon(i,j,k)* &
2174!^ & (tl_t(i-1,j,k,3,itrc)+ &
2175!^ & tl_t(i ,j,k,3,itrc)))- &
2176!^ & cff1* &
2177!^ & (tl_curv(i-1,j)*MAX(Huon(i,j,k),0.0_r8)+ &
2178!^ & curv(i-1,j)* &
2179!^ & (0.5_r8+SIGN(0.5_r8, Huon(i,j,k)))* &
2180!^ & tl_Huon(i,j,k)+ &
2181!^ & tl_curv(i ,j)*MIN(Huon(i,j,k),0.0_r8)+ &
2182!^ & curv(i ,j)* &
2183!^ & (0.5_r8+SIGN(0.5_r8,-Huon(i,j,k)))* &
2184!^ & tl_Huon(i,j,k))
2185!^
2186 adfac=0.5_r8*ad_fx(i,j)
2187 adfac1=adfac*huon(i,j,k)
2188 adfac2=cff1*ad_fx(i,j)
2189 ad_huon(i,j,k)=ad_huon(i,j,k)+ &
2190 & (t(i-1,j,k,3,itrc)+ &
2191 & t(i ,j,k,3,itrc))*adfac- &
2192 & (curv(i-1,j)* &
2193 & (0.5_r8+sign(0.5_r8, huon(i,j,k)))+ &
2194 & curv(i ,j)* &
2195 & (0.5_r8+sign(0.5_r8,-huon(i,j,k))))* &
2196 & adfac2
2197 ad_t(i-1,j,k,3,itrc)=ad_t(i-1,j,k,3,itrc)+adfac1
2198 ad_t(i ,j,k,3,itrc)=ad_t(i ,j,k,3,itrc)+adfac1
2199 ad_curv(i-1,j)=ad_curv(i-1,j)- &
2200 & max(huon(i,j,k),0.0_r8)*adfac2
2201 ad_curv(i ,j)=ad_curv(i ,j)- &
2202 & min(huon(i,j,k),0.0_r8)*adfac2
2203 ad_fx(i,j)=0.0_r8
2204 ELSE IF ((ad_hadvection(itrc,ng)%AKIMA4).or. &
2205 & (ad_hadvection(itrc,ng)%CENTERED4).or. &
2206 & (ad_hadvection(itrc,ng)%SPLIT_U3)) THEN
2207!^ tl_FX(i,j)=0.5_r8* &
2208!^ & (tl_Huon(i,j,k)* &
2209!^ & (t(i-1,j,k,3,itrc)+ &
2210!^ & t(i ,j,k,3,itrc)- &
2211!^ & cff2*(grad(i ,j)- &
2212!^ & grad(i-1,j)))+ &
2213!^ & Huon(i,j,k)* &
2214!^ & (tl_t(i-1,j,k,3,itrc)+ &
2215!^ & tl_t(i ,j,k,3,itrc)- &
2216!^ & cff2*(tl_grad(i ,j)- &
2217!^ & tl_grad(i-1,j))))
2218!^
2219 adfac=0.5_r8*ad_fx(i,j)
2220 adfac1=adfac*huon(i,j,k)
2221 adfac2=adfac1*cff2
2222 ad_huon(i,j,k)=ad_huon(i,j,k)+ &
2223 & adfac*(t(i-1,j,k,3,itrc)+ &
2224 & t(i ,j,k,3,itrc)- &
2225 & cff2*(grad(i ,j)- &
2226 & grad(i-1,j)))
2227 ad_t(i-1,j,k,3,itrc)=ad_t(i-1,j,k,3,itrc)+adfac1
2228 ad_t(i ,j,k,3,itrc)=ad_t(i ,j,k,3,itrc)+adfac1
2229 ad_grad(i-1,j)=ad_grad(i-1,j)+adfac2
2230 ad_grad(i ,j)=ad_grad(i ,j)-adfac2
2231 ad_fx(i,j)=0.0_r8
2232 END IF
2233 END DO
2234 END DO
2235!
2236 DO j=jstr,jend
2237 DO i=istr-1,iend+1
2238 IF (ad_hadvection(itrc,ng)%UPSTREAM3) THEN
2239!^ tl_curv(i,j)=tl_FX(i+1,j)-tl_FX(i,j)
2240!^
2241 ad_fx(i ,j)=ad_fx(i ,j)-ad_curv(i,j)
2242 ad_fx(i+1,j)=ad_fx(i+1,j)+ad_curv(i,j)
2243 ad_curv(i,j)=0.0_r8
2244 ELSE IF (ad_hadvection(itrc,ng)%AKIMA4) THEN
2245 cff=2.0_r8*fx(i+1,j)*fx(i,j)
2246 IF (cff.gt.eps) THEN
2247!^ tl_grad(i,j)=((FX(i+1,j)+FX(i,j))*tl_cff- &
2248!^ & cff*(tl_FX(i+1,j)+tl_FX(i,j)))/ &
2249!^ & ((FX(i+1,j)+FX(i,j))* &
2250!^ & (FX(i+1,j)+FX(i,j)))
2251!^
2252 adfac=ad_grad(i,j)/ &
2253 & ((fx(i+1,j)+fx(i,j))*(fx(i+1,j)+fx(i,j)))
2254 adfac1=adfac*cff
2255 ad_fx(i ,j)=ad_fx(i ,j)-adfac1
2256 ad_fx(i+1,j)=ad_fx(i+1,j)-adfac1
2257 ad_cff=ad_cff+(fx(i+1,j)+fx(i,j))*adfac
2258 ad_grad(i,j)=0.0_r8
2259 ELSE
2260!^ tl_grad(i,j)=0.0_r8
2261!^
2262 ad_grad(i,j)=0.0_r8
2263 END IF
2264 ELSE IF ((ad_hadvection(itrc,ng)%CENTERED4).or. &
2265 & (ad_hadvection(itrc,ng)%SPLIT_U3)) THEN
2266!^ tl_grad(i,j)=0.5_r8*(tl_FX(i+1,j)+tl_FX(i,j))
2267!^
2268 adfac=0.5_r8*ad_grad(i,j)
2269 ad_fx(i ,j)=ad_fx(i ,j)+adfac
2270 ad_fx(i+1,j)=ad_fx(i+1,j)+adfac
2271 ad_grad(i,j)=0.0_r8
2272 END IF
2273!^ tl_cff=2.0_r8*(tl_FX(i+1,j)*FX(i,j)+ &
2274!^ & FX(i+1,j)*tl_FX(i,j))
2275!^
2276 adfac=2.0_r8*ad_cff
2277 ad_fx(i ,j)=ad_fx(i ,j)+fx(i+1,j)*adfac
2278 ad_fx(i+1,j)=ad_fx(i+1,j)+fx(i ,j)*adfac
2279 ad_cff=0.0_r8
2280 END DO
2281 END DO
2282 IF (.not.(compositegrid(ieast,ng).or.ewperiodic(ng))) THEN
2283 IF (domain(ng)%Eastern_Edge(tile)) THEN
2284 DO j=jstr,jend
2285!^ tl_FX(Iend+2,j)=tl_FX(Iend+1,j)
2286!^
2287 ad_fx(iend+1,j)=ad_fx(iend+1,j)+ad_fx(iend+2,j)
2288 ad_fx(iend+2,j)=0.0_r8
2289 END DO
2290 END IF
2291 END IF
2292 IF (.not.(compositegrid(iwest,ng).or.ewperiodic(ng))) THEN
2293 IF (domain(ng)%Western_Edge(tile)) THEN
2294 DO j=jstr,jend
2295!^ tl_FX(Istr-1,j)=tl_FX(Istr,j)
2296!^
2297 ad_fx(istr,j)=ad_fx(istr,j)+ad_fx(istr-1,j)
2298 ad_fx(istr-1,j)=0.0_r8
2299 END DO
2300 END IF
2301 END IF
2302!
2303 DO j=jstr,jend
2304 DO i=istrm1,iendp2
2305# ifdef MASKING
2306!^ tl_FX(i,j)=tl_FX(i,j)*umask(i,j)
2307!^
2308 ad_fx(i,j)=ad_fx(i,j)*umask(i,j)
2309# endif
2310!^ tl_FX(i,j)=tl_t(i ,j,k,nstp,itrc)- &
2311!^ & tl_t(i-1,j,k,nstp,itrc)
2312!^
2313 ad_t(i-1,j,k,3,itrc)=ad_t(i-1,j,k,3,itrc)-ad_fx(i,j)
2314 ad_t(i ,j,k,3,itrc)=ad_t(i ,j,k,3,itrc)+ad_fx(i,j)
2315 ad_fx(i,j)=0.0_r8
2316 END DO
2317 END DO
2318 END IF hadv_flux
2319 END DO k_loop
2320
2321# ifdef AD_SUPPORTED
2322!
2323! The MPDATA algorithm requires a three-point footprint, so exchange
2324! boundary data on t(:,:,:,nnew,:) so other processes computed earlier
2325! (horizontal diffusion, biology, or sediment) are accounted.
2326!
2327 IF ((ad_hadvection(itrc,ng)%MPDATA).or. &
2328 & (ad_hadvection(itrc,ng)%HSIMT)) THEN
2329# ifdef DISTRIBUTE
2330! > CALL mp_exchange3d (ng, tile, iNLM, 1, &
2331!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
2332!^ & NghostPoints, &
2333!^ & EWperiodic(ng), NSperiodic(ng), &
2334!^ & tl_t(:,:,:,nnew,itrc))
2335!^
2336 CALL ad_mp_exchange3d (ng, tile, iadm, 1, &
2337 & lbi, ubi, lbj, ubj, 1, n(ng), &
2338 & nghostpoints, &
2339 & ewperiodic(ng), nsperiodic(ng), &
2340 & ad_t(:,:,:,nnew,itrc))
2341# endif
2342 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
2343!^ CALL exchange_r3d_tile (ng, tile, &
2344!^ & LBi, UBi, LBj, UBj, 1, N(ng), &
2345!^ & tl_t(:,:,:,nnew,itrc))
2346!^
2347 CALL ad_exchange_r3d_tile (ng, tile, &
2348 & lbi, ubi, lbj, ubj, 1, n(ng), &
2349 & ad_t(:,:,:,nnew,itrc))
2350 END IF
2351 END IF
2352# endif
2353
2354 END DO t_loop1
2355!
2356! Compute adjoint inverse thickness.
2357!
2358 IF (lhsimt) THEN
2359 DO k=1,n(ng)
2360 DO j=jstrm2,jendp2
2361 DO i=istrm2,iendp2
2362 ohz(i,j,k)=1.0_r8/hz(i,j,k)
2363!^ tl_oHz(i,j,k)=-oHz(i,j,k)*oHz(i,j,k)*tl_Hz(i,j,k)
2364!^
2365 ad_hz(i,j,k)=ad_hz(i,j,k)- &
2366 & ohz(i,j,k)*ohz(i,j,k)*ad_ohz(i,j,k)
2367 ad_ohz(i,j,k)=0.0_r8
2368 END DO
2369 END DO
2370 END DO
2371 ELSE
2372 DO k=1,n(ng)
2373 DO j=jstr,jend
2374 DO i=istr,iend
2375 ohz(i,j,k)=1.0_r8/hz(i,j,k)
2376!^ tl_oHz(i,j,k)=-oHz(i,j,k)*oHz(i,j,k)*tl_Hz(i,j,k)
2377!^
2378 ad_hz(i,j,k)=ad_hz(i,j,k)- &
2379 & ohz(i,j,k)*ohz(i,j,k)*ad_ohz(i,j,k)
2380 ad_ohz(i,j,k)=0.0_r8
2381 END DO
2382 END DO
2383 END DO
2384 END IF
2385!
2386 RETURN
subroutine ad_exchange_r3d_tile(ng, tile, lbi, ubi, lbj, ubj, lbk, ubk, ad_a)
subroutine, public ad_t3dbc_tile(ng, tile, itrc, ic, lbi, ubi, lbj, ubj, ubk, ubt, imins, imaxs, jmins, jmaxs, nstp, nout, ad_t)
Definition ad_t3dbc_im.F:57
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
integer nat
Definition mod_param.F:499
type(t_adv), dimension(:,:), allocatable ad_hadvection
Definition mod_param.F:407
integer, parameter inlm
Definition mod_param.F:662
integer, dimension(:), allocatable n
Definition mod_param.F:479
integer nghostpoints
Definition mod_param.F:710
type(t_domain), dimension(:), allocatable domain
Definition mod_param.F:329
integer, dimension(:), allocatable nt
Definition mod_param.F:489
type(t_adv), dimension(:,:), allocatable ad_vadvection
Definition mod_param.F:408
integer npt
Definition mod_param.F:505
logical, dimension(:), allocatable luvsrc
logical, dimension(:,:), allocatable ltracersrc
real(dp), dimension(:), allocatable dt
real(r8) lambda
logical, dimension(:), allocatable ewperiodic
integer, parameter iwest
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable lwsrc
logical, dimension(:,:), allocatable compositegrid
integer, parameter isouth
integer, dimension(:), pointer inert
integer, parameter ieast
logical, dimension(:,:), allocatable ltracerclm
integer, parameter inorth
logical, dimension(:,:), allocatable lnudgetclm
type(t_sources), dimension(:), allocatable sources
Definition mod_sources.F:90
integer, dimension(:), allocatable nsrc
Definition mod_sources.F:97
subroutine ad_mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, ad_a, ad_b, ad_c, ad_d)
subroutine ad_mp_exchange4d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, lbt, ubt, nghost, ew_periodic, ns_periodic, ad_a, ad_b, ad_c)
subroutine mp_exchange3d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, lbk, ubk, nghost, ew_periodic, ns_periodic, a, b, c, d)

References ad_exchange_3d_mod::ad_exchange_r3d_tile(), mod_param::ad_hadvection, mp_exchange_mod::ad_mp_exchange3d(), mp_exchange_mod::ad_mp_exchange4d(), ad_t3dbc_mod::ad_t3dbc_tile(), mod_param::ad_vadvection, mod_clima::clima, mod_scalars::compositegrid, mod_param::domain, mod_scalars::dt, mod_scalars::ewperiodic, mod_param::iadm, mod_scalars::ieast, mod_scalars::inert, mod_param::inlm, mod_scalars::inorth, mod_scalars::isouth, mod_scalars::iwest, mod_scalars::lambda, mod_scalars::lnudgetclm, mod_scalars::ltracerclm, mod_scalars::ltracersrc, mod_scalars::luvsrc, mod_scalars::lwsrc, mp_exchange_mod::mp_exchange3d(), mod_param::nghostpoints, mod_param::npt, mod_scalars::nsperiodic, mod_sources::nsrc, and mod_sources::sources.

Referenced by ad_step3d_t().

Here is the call graph for this function:
Here is the caller graph for this function: