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

Functions/Subroutines

subroutine, public ad_uv3dmix2 (ng, tile)
 
subroutine ad_uv3dmix2_geo_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nnew, pmask, rmask, umask, vmask, om_p, om_r, om_u, om_v, on_p, on_r, on_u, on_v, pm, pn, hz, ad_hz, z_r, ad_z_r, visc3d_r, ad_visc3d_r, visc2_p, visc2_r, u, v, ad_u, ad_v, ad_rufrc, ad_rvfrc)
 
subroutine ad_uv3dmix2_s_tile (ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, nnew, pmask, hz, ad_hz, om_p, om_r, on_p, on_r, pm, pmon_p, pmon_r, pn, pnom_p, pnom_r, visc2_p, visc2_r, u, v, ad_rufrc, ad_rvfrc, ad_u, ad_v)
 

Function/Subroutine Documentation

◆ ad_uv3dmix2()

subroutine public ad_uv3dmix2_mod::ad_uv3dmix2 ( integer, intent(in) ng,
integer, intent(in) tile )

Definition at line 42 of file ad_uv3dmix2_geo.h.

43!***********************************************************************
44!
45 USE mod_param
46 USE mod_coupling
47#ifdef DIAGNOSTICS_UV
48!! USE mod_diags
49#endif
50 USE mod_grid
51 USE mod_mixing
52 USE mod_ocean
53 USE mod_stepping
54!
55! Imported variable declarations.
56!
57 integer, intent(in) :: ng, tile
58!
59! Local variable declarations.
60!
61 character (len=*), parameter :: MyFile = &
62 & __FILE__
63!
64#include "tile.h"
65!
66#ifdef PROFILE
67 CALL wclock_on (ng, iadm, 31, __line__, myfile)
68#endif
69 CALL ad_uv3dmix2_geo_tile (ng, tile, &
70 & lbi, ubi, lbj, ubj, &
71 & imins, imaxs, jmins, jmaxs, &
72 & nrhs(ng), nnew(ng), &
73#ifdef MASKING
74 & grid(ng) % pmask, &
75 & grid(ng) % rmask, &
76 & grid(ng) % umask, &
77 & grid(ng) % vmask, &
78#endif
79 & grid(ng) % om_p, &
80 & grid(ng) % om_r, &
81 & grid(ng) % om_u, &
82 & grid(ng) % om_v, &
83 & grid(ng) % on_p, &
84 & grid(ng) % on_r, &
85 & grid(ng) % on_u, &
86 & grid(ng) % on_v, &
87 & grid(ng) % pm, &
88 & grid(ng) % pn, &
89 & grid(ng) % Hz, &
90 & grid(ng) % ad_Hz, &
91 & grid(ng) % z_r, &
92 & grid(ng) % ad_z_r, &
93#ifdef VISC_3DCOEF
94 & mixing(ng) % visc3d_r, &
95 & mixing(ng) % ad_visc3d_r, &
96#else
97 & mixing(ng) % visc2_p, &
98 & mixing(ng) % visc2_r, &
99#endif
100#ifdef DIAGNOSTICS_UV
101!! & DIAGS(ng) % DiaRUfrc, &
102!! & DIAGS(ng) % DiaRVfrc, &
103!! & DIAGS(ng) % DiaU3wrk, &
104!! & DIAGS(ng) % DiaV3wrk, &
105#endif
106 & ocean(ng) % u, &
107 & ocean(ng) % v, &
108 & ocean(ng) % ad_u, &
109 & ocean(ng) % ad_v, &
110 & coupling(ng) % ad_rufrc, &
111 & coupling(ng) % ad_rvfrc)
112#ifdef PROFILE
113 CALL wclock_off (ng, iadm, 31, __line__, myfile)
114#endif
115!
116 RETURN
type(t_coupling), dimension(:), allocatable coupling
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
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_uv3dmix2_geo_tile(), mod_coupling::coupling, mod_grid::grid, mod_param::iadm, mod_mixing::mixing, mod_stepping::nnew, mod_stepping::nrhs, mod_ocean::ocean, wclock_off(), and wclock_on().

Referenced by ad_rhs3d_mod::ad_rhs3d().

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

◆ ad_uv3dmix2_geo_tile()

subroutine ad_uv3dmix2_mod::ad_uv3dmix2_geo_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) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) rmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) umask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) vmask,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_u,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) hz,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(inout) ad_hz,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) z_r,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(inout) ad_z_r,
real(r8), dimension(lbi:,lbj:,:), intent(in) visc3d_r,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(inout) ad_visc3d_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_r,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(in) u,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(in) v,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(inout) ad_u,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(inout) ad_v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) ad_rufrc,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) ad_rvfrc )
private

Definition at line 121 of file ad_uv3dmix2_geo.h.

144!***********************************************************************
145!
146 USE mod_param
147 USE mod_scalars
148!
149! Imported variable declarations.
150!
151 integer, intent(in) :: ng, tile
152 integer, intent(in) :: LBi, UBi, LBj, UBj
153 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
154 integer, intent(in) :: nrhs, nnew
155
156#ifdef ASSUMED_SHAPE
157# ifdef MASKING
158 real(r8), intent(in) :: pmask(LBi:,LBj:)
159 real(r8), intent(in) :: rmask(LBi:,LBj:)
160 real(r8), intent(in) :: umask(LBi:,LBj:)
161 real(r8), intent(in) :: vmask(LBi:,LBj:)
162# endif
163 real(r8), intent(in) :: om_p(LBi:,LBj:)
164 real(r8), intent(in) :: om_r(LBi:,LBj:)
165 real(r8), intent(in) :: om_u(LBi:,LBj:)
166 real(r8), intent(in) :: om_v(LBi:,LBj:)
167 real(r8), intent(in) :: on_p(LBi:,LBj:)
168 real(r8), intent(in) :: on_r(LBi:,LBj:)
169 real(r8), intent(in) :: on_u(LBi:,LBj:)
170 real(r8), intent(in) :: on_v(LBi:,LBj:)
171 real(r8), intent(in) :: pm(LBi:,LBj:)
172 real(r8), intent(in) :: pn(LBi:,LBj:)
173 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
174 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
175# ifdef VISC_3DCOEF
176 real(r8), intent(in) :: visc3d_r(LBi:,LBj:,:)
177# else
178 real(r8), intent(in) :: visc2_p(LBi:,LBj:)
179 real(r8), intent(in) :: visc2_r(LBi:,LBj:)
180# endif
181 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
182 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
183
184# ifdef DIAGNOSTICS_UV
185!! real(r8), intent(inout) :: DiaRUfrc(LBi:,LBj:,:,:)
186!! real(r8), intent(inout) :: DiaRVfrc(LBi:,LBj:,:,:)
187!! real(r8), intent(inout) :: DiaU3wrk(LBi:,LBj:,:,:)
188!! real(r8), intent(inout) :: DiaV3wrk(LBi:,LBj:,:,:)
189# endif
190 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
191 real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:)
192 real(r8), intent(inout) :: ad_rufrc(LBi:,LBj:)
193 real(r8), intent(inout) :: ad_rvfrc(LBi:,LBj:)
194 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
195 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
196# ifdef VISC_3DCOEF
197 real(r8), intent(inout) :: ad_visc3d_r(LBi:,LBj:,:)
198# endif
199#else
200# ifdef MASKING
201 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
202 real(r8), intent(in) :: rmask(LBi:UBi,LBj:UBj)
203 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
204 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
205# endif
206 real(r8), intent(in) :: om_p(LBi:UBi,LBj:UBj)
207 real(r8), intent(in) :: om_r(LBi:UBi,LBj:UBj)
208 real(r8), intent(in) :: om_u(LBi:UBi,LBj:UBj)
209 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
210 real(r8), intent(in) :: on_p(LBi:UBi,LBj:UBj)
211 real(r8), intent(in) :: on_r(LBi:UBi,LBj:UBj)
212 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
213 real(r8), intent(in) :: on_v(LBi:UBi,LBj:UBj)
214 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
215 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
216 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
217 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
218 real(r8), intent(in) :: visc2_p(LBi:UBi,LBj:UBj)
219 real(r8), intent(in) :: visc2_r(LBi:UBi,LBj:UBj)
220
221 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
222 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
223
224# ifdef DIAGNOSTICS_UV
225!! real(r8), intent(inout) :: DiaRUfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
226!! real(r8), intent(inout) :: DiaRVfrc(LBi:UBi,LBj:UBj,3,NDM2d-1)
227!! real(r8), intent(inout) :: DiaU3wrk(LBi:UBi,LBj:UBj,N(ng),NDM3d)
228!! real(r8), intent(inout) :: DiaV3wrk(LBi:UBi,LBj:UBj,N(ng),NDM3d)
229# endif
230
231 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
232 real(r8), intent(inout) :: ad_z_r(LBi:UBi,LBj:UBj,N(ng))
233 real(r8), intent(inout) :: ad_rufrc(LBi:UBi,LBj:UBj)
234 real(r8), intent(inout) :: ad_rvfrc(LBi:UBi,LBj:UBj)
235 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
236 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
237# ifdef VISC_3DCOEF
238 real(r8), intent(inout) :: ad_visc3d_r(LBi:UBi,LBj:UBj,N(ng))
239# endif
240
241#endif
242!
243! Local variable declarations.
244!
245 integer :: i, j, k, kk, kt, k1, k1b, k2, k2b
246
247 real(r8) :: cff, fac1, fac2, pm_p, pn_p
248 real(r8) :: cff1, cff2, cff3, cff4
249 real(r8) :: cff5, cff6, cff7, cff8
250 real(r8) :: dmUdz, dnUdz, dmVdz, dnVdz
251#ifdef VISC_3DCOEF
252 real(r8) :: visc_p
253 real(r8) :: ad_fac1, ad_fac2, ad_visc_p
254#endif
255 real(r8) :: adfac, ad_cff
256 real(r8) :: adfac1, adfac2, adfac3, adfac4, adfac5, adfac6
257 real(r8) :: ad_cff1, ad_cff2, ad_cff3, ad_cff4
258 real(r8) :: ad_cff5, ad_cff6, ad_cff7, ad_cff8
259 real(r8) :: ad_dmUdz, ad_dnUdz, ad_dmVdz, ad_dnVdz
260
261 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: UFx
262 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: VFe
263
264 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFe
265 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFx
266 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFe
267 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFx
268
269 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dmUde
270 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dmVde
271 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dnUdx
272 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dnVdx
273 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dUdz
274 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dVdz
275 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZde_p
276 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZde_r
277 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZdx_p
278 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: dZdx_r
279
280 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_UFse
281 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_UFsx
282 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_VFse
283 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_VFsx
284 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dmUde
285 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dmVde
286 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dnUdx
287 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dnVdx
288 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dUdz
289 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dVdz
290 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dZde_p
291 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dZde_r
292 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dZdx_p
293 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,2) :: ad_dZdx_r
294
295#include "set_bounds.h"
296!
297!-----------------------------------------------------------------------
298! Initialize private adjoint variables and arrays.
299!-----------------------------------------------------------------------
300!
301 ad_cff=0.0_r8
302 ad_cff1=0.0_r8
303 ad_cff2=0.0_r8
304 ad_cff3=0.0_r8
305 ad_cff4=0.0_r8
306 ad_cff5=0.0_r8
307 ad_cff6=0.0_r8
308 ad_cff7=0.0_r8
309 ad_cff8=0.0_r8
310
311#ifdef VISC_3DCOEF
312 ad_fac1=0.0_r8
313 ad_fac2=0.0_r8
314 ad_visc_p=0.0_r8
315#endif
316
317 ad_dmudz=0.0_r8
318 ad_dnudz=0.0_r8
319 ad_dmvdz=0.0_r8
320 ad_dnvdz=0.0_r8
321
322 ad_ufe(imins:imaxs,jmins:jmaxs)=0.0_r8
323 ad_ufx(imins:imaxs,jmins:jmaxs)=0.0_r8
324 ad_vfe(imins:imaxs,jmins:jmaxs)=0.0_r8
325 ad_vfx(imins:imaxs,jmins:jmaxs)=0.0_r8
326
327 ad_ufse(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
328 ad_ufsx(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
329 ad_vfse(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
330 ad_vfsx(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
331
332 ad_dmude(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
333 ad_dmvde(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
334 ad_dnudx(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
335 ad_dnvdx(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
336
337 ad_dudz(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
338 ad_dvdz(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
339
340 ad_dzde_p(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
341 ad_dzde_r(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
342 ad_dzdx_p(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
343 ad_dzdx_r(imins:imaxs,jmins:jmaxs,1:2)=0.0_r8
344!
345!----------------------------------------------------------------------
346! Compute horizontal harmonic viscosity along geopotential surfaces.
347!----------------------------------------------------------------------
348!
349! Compute basic state horizontal and vertical gradients. Notice the
350! recursive blocking sequence. The vertical placement of the
351! gradients is:
352!
353! dZdx_r, dZde_r, dnUdx, dmVde(:,:,k1) k rho-points
354! dZdx_r, dZde_r, dnUdx, dmVde(:,:,k2) k+1 rho-points
355! dZdx_p, dZde_p, dnVdx, dmUde(:,:,k1) k psi-points
356! dZdx_p, dZde_p, dnVdx, dmUde(:,:,k2) k+1 psi-points
357! UFse, UFsx, dUdz(:,:,k1) k-1/2 WU-points
358! UFse, UFsx, dUdz(:,:,k2) k+1/2 WU-points
359! VFse, VFsx, dVdz(:,:,k1) k-1/2 WV-points
360! VFse, VFsx, dVdz(:,:,k2) k+1/2 WV-points
361!
362! Compute adjoint starting values of k1 and k2. We need to compute
363! the adjoint equivalent of kt=k1, k1=k2, k2=k1.
364!
365 k1=2
366 k2=1
367 DO k=0,n(ng)
368 k1=k2
369 k2=3-k1
370 END DO
371!
372! Compute required BASIC STATE fields. Need to look forward in "kk"
373! index.
374!
375 k_loop : DO k=n(ng),0,-1
376 k2b=1
377 DO kk=0,k
378 k1b=k2b
379 k2b=3-k1b
380 IF (kk.lt.n(ng)) THEN
381!
382! Compute slopes (nondimensional) at RHO- and PSI-points.
383!
384 DO j=jstr-1,jend+1
385 DO i=istru-1,iend+1
386 cff=0.5_r8*(pm(i-1,j)+pm(i,j))
387#ifdef MASKING
388 cff=cff*umask(i,j)
389#endif
390 ufx(i,j)=cff*(z_r(i ,j,kk+1)- &
391 & z_r(i-1,j,kk+1))
392 END DO
393 END DO
394 DO j=jstrv-1,jend+1
395 DO i=istr-1,iend+1
396 cff=0.5_r8*(pn(i,j-1)+pn(i,j))
397#ifdef MASKING
398 cff=cff*vmask(i,j)
399#endif
400 vfe(i,j)=cff*(z_r(i,j ,kk+1)- &
401 & z_r(i,j-1,kk+1))
402 END DO
403 END DO
404!
405 DO j=jstr,jend+1
406 DO i=istr,iend+1
407 dzdx_p(i,j,k2b)=0.5_r8*(ufx(i,j-1)+ &
408 & ufx(i,j ))
409 dzde_p(i,j,k2b)=0.5_r8*(vfe(i-1,j)+ &
410 & vfe(i ,j))
411 END DO
412 END DO
413 DO j=jstrv-1,jend
414 DO i=istru-1,iend
415 dzdx_r(i,j,k2b)=0.5_r8*(ufx(i ,j)+ &
416 & ufx(i+1,j))
417 dzde_r(i,j,k2b)=0.5_r8*(vfe(i,j )+ &
418 & vfe(i,j+1))
419 END DO
420 END DO
421
422 IF (kk.eq.0) THEN
423 DO j=jstr,jend+1
424 DO i=istr,iend+1
425 dzdx_p(i,j,k1b)=0.0_r8
426 dzde_p(i,j,k1b)=0.0_r8
427 END DO
428 END DO
429 DO j=jstrv-1,jend
430 DO i=istru-1,iend
431 dzdx_r(i,j,k1b)=0.0_r8
432 dzde_r(i,j,k1b)=0.0_r8
433 END DO
434 END DO
435 END IF
436!
437! Compute momentum horizontal (1/m/s) and vertical (1/s) gradients.
438!
439 DO j=jstrv-1,jend
440 DO i=istru-1,iend
441 cff=0.5_r8*pm(i,j)
442#ifdef MASKING
443 cff=cff*rmask(i,j)
444#endif
445 dnudx(i,j,k2b)=cff*((pn(i ,j)+pn(i+1,j))* &
446 & u(i+1,j,kk+1,nrhs)- &
447 & (pn(i-1,j)+pn(i ,j))* &
448 & u(i ,j,kk+1,nrhs))
449 END DO
450 END DO
451
452 DO j=jstr,jend+1
453 DO i=istr,iend+1
454 cff=0.125_r8*(pn(i-1,j )+pn(i,j )+ &
455 & pn(i-1,j-1)+pn(i,j-1))
456#ifdef MASKING
457 cff=cff*pmask(i,j)
458#endif
459 dmude(i,j,k2b)=cff*((pm(i-1,j )+pm(i,j ))* &
460 & u(i,j ,kk+1,nrhs)- &
461 & (pm(i-1,j-1)+pm(i,j-1))* &
462 & u(i,j-1,kk+1,nrhs))
463 END DO
464 END DO
465
466 DO j=jstr,jend+1
467 DO i=istr,iend+1
468 cff=0.125_r8*(pm(i-1,j )+pm(i,j )+ &
469 & pm(i-1,j-1)+pm(i,j-1))
470#ifdef MASKING
471 cff=cff*pmask(i,j)
472#endif
473 dnvdx(i,j,k2b)=cff*((pn(i ,j-1)+pn(i ,j))* &
474 & v(i ,j,kk+1,nrhs)- &
475 & (pn(i-1,j-1)+pn(i-1,j))* &
476 & v(i-1,j,kk+1,nrhs))
477 END DO
478 END DO
479
480 DO j=jstrv-1,jend
481 DO i=istru-1,iend
482 cff=0.5_r8*pn(i,j)
483#ifdef MASKING
484 cff=cff*rmask(i,j)
485#endif
486 dmvde(i,j,k2b)=cff*((pm(i,j )+pm(i,j+1))* &
487 & v(i,j+1,kk+1,nrhs)- &
488 & (pm(i,j-1)+pm(i,j ))* &
489 & v(i,j ,kk+1,nrhs))
490 END DO
491 END DO
492
493 IF (kk.eq.0) THEN
494 DO j=jstrv-1,jend
495 DO i=istru-1,iend
496 dnudx(i,j,k1b)=0.0_r8
497 END DO
498 END DO
499 DO j=jstr,jend+1
500 DO i=istr,iend+1
501 dmude(i,j,k1b)=0.0_r8
502 END DO
503 END DO
504 DO j=jstr,jend+1
505 DO i=istr,iend+1
506 dnvdx(i,j,k1b)=0.0_r8
507 END DO
508 END DO
509 DO j=jstrv-1,jend
510 DO i=istru-1,iend
511 dmvde(i,j,k1b)=0.0_r8
512 END DO
513 END DO
514 END IF
515 END IF
516
517 IF ((kk.eq.0).or.(kk.eq.n(ng))) THEN
518 DO j=jstr-1,jend+1
519 DO i=istru-1,iend+1
520 dudz(i,j,k2b)=0.0_r8
521 END DO
522 END DO
523 DO j=jstrv-1,jend+1
524 DO i=istr-1,iend+1
525 dvdz(i,j,k2b)=0.0_r8
526 END DO
527 END DO
528 IF (kk.eq.0) THEN
529 DO j=jstr-1,jend+1
530 DO i=istru-1,iend+1
531 dudz(i,j,k1b)=0.0_r8
532 END DO
533 END DO
534 DO j=jstrv-1,jend+1
535 DO i=istr-1,iend+1
536 dvdz(i,j,k1b)=0.0_r8
537 END DO
538 END DO
539 END IF
540 ELSE
541 DO j=jstr-1,jend+1
542 DO i=istru-1,iend+1
543 cff=1.0_r8/(0.5_r8*(z_r(i-1,j,kk+1)- &
544 & z_r(i-1,j,kk )+ &
545 & z_r(i ,j,kk+1)- &
546 & z_r(i ,j,kk )))
547 dudz(i,j,k2b)=cff*(u(i,j,kk+1,nrhs)- &
548 & u(i,j,kk ,nrhs))
549 END DO
550 END DO
551
552 DO j=jstrv-1,jend+1
553 DO i=istr-1,iend+1
554 cff=1.0_r8/(0.5_r8*(z_r(i,j-1,kk+1)- &
555 & z_r(i,j-1,kk )+ &
556 & z_r(i,j ,kk+1)- &
557 & z_r(i,j ,kk )))
558 dvdz(i,j,k2b)=cff*(v(i,j,kk+1,nrhs)- &
559 & v(i,j,kk ,nrhs))
560 END DO
561 END DO
562 END IF
563 END DO
564!
565! Time-step harmonic, geopotential viscosity term. Notice that
566! momentum at this stage is HzU and HzV and has m2/s units. Add
567! contribution for barotropic forcing terms.
568!
569 above_bottom : IF (k.gt.0) THEN
570 DO j=jstrv,jend
571 DO i=istr,iend
572 cff=dt(ng)*0.25_r8*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
573#ifdef DIAGNOSTICS_UV
574!! DiaV3wrk(i,j,k,M3yvis)=-cff*cff2+dt(ng)*cff4
575!! DiaV3wrk(i,j,k,M3xvis)= cff*cff1+dt(ng)*cff3
576!! DiaV3wrk(i,j,k,M3hvis)=cff5+cff6
577!! DiaRVfrc(i,j,3,M2yvis)=DiaRVfrc(i,j,3,M2yvis)-cff2+cff4
578!! DiaRVfrc(i,j,3,M2xvis)=DiaRVfrc(i,j,3,M2xvis)+cff1+cff3
579!! DiaRVfrc(i,j,3,M2hvis)=DiaRVfrc(i,j,3,M2hvis)+cff1-cff2+ &
580!! & cff3+cff4
581#endif
582!^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)+tl_cff5+tl_cff6
583!^
584 ad_cff5=ad_cff5+ad_v(i,j,k,nnew)
585 ad_cff6=ad_cff6+ad_v(i,j,k,nnew)
586!^ tl_rvfrc(i,j)=tl_rvfrc(i,j)+ &
587!^ & tl_cff1-tl_cff2+tl_cff3+tl_cff4
588!^
589 ad_cff1=ad_cff1+ad_rvfrc(i,j)
590 ad_cff2=ad_cff2-ad_rvfrc(i,j)
591 ad_cff3=ad_cff3+ad_rvfrc(i,j)
592 ad_cff4=ad_cff4+ad_rvfrc(i,j)
593!^ tl_cff6=dt(ng)*(tl_cff3+tl_cff4)
594!^
595 adfac=dt(ng)*ad_cff6
596 ad_cff3=ad_cff3+adfac
597 ad_cff4=ad_cff4+adfac
598 ad_cff6=0.0_r8
599!^ tl_cff5=cff*(tl_cff1-tl_cff2)
600!^
601 adfac=cff*ad_cff5
602 ad_cff1=ad_cff1+adfac
603 ad_cff2=ad_cff2-adfac
604 ad_cff5=0.0_r8
605!^ tl_cff4=tl_VFse(i,j,k2)-tl_VFse(i,j,k1)
606!^
607 ad_vfse(i,j,k1)=ad_vfse(i,j,k1)-ad_cff4
608 ad_vfse(i,j,k2)=ad_vfse(i,j,k2)+ad_cff4
609 ad_cff4=0.0_r8
610!^ tl_cff3=tl_VFsx(i,j,k2)-tl_VFsx(i,j,k1)
611!^
612 ad_vfsx(i,j,k1)=ad_vfsx(i,j,k1)-ad_cff3
613 ad_vfsx(i,j,k2)=ad_vfsx(i,j,k2)+ad_cff3
614 ad_cff3=0.0_r8
615!^ tl_cff2=0.5_r8*(pm(i,j-1)+pm(i,j))* &
616!^ & (tl_VFe(i ,j)-tl_VFe(i,j-1))
617!^
618 adfac=0.5_r8*(pm(i,j-1)+pm(i,j))*ad_cff2
619 ad_vfe(i,j-1)=ad_vfe(i,j-1)-adfac
620 ad_vfe(i,j )=ad_vfe(i,j )+adfac
621 ad_cff2=0.0_r8
622!^ tl_cff1=0.5_r8*(pn(i,j-1)+pn(i,j))* &
623!^ & (tl_VFx(i+1,j)-tl_VFx(i,j ))
624!^
625 adfac=0.5_r8*(pn(i,j-1)+pn(i,j))*ad_cff1
626 ad_vfx(i ,j)=ad_vfx(i ,j)-adfac
627 ad_vfx(i+1,j)=ad_vfx(i+1,j)+adfac
628 ad_cff1=0.0_r8
629 END DO
630 END DO
631
632 DO j=jstr,jend
633 DO i=istru,iend
634 cff=dt(ng)*0.25_r8*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
635#ifdef DIAGNOSTICS_UV
636!! DiaU3wrk(i,j,k,M3yvis)=cff*cff2+dt(ng)*cff4
637!! DiaU3wrk(i,j,k,M3xvis)=cff*cff1+dt(ng)*cff3
638!! DiaU3wrk(i,j,k,M3hvis)=cff5+cff6
639!! DiaRUfrc(i,j,3,M2yvis)=DiaRUfrc(i,j,3,M2yvis)+cff2+cff4
640!! DiaRUfrc(i,j,3,M2xvis)=DiaRUfrc(i,j,3,M2xvis)+cff1+cff3
641!! DiaRUfrc(i,j,3,M2hvis)=DiaRUfrc(i,j,3,M2hvis)+cff1+cff2+ &
642!! & cff3+cff4
643#endif
644!^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)+tl_cff5+tl_cff6
645!^
646 ad_cff5=ad_cff5+ad_u(i,j,k,nnew)
647 ad_cff6=ad_cff6+ad_u(i,j,k,nnew)
648!^ tl_rufrc(i,j)=tl_rufrc(i,j)+ &
649!^ & tl_cff1+tl_cff2+tl_cff3+tl_cff4
650!^
651 ad_cff1=ad_cff1+ad_rufrc(i,j)
652 ad_cff2=ad_cff2+ad_rufrc(i,j)
653 ad_cff3=ad_cff3+ad_rufrc(i,j)
654 ad_cff4=ad_cff4+ad_rufrc(i,j)
655!^ tl_cff6=dt(ng)*(tl_cff3+tl_cff4)
656!^
657 adfac=dt(ng)*ad_cff6
658 ad_cff3=ad_cff3+adfac
659 ad_cff4=ad_cff4+adfac
660 ad_cff6=0.0_r8
661!^ tl_cff5=cff*(tl_cff1+tl_cff2)
662!^
663 adfac=cff*ad_cff5
664 ad_cff1=ad_cff1+adfac
665 ad_cff2=ad_cff2+adfac
666 ad_cff5=0.0_r8
667!^ tl_cff4=tl_UFse(i,j,k2)-tl_UFse(i,j,k1)
668!^
669 ad_ufse(i,j,k1)=ad_ufse(i,j,k1)-ad_cff4
670 ad_ufse(i,j,k2)=ad_ufse(i,j,k2)+ad_cff4
671 ad_cff4=0.0_r8
672!^ tl_cff3=tl_UFsx(i,j,k2)-tl_UFsx(i,j,k1)
673!^
674 ad_ufsx(i,j,k1)=ad_ufsx(i,j,k1)-ad_cff3
675 ad_ufsx(i,j,k2)=ad_ufsx(i,j,k2)+ad_cff3
676 ad_cff3=0.0_r8
677!^ tl_cff2=0.5_r8*(pm(i-1,j)+pm(i,j))* &
678!^ & (tl_UFe(i,j+1)-tl_UFe(i ,j))
679!^
680 adfac=0.5_r8*(pm(i-1,j)+pm(i,j))*ad_cff2
681 ad_ufe(i,j )=ad_ufe(i,j )-adfac
682 ad_ufe(i,j+1)=ad_ufe(i,j+1)+adfac
683 ad_cff2=0.0_r8
684!^ tl_cff1=0.5_r8*(pn(i-1,j)+pn(i,j))* &
685!^ & (tl_UFx(i,j )-tl_UFx(i-1,j))
686!^
687 adfac=0.5_r8*(pn(i-1,j)+pn(i,j))*ad_cff1
688 ad_ufx(i-1,j)=ad_ufx(i-1,j)-adfac
689 ad_ufx(i ,j)=ad_ufx(i ,j)+adfac
690 ad_cff1=0.0_r8
691 END DO
692 END DO
693!
694! Compute adjoint vertical flux (m2/s2) due to sloping
695! terrain-following surfaces.
696!
697 below_surface : IF (k.lt.n(ng)) THEN
698 DO j=jstrv,jend
699 DO i=istr,iend
700#ifdef VISC_3DCOEF
701 cff=0.125_r8* &
702 & (visc3d_r(i,j-1,k )+visc3d_r(i,j,k )+ &
703 & visc3d_r(i,j-1,k+1)+visc3d_r(i,j,k+1))
704 fac1=cff*on_v(i,j)
705 fac2=cff*om_v(i,j)
706#else
707 cff=0.25_r8*(visc2_r(i,j-1)+visc2_r(i,j))
708 fac1=cff*on_v(i,j)
709 fac2=cff*om_v(i,j)
710#endif
711 cff=0.5_r8*(pn(i,j-1)+pn(i,j))
712 dnudz=cff*0.25_r8*(dudz(i ,j ,k2)+ &
713 & dudz(i+1,j ,k2)+ &
714 & dudz(i ,j-1,k2)+ &
715 & dudz(i+1,j-1,k2))
716 dnvdz=cff*dvdz(i,j,k2)
717 cff=0.5_r8*(pm(i,j-1)+pm(i,j))
718 dmudz=cff*0.25_r8*(dudz(i ,j ,k2)+ &
719 & dudz(i+1,j ,k2)+ &
720 & dudz(i ,j-1,k2)+ &
721 & dudz(i+1,j-1,k2))
722 dmvdz=cff*dvdz(i,j,k2)
723!
724 cff1=min(dzdx_p(i ,j,k1),0.0_r8)
725 cff2=min(dzdx_p(i+1,j,k2),0.0_r8)
726 cff3=max(dzdx_p(i ,j,k2),0.0_r8)
727 cff4=max(dzdx_p(i+1,j,k1),0.0_r8)
728 cff5=min(dzde_p(i ,j,k1),0.0_r8)
729 cff6=min(dzde_p(i+1,j,k2),0.0_r8)
730 cff7=max(dzde_p(i ,j,k2),0.0_r8)
731 cff8=max(dzde_p(i+1,j,k1),0.0_r8)
732#ifdef VISC_3DCOEF
733!^ tl_VFse(i,j,k2)=tl_VFse(i,j,k2)+ &
734!^ & tl_fac2* &
735!^ & (cff1*(cff5*dmUdz-dmUde(i ,j,k1))+ &
736!^ & cff2*(cff6*dmUdz-dmUde(i+1,j,k2))+ &
737!^ & cff3*(cff7*dmUdz-dmUde(i ,j,k2))+ &
738!^ & cff4*(cff8*dmUdz-dmUde(i+1,j,k1)))
739!^
740 ad_fac2=ad_fac2+ &
741 & (cff1*(cff5*dmudz-dmude(i ,j,k1))+ &
742 & cff2*(cff6*dmudz-dmude(i+1,j,k2))+ &
743 & cff3*(cff7*dmudz-dmude(i ,j,k2))+ &
744 & cff4*(cff8*dmudz-dmude(i+1,j,k1)))* &
745 & ad_vfse(i,j,k2)
746#endif
747!^ tl_VFse(i,j,k2)=tl_VFse(i,j,k2)+ &
748!^ & fac2* &
749!^ & (tl_cff1*(cff5*dmUdz-dmUde(i ,j,k1))+ &
750!^ & tl_cff2*(cff6*dmUdz-dmUde(i+1,j,k2))+ &
751!^ & tl_cff3*(cff7*dmUdz-dmUde(i ,j,k2))+ &
752!^ & tl_cff4*(cff8*dmUdz-dmUde(i+1,j,k1))+ &
753!^ & cff1*(tl_cff5*dmUdz+cff5*tl_dmUdz- &
754!^ & tl_dmUde(i ,j,k1))+ &
755!^ & cff2*(tl_cff6*dmUdz+cff6*tl_dmUdz- &
756!^ & tl_dmUde(i+1,j,k2))+ &
757!^ & cff3*(tl_cff7*dmUdz+cff7*tl_dmUdz- &
758!^ & tl_dmUde(i ,j,k2))+ &
759!^ & cff4*(tl_cff8*dmUdz+cff8*tl_dmUdz- &
760!^ & tl_dmUde(i+1,j,k1)))
761!^
762 adfac=fac2*ad_vfse(i,j,k2)
763 adfac1=adfac*dmudz
764 ad_cff1=ad_cff1+(cff5*dmudz-dmude(i ,j,k1))*adfac
765 ad_cff2=ad_cff2+(cff6*dmudz-dmude(i+1,j,k2))*adfac
766 ad_cff3=ad_cff3+(cff7*dmudz-dmude(i ,j,k2))*adfac
767 ad_cff4=ad_cff4+(cff8*dmudz-dmude(i+1,j,k1))*adfac
768 ad_cff5=ad_cff5+cff1*adfac1
769 ad_cff6=ad_cff6+cff2*adfac1
770 ad_cff7=ad_cff7+cff3*adfac1
771 ad_cff8=ad_cff8+cff4*adfac1
772 ad_dmudz=ad_dmudz+ &
773 & (cff1*cff5+cff2*cff6+cff3*cff7+cff4*cff8)* &
774 & adfac
775 ad_dmude(i ,j,k1)=ad_dmude(i ,j,k1)-cff1*adfac
776 ad_dmude(i+1,j,k2)=ad_dmude(i+1,j,k2)-cff2*adfac
777 ad_dmude(i ,j,k2)=ad_dmude(i ,j,k2)-cff3*adfac
778 ad_dmude(i+1,j,k1)=ad_dmude(i+1,j,k1)-cff4*adfac
779!^ tl_cff8=(0.5_r8+SIGN(0.5_r8, dZde_p(i+1,j,k1)))* &
780!^ & tl_dZde_p(i+1,j,k1)
781!^
782 ad_dzde_p(i+1,j,k1)=ad_dzde_p(i+1,j,k1)+ &
783 & (0.5_r8+ &
784 & sign(0.5_r8, dzde_p(i+1,j,k1)))* &
785 & ad_cff8
786 ad_cff8=0.0_r8
787!^ tl_cff7=(0.5_r8+SIGN(0.5_r8, dZde_p(i ,j,k2)))* &
788!^ & tl_dZde_p(i ,j,k2)
789!^
790 ad_dzde_p(i ,j,k2)=ad_dzde_p(i ,j,k2)+ &
791 & (0.5_r8+ &
792 & sign(0.5_r8, dzde_p(i ,j,k2)))* &
793 & ad_cff7
794 ad_cff7=0.0_r8
795!^ tl_cff6=(0.5_r8+SIGN(0.5_r8,-dZde_p(i+1,j,k2)))* &
796!^ & tl_dZde_p(i+1,j,k2)
797!^
798 ad_dzde_p(i+1,j,k2)=ad_dzde_p(i+1,j,k2)+ &
799 & (0.5_r8+ &
800 & sign(0.5_r8,-dzde_p(i+1,j,k2)))* &
801 & ad_cff6
802 ad_cff6=0.0_r8
803!^ tl_cff5=(0.5_r8+SIGN(0.5_r8,-dZde_p(i ,j,k1)))* &
804!^ & tl_dZde_p(i ,j,k1)
805!^
806 ad_dzde_p(i ,j,k1)=ad_dzde_p(i ,j,k1)+ &
807 & (0.5_r8+ &
808 & sign(0.5_r8,-dzde_p(i ,j,k1)))* &
809 & ad_cff5
810 ad_cff5=0.0_r8
811!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZdx_p(i+1,j,k1)))* &
812!^ & tl_dZdx_p(i+1,j,k1)
813!^
814 ad_dzdx_p(i+1,j,k1)=ad_dzdx_p(i+1,j,k1)+ &
815 & (0.5_r8+ &
816 & sign(0.5_r8, dzdx_p(i+1,j,k1)))* &
817 & ad_cff4
818 ad_cff4=0.0_r8
819!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZdx_p(i ,j,k2)))* &
820!^ & tl_dZdx_p(i ,j,k2)
821!^
822 ad_dzdx_p(i ,j,k2)=ad_dzdx_p(i ,j,k2)+ &
823 & (0.5_r8+ &
824 & sign(0.5_r8, dzdx_p(i ,j,k2)))* &
825 & ad_cff3
826 ad_cff3=0.0_r8
827!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZdx_p(i+1,j,k2)))* &
828!^ & tl_dZdx_p(i+1,j,k2)
829!^
830 ad_dzdx_p(i+1,j,k2)=ad_dzdx_p(i+1,j,k2)+ &
831 & (0.5_r8+ &
832 & sign(0.5_r8,-dzdx_p(i+1,j,k2)))* &
833 & ad_cff2
834 ad_cff2=0.0_r8
835!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZdx_p(i ,j,k1)))* &
836!^ & tl_dZdx_p(i ,j,k1)
837!^
838 ad_dzdx_p(i ,j,k1)=ad_dzdx_p(i ,j,k1)+ &
839 & (0.5_r8+ &
840 & sign(0.5_r8,-dzdx_p(i ,j,k1)))* &
841 & ad_cff1
842 ad_cff1=0.0_r8
843!
844 cff1=min(dzde_r(i,j-1,k1),0.0_r8)
845 cff2=min(dzde_r(i,j ,k2),0.0_r8)
846 cff3=max(dzde_r(i,j-1,k2),0.0_r8)
847 cff4=max(dzde_r(i,j ,k1),0.0_r8)
848 cff5=min(dzdx_r(i,j-1,k1),0.0_r8)
849 cff6=min(dzdx_r(i,j ,k2),0.0_r8)
850 cff7=max(dzdx_r(i,j-1,k2),0.0_r8)
851 cff8=max(dzdx_r(i,j ,k1),0.0_r8)
852#ifdef VISC_3DCOEF
853!^ tl_VFsx(i,j,k2)=tl_VFsx(i,j,k2)- &
854!^ & tl_fac1* &
855!^ & (cff1*(cff5*dnUdz-dnUdx(i,j-1,k1))+ &
856!^ & cff2*(cff6*dnUdz-dnUdx(i,j ,k2))+ &
857!^ & cff3*(cff7*dnUdz-dnUdx(i,j-1,k2))+ &
858!^ & cff4*(cff8*dnUdz-dnUdx(i,j ,k1)))
859!^
860 ad_fac1=ad_fac1- &
861 & (cff1*(cff5*dnudz-dnudx(i,j-1,k1))+ &
862 & cff2*(cff6*dnudz-dnudx(i,j ,k2))+ &
863 & cff3*(cff7*dnudz-dnudx(i,j-1,k2))+ &
864 & cff4*(cff8*dnudz-dnudx(i,j ,k1)))* &
865 & ad_vfsx(i,j,k2)
866#endif
867!^ tl_VFsx(i,j,k2)=tl_VFsx(i,j,k2)- &
868!^ & fac1* &
869!^ & (tl_cff1*(cff5*dnUdz-dnUdx(i,j-1,k1))+ &
870!^ & tl_cff2*(cff6*dnUdz-dnUdx(i,j ,k2))+ &
871!^ & tl_cff3*(cff7*dnUdz-dnUdx(i,j-1,k2))+ &
872!^ & tl_cff4*(cff8*dnUdz-dnUdx(i,j ,k1))+ &
873!^ & cff1*(tl_cff5*dnUdz+cff5*tl_dnUdz- &
874!^ & tl_dnUdx(i,j-1,k1))+ &
875!^ & cff2*(tl_cff6*dnUdz+cff6*tl_dnUdz- &
876!^ & tl_dnUdx(i,j ,k2))+ &
877!^ & cff3*(tl_cff7*dnUdz+cff7*tl_dnUdz- &
878!^ & tl_dnUdx(i,j-1,k2))+ &
879!^ & cff4*(tl_cff8*dnUdz+cff8*tl_dnUdz- &
880!^ & tl_dnUdx(i,j ,k1)))
881!^
882 adfac=fac1*ad_vfsx(i,j,k2)
883 adfac1=adfac*dnudz
884 ad_cff1=ad_cff1-(cff5*dnudz-dnudx(i,j-1,k1))*adfac
885 ad_cff2=ad_cff2-(cff6*dnudz-dnudx(i,j ,k2))*adfac
886 ad_cff3=ad_cff3-(cff7*dnudz-dnudx(i,j-1,k2))*adfac
887 ad_cff4=ad_cff4-(cff8*dnudz-dnudx(i,j ,k1))*adfac
888 ad_cff5=ad_cff5-cff1*adfac1
889 ad_cff6=ad_cff6-cff2*adfac1
890 ad_cff7=ad_cff7-cff3*adfac1
891 ad_cff8=ad_cff8-cff4*adfac1
892 ad_dnudz=ad_dnudz- &
893 & (cff1*cff5+cff2*cff6+cff3*cff7+cff4*cff8)* &
894 & adfac
895 ad_dnudx(i,j-1,k1)=ad_dnudx(i,j-1,k1)+cff1*adfac
896 ad_dnudx(i,j ,k2)=ad_dnudx(i,j ,k2)+cff2*adfac
897 ad_dnudx(i,j-1,k2)=ad_dnudx(i,j-1,k2)+cff3*adfac
898 ad_dnudx(i,j ,k1)=ad_dnudx(i,j ,k1)+cff4*adfac
899!^ tl_cff8=(0.5_r8+SIGN(0.5_r8, dZdx_r(i,j ,k1)))* &
900!^ & tl_dZdx_r(i,j ,k1)
901!^
902 ad_dzdx_r(i,j ,k1)=ad_dzdx_r(i,j ,k1)+ &
903 & (0.5_r8+ &
904 & sign(0.5_r8, dzdx_r(i,j ,k1)))* &
905 & ad_cff8
906 ad_cff8=0.0_r8
907!^ tl_cff7=(0.5_r8+SIGN(0.5_r8, dZdx_r(i,j-1,k2)))* &
908!^ & tl_dZdx_r(i,j-1,k2)
909!^
910 ad_dzdx_r(i,j-1,k2)=ad_dzdx_r(i,j-1,k2)+ &
911 & (0.5_r8+ &
912 & sign(0.5_r8, dzdx_r(i,j-1,k2)))* &
913 & ad_cff7
914 ad_cff7=0.0_r8
915!^ tl_cff6=(0.5_r8+SIGN(0.5_r8,-dZdx_r(i,j ,k2)))* &
916!^ & tl_dZdx_r(i,j ,k2)
917!^
918 ad_dzdx_r(i,j ,k2)=ad_dzdx_r(i,j ,k2)+ &
919 & (0.5_r8+ &
920 & sign(0.5_r8,-dzdx_r(i,j ,k2)))* &
921 & ad_cff6
922 ad_cff6=0.0_r8
923!^ tl_cff5=(0.5_r8+SIGN(0.5_r8,-dZdx_r(i,j-1,k1)))* &
924!^ & tl_dZdx_r(i,j-1,k1)
925!^
926 ad_dzdx_r(i,j-1,k1)=ad_dzdx_r(i,j-1,k1)+ &
927 & (0.5_r8+ &
928 & sign(0.5_r8,-dzdx_r(i,j-1,k1)))* &
929 & ad_cff5
930 ad_cff5=0.0_r8
931!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZde_r(i,j ,k1)))* &
932!^ & tl_dZde_r(i,j ,k1)
933!^
934 ad_dzde_r(i,j ,k1)=ad_dzde_r(i,j ,k1)+ &
935 & (0.5_r8+ &
936 & sign(0.5_r8, dzde_r(i,j ,k1)))* &
937 & ad_cff4
938 ad_cff4=0.0_r8
939!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZde_r(i,j-1,k2)))* &
940!^ & tl_dZde_r(i,j-1,k2)
941!^
942 ad_dzde_r(i,j-1,k2)=ad_dzde_r(i,j-1,k2)+ &
943 & (0.5_r8+ &
944 & sign(0.5_r8, dzde_r(i,j-1,k2)))* &
945 & ad_cff3
946 ad_cff3=0.0_r8
947!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZde_r(i,j ,k2)))* &
948!^ & tl_dZde_r(i,j ,k2)
949!^
950 ad_dzde_r(i,j ,k2)=ad_dzde_r(i,j ,k2)+ &
951 & (0.5_r8+ &
952 & sign(0.5_r8,-dzde_r(i,j ,k2)))* &
953 & ad_cff2
954 ad_cff2=0.0_r8
955!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZde_r(i,j-1,k1)))* &
956!^ & tl_dZde_r(i,j-1,k1)
957!^
958 ad_dzde_r(i,j-1,k1)=ad_dzde_r(i,j-1,k1)+ &
959 & (0.5_r8+ &
960 & sign(0.5_r8,-dzde_r(i,j-1,k1)))* &
961 & ad_cff1
962 ad_cff1=0.0_r8
963!
964 cff1=min(dzde_r(i,j-1,k1),0.0_r8)
965 cff2=min(dzde_r(i,j ,k2),0.0_r8)
966 cff3=max(dzde_r(i,j-1,k2),0.0_r8)
967 cff4=max(dzde_r(i,j ,k1),0.0_r8)
968#ifdef VISC_3DCOEF
969!^ tl_VFse(i,j,k2)=tl_VFse(i,j,k2)+ &
970!^ & tl_fac2* &
971!^ & (cff1*(cff1*dmVdz-dmVde(i,j-1,k1))+ &
972!^ & cff2*(cff2*dmVdz-dmVde(i,j ,k2))+ &
973!^ & cff3*(cff3*dmVdz-dmVde(i,j-1,k2))+ &
974!^ & cff4*(cff4*dmVdz-dmVde(i,j ,k1)))
975!^
976 ad_fac2=ad_fac2+ &
977 & (cff1*(cff1*dmvdz-dmvde(i,j-1,k1))+ &
978 & cff2*(cff2*dmvdz-dmvde(i,j ,k2))+ &
979 & cff3*(cff3*dmvdz-dmvde(i,j-1,k2))+ &
980 & cff4*(cff4*dmvdz-dmvde(i,j ,k1)))* &
981 & ad_vfse(i,j,k2)
982#endif
983!^ tl_VFse(i,j,k2)=fac2* &
984!^ & (tl_cff1*(cff1*dmVdz-dmVde(i,j-1,k1))+ &
985!^ & tl_cff2*(cff2*dmVdz-dmVde(i,j ,k2))+ &
986!^ & tl_cff3*(cff3*dmVdz-dmVde(i,j-1,k2))+ &
987!^ & tl_cff4*(cff4*dmVdz-dmVde(i,j ,k1))+ &
988!^ & cff1*(tl_cff1*dmVdz+cff1*tl_dmVdz- &
989!^ & tl_dmVde(i,j-1,k1))+ &
990!^ & cff2*(tl_cff2*dmVdz+cff2*tl_dmVdz- &
991!^ & tl_dmVde(i,j ,k2))+ &
992!^ & cff3*(tl_cff3*dmVdz+cff3*tl_dmVdz- &
993!^ & tl_dmVde(i,j-1,k2))+ &
994!^ & cff4*(tl_cff4*dmVdz+cff4*tl_dmVdz- &
995!^ & tl_dmVde(i,j ,k1)))
996!^
997 cff=2.0_r8*dmvdz
998 adfac=fac2*ad_vfse(i,j,k2)
999 ad_cff1=ad_cff1+(cff1*cff-dmvde(i,j-1,k1))*adfac
1000 ad_cff2=ad_cff2+(cff2*cff-dmvde(i,j ,k2))*adfac
1001 ad_cff3=ad_cff3+(cff3*cff-dmvde(i,j-1,k2))*adfac
1002 ad_cff4=ad_cff4+(cff4*cff-dmvde(i,j ,k1))*adfac
1003 ad_dmvdz=ad_dmvdz+ &
1004 & (cff1*cff1+cff2*cff2+cff3*cff3+cff4*cff4)* &
1005 & adfac
1006 ad_dmvde(i,j-1,k1)=ad_dmvde(i,j-1,k1)-cff1*adfac
1007 ad_dmvde(i,j ,k2)=ad_dmvde(i,j ,k2)-cff2*adfac
1008 ad_dmvde(i,j-1,k2)=ad_dmvde(i,j-1,k2)-cff3*adfac
1009 ad_dmvde(i,j ,k1)=ad_dmvde(i,j ,k1)-cff4*adfac
1010 ad_vfse(i,j,k2)=0.0_r8
1011!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZde_r(i,j ,k1)))* &
1012!^ & tl_dZde_r(i,j ,k1)
1013!^
1014 ad_dzde_r(i,j ,k1)=ad_dzde_r(i,j ,k1)+ &
1015 & (0.5_r8+ &
1016 & sign(0.5_r8, dzde_r(i,j ,k1)))* &
1017 & ad_cff4
1018 ad_cff4=0.0_r8
1019!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZde_r(i,j-1,k2)))* &
1020!^ & tl_dZde_r(i,j-1,k2)
1021!^
1022 ad_dzde_r(i,j-1,k2)=ad_dzde_r(i,j-1,k2)+ &
1023 & (0.5_r8+ &
1024 & sign(0.5_r8, dzde_r(i,j-1,k2)))* &
1025 & ad_cff3
1026 ad_cff3=0.0_r8
1027!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZde_r(i,j ,k2)))* &
1028!^ & tl_dZde_r(i,j ,k2)
1029!^
1030 ad_dzde_r(i,j ,k2)=ad_dzde_r(i,j ,k2)+ &
1031 & (0.5_r8+ &
1032 & sign(0.5_r8,-dzde_r(i,j ,k2)))* &
1033 & ad_cff2
1034 ad_cff2=0.0_r8
1035!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZde_r(i,j-1,k1)))* &
1036!^ & tl_dZde_r(i,j-1,k1)
1037!^
1038 ad_dzde_r(i,j-1,k1)=ad_dzde_r(i,j-1,k1)+ &
1039 & (0.5_r8+ &
1040 & sign(0.5_r8,-dzde_r(i,j-1,k1)))* &
1041 & ad_cff1
1042 ad_cff1=0.0_r8
1043!
1044 cff1=min(dzdx_p(i ,j,k1),0.0_r8)
1045 cff2=min(dzdx_p(i+1,j,k2),0.0_r8)
1046 cff3=max(dzdx_p(i ,j,k2),0.0_r8)
1047 cff4=max(dzdx_p(i+1,j,k1),0.0_r8)
1048#ifdef VISC_3DCOEF
1049!^ tl_VFsx(i,j,k2)=tl_VFsx(i,j,k2)+ &
1050!^ & tl_fac1* &
1051!^ & (cff1*(cff1*dnVdz-dnVdx(i ,j,k1))+ &
1052!^ & cff2*(cff2*dnVdz-dnVdx(i+1,j,k2))+ &
1053!^ & cff3*(cff3*dnVdz-dnVdx(i ,j,k2))+ &
1054!^ & cff4*(cff4*dnVdz-dnVdx(i+1,j,k1)))
1055!^
1056 ad_fac1=ad_fac1+ &
1057 & (cff1*(cff1*dnvdz-dnvdx(i ,j,k1))+ &
1058 & cff2*(cff2*dnvdz-dnvdx(i+1,j,k2))+ &
1059 & cff3*(cff3*dnvdz-dnvdx(i ,j,k2))+ &
1060 & cff4*(cff4*dnvdz-dnvdx(i+1,j,k1)))* &
1061 & ad_vfsx(i,j,k2)
1062#endif
1063!^ tl_VFsx(i,j,k2)=fac1* &
1064!^ & (tl_cff1*(cff1*dnVdz-dnVdx(i ,j,k1))+ &
1065!^ & tl_cff2*(cff2*dnVdz-dnVdx(i+1,j,k2))+ &
1066!^ & tl_cff3*(cff3*dnVdz-dnVdx(i ,j,k2))+ &
1067!^ & tl_cff4*(cff4*dnVdz-dnVdx(i+1,j,k1))+ &
1068!^ & cff1*(tl_cff1*dnVdz+cff1*tl_dnVdz- &
1069!^ & tl_dnVdx(i ,j,k1))+ &
1070!^ & cff2*(tl_cff2*dnVdz+cff2*tl_dnVdz- &
1071!^ & tl_dnVdx(i+1,j,k2))+ &
1072!^ & cff3*(tl_cff3*dnVdz+cff3*tl_dnVdz- &
1073!^ & tl_dnVdx(i ,j,k2))+ &
1074!^ & cff4*(tl_cff4*dnVdz+cff4*tl_dnVdz- &
1075!^ & tl_dnVdx(i+1,j,k1)))
1076!^
1077 cff=2.0_r8*dnvdz
1078 adfac=fac1*ad_vfsx(i,j,k2)
1079 ad_cff1=ad_cff1+(cff1*cff-dnvdx(i ,j,k1))*adfac
1080 ad_cff2=ad_cff2+(cff2*cff-dnvdx(i+1,j,k2))*adfac
1081 ad_cff3=ad_cff3+(cff3*cff-dnvdx(i ,j,k2))*adfac
1082 ad_cff4=ad_cff4+(cff4*cff-dnvdx(i+1,j,k1))*adfac
1083 ad_dnvdz=ad_dnvdz+ &
1084 & (cff1*cff1+cff2*cff2+cff3*cff3+cff4*cff4)* &
1085 & adfac
1086 ad_dnvdx(i ,j,k1)=ad_dnvdx(i ,j,k1)-cff1*adfac
1087 ad_dnvdx(i+1,j,k2)=ad_dnvdx(i+1,j,k2)-cff2*adfac
1088 ad_dnvdx(i ,j,k2)=ad_dnvdx(i ,j,k2)-cff3*adfac
1089 ad_dnvdx(i+1,j,k1)=ad_dnvdx(i+1,j,k1)-cff4*adfac
1090 ad_vfsx(i,j,k2)=0.0_r8
1091!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZdx_p(i+1,j,k1)))* &
1092!^ & tl_dZdx_p(i+1,j,k1)
1093!^
1094 ad_dzdx_p(i+1,j,k1)=ad_dzdx_p(i+1,j,k1)+ &
1095 & (0.5_r8+ &
1096 & sign(0.5_r8, dzdx_p(i+1,j,k1)))* &
1097 & ad_cff4
1098 ad_cff4=0.0_r8
1099!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZdx_p(i ,j,k2)))* &
1100!^ & tl_dZdx_p(i ,j,k2)
1101!^
1102 ad_dzdx_p(i ,j,k2)=ad_dzdx_p(i ,j,k2)+ &
1103 & (0.5_r8+ &
1104 & sign(0.5_r8, dzdx_p(i ,j,k2)))* &
1105 & ad_cff3
1106 ad_cff3=0.0_r8
1107!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZdx_p(i+1,j,k2)))* &
1108!^ & tl_dZdx_p(i+1,j,k2)
1109!^
1110 ad_dzdx_p(i+1,j,k2)=ad_dzdx_p(i+1,j,k2)+ &
1111 & (0.5_r8+ &
1112 & sign(0.5_r8,-dzdx_p(i+1,j,k2)))* &
1113 & ad_cff2
1114 ad_cff2=0.0_r8
1115!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZdx_p(i ,j,k1)))* &
1116!^ & tl_dZdx_p(i ,j,k1)
1117!^
1118 ad_dzdx_p(i ,j,k1)=ad_dzdx_p(i ,j,k1)+ &
1119 & (0.5_r8+ &
1120 & sign(0.5_r8,-dzdx_p(i ,j,k1)))* &
1121 & ad_cff1
1122 ad_cff1=0.0_r8
1123!
1124 cff=0.5_r8*(pm(i,j-1)+pm(i,j))
1125!^ tl_dmVdz=cff*tl_dVdz(i,j,k2)
1126!^
1127 ad_dvdz(i,j,k2)=ad_dvdz(i,j,k2)+cff*ad_dmvdz
1128 ad_dmvdz=0.0_r8
1129!^ tl_dmUdz=cff*0.25_r8*(tl_dUdz(i ,j ,k2)+ &
1130!^ & tl_dUdz(i+1,j ,k2)+ &
1131!^ & tl_dUdz(i ,j-1,k2)+ &
1132!^ & tl_dUdz(i+1,j-1,k2))
1133!^
1134 adfac=cff*0.25_r8*ad_dmudz
1135 ad_dudz(i ,j-1,k2)=ad_dudz(i ,j-1,k2)+adfac
1136 ad_dudz(i+1,j-1,k2)=ad_dudz(i+1,j-1,k2)+adfac
1137 ad_dudz(i ,j ,k2)=ad_dudz(i ,j ,k2)+adfac
1138 ad_dudz(i+1,j ,k2)=ad_dudz(i+1,j ,k2)+adfac
1139 ad_dmudz=0.0_r8
1140!
1141 cff=0.5_r8*(pn(i,j-1)+pn(i,j))
1142!^ tl_dnVdz=cff*tl_dVdz(i,j,k2)
1143!^
1144 ad_dvdz(i,j,k2)=ad_dvdz(i,j,k2)+cff*ad_dnvdz
1145 ad_dnvdz=0.0_r8
1146!^ tl_dnUdz=cff*0.25_r8*(tl_dUdz(i ,j ,k2)+ &
1147!^ & tl_dUdz(i+1,j ,k2)+ &
1148!^ & tl_dUdz(i ,j-1,k2)+ &
1149!^ & tl_dUdz(i+1,j-1,k2))
1150!^
1151 adfac=cff*0.25_r8*ad_dnudz
1152 ad_dudz(i ,j-1,k2)=ad_dudz(i ,j-1,k2)+adfac
1153 ad_dudz(i+1,j-1,k2)=ad_dudz(i+1,j-1,k2)+adfac
1154 ad_dudz(i ,j ,k2)=ad_dudz(i ,j ,k2)+adfac
1155 ad_dudz(i+1,j ,k2)=ad_dudz(i+1,j ,k2)+adfac
1156 ad_dnudz=0.0_r8
1157#ifdef VISC_3DCOEF
1158!^ tl_fac2=tl_cff*om_v(i,j)
1159!^ tl_fac1=tl_cff*on_v(i,j)
1160!^
1161 ad_cff=ad_cff+ &
1162 & on_v(i,j)*ad_fac1+om_v(i,j)*ad_fac2
1163 ad_fac1=0.0_r8
1164 ad_fac2=0.0_r8
1165!^ tl_cff=0.125_r8* &
1166!^ & (tl_visc3d_r(i,j-1,k )+tl_visc3d_r(i,j,k )+ &
1167!^ & tl_visc3d_r(i,j-1,k+1)+tl_visc3d_r(i,j,k+1))
1168!^
1169 adfac=0.125_r8*ad_cff
1170 ad_visc3d_r(i,j-1,k )=ad_visc3d_r(i,j-1,k )+adfac
1171 ad_visc3d_r(i,j ,k )=ad_visc3d_r(i,j ,k )+adfac
1172 ad_visc3d_r(i,j-1,k+1)=ad_visc3d_r(i,j-1,k+1)+adfac
1173 ad_visc3d_r(i,j ,k+1)=ad_visc3d_r(i,j ,k+1)+adfac
1174 ad_cff=0.0_r8
1175#endif
1176 END DO
1177 END DO
1178!
1179 DO j=jstr,jend
1180 DO i=istru,iend
1181#ifdef VISC_3DCOEF
1182 cff=0.125_r8* &
1183 & (visc3d_r(i-1,j,k )+visc3d_r(i,j,k )+ &
1184 & visc3d_r(i-1,j,k+1)+visc3d_r(i,j,k+1))
1185 fac1=cff*on_u(i,j)
1186 fac2=cff*om_u(i,j)
1187#else
1188 cff=0.25_r8*(visc2_r(i-1,j)+visc2_r(i,j))
1189 fac1=cff*on_u(i,j)
1190 fac2=cff*om_u(i,j)
1191#endif
1192 cff=0.5_r8*(pn(i-1,j)+pn(i,j))
1193 dnudz=cff*dudz(i,j,k2)
1194 dnvdz=cff*0.25_r8*(dvdz(i-1,j+1,k2)+ &
1195 & dvdz(i ,j+1,k2)+ &
1196 & dvdz(i-1,j ,k2)+ &
1197 & dvdz(i ,j ,k2))
1198 cff=0.5_r8*(pm(i-1,j)+pm(i,j))
1199 dmudz=cff*dudz(i,j,k2)
1200 dmvdz=cff*0.25_r8*(dvdz(i-1,j+1,k2)+ &
1201 & dvdz(i ,j+1,k2)+ &
1202 & dvdz(i-1,j ,k2)+ &
1203 & dvdz(i ,j ,k2))
1204!
1205 cff1=min(dzdx_r(i-1,j,k1),0.0_r8)
1206 cff2=min(dzdx_r(i ,j,k2),0.0_r8)
1207 cff3=max(dzdx_r(i-1,j,k2),0.0_r8)
1208 cff4=max(dzdx_r(i ,j,k1),0.0_r8)
1209 cff5=min(dzde_r(i-1,j,k1),0.0_r8)
1210 cff6=min(dzde_r(i ,j,k2),0.0_r8)
1211 cff7=max(dzde_r(i-1,j,k2),0.0_r8)
1212 cff8=max(dzde_r(i ,j,k1),0.0_r8)
1213#ifdef VISC_3DCOEF
1214!^ tl_UFse(i,j,k2)=tl_UFse(i,j,k2)- &
1215!^ & tl_fac2* &
1216!^ & (cff1*(cff5*dmVdz-dmVde(i-1,j,k1))+ &
1217!^ & cff2*(cff6*dmVdz-dmVde(i ,j,k2))+ &
1218!^ & cff3*(cff7*dmVdz-dmVde(i-1,j,k2))+ &
1219!^ & cff4*(cff8*dmVdz-dmVde(i ,j,k1)))
1220!^
1221 ad_fac2=ad_fac2- &
1222 & (cff1*(cff5*dmvdz-dmvde(i-1,j,k1))+ &
1223 & cff2*(cff6*dmvdz-dmvde(i ,j,k2))+ &
1224 & cff3*(cff7*dmvdz-dmvde(i-1,j,k2))+ &
1225 & cff4*(cff8*dmvdz-dmvde(i ,j,k1)))* &
1226 & ad_ufse(i,j,k2)
1227#endif
1228!^ tl_UFse(i,j,k2)=tl_UFse(i,j,k2)- &
1229!^ & fac2* &
1230!^ & (tl_cff1*(cff5*dmVdz-dmVde(i-1,j,k1))+ &
1231!^ & tl_cff2*(cff6*dmVdz-dmVde(i ,j,k2))+ &
1232!^ & tl_cff3*(cff7*dmVdz-dmVde(i-1,j,k2))+ &
1233!^ & tl_cff4*(cff8*dmVdz-dmVde(i ,j,k1))+ &
1234!^ & cff1*(tl_cff5*dmVdz+cff5*tl_dmVdz- &
1235!^ & tl_dmVde(i-1,j,k1))+ &
1236!^ & cff2*(tl_cff6*dmVdz+cff6*tl_dmVdz- &
1237!^ & tl_dmVde(i ,j,k2))+ &
1238!^ & cff3*(tl_cff7*dmVdz+cff7*tl_dmVdz- &
1239!^ & tl_dmVde(i-1,j,k2))+ &
1240!^ & cff4*(tl_cff8*dmVdz+cff8*tl_dmVdz- &
1241!^ & tl_dmVde(i ,j,k1)))
1242!^
1243 adfac=fac2*ad_ufse(i,j,k2)
1244 adfac1=adfac*dmvdz
1245 ad_cff1=ad_cff1-(cff5*dmvdz-dmvde(i-1,j,k1))*adfac
1246 ad_cff2=ad_cff2-(cff6*dmvdz-dmvde(i ,j,k2))*adfac
1247 ad_cff3=ad_cff3-(cff7*dmvdz-dmvde(i-1,j,k2))*adfac
1248 ad_cff4=ad_cff4-(cff8*dmvdz-dmvde(i ,j,k1))*adfac
1249 ad_cff5=ad_cff5-cff1*adfac1
1250 ad_cff6=ad_cff6-cff2*adfac1
1251 ad_cff7=ad_cff7-cff3*adfac1
1252 ad_cff8=ad_cff8-cff4*adfac1
1253 ad_dmvdz=ad_dmvdz- &
1254 & (cff1*cff5+cff2*cff6+cff3*cff7+cff4*cff8)* &
1255 & adfac
1256 ad_dmvde(i-1,j,k1)=ad_dmvde(i-1,j,k1)+cff1*adfac
1257 ad_dmvde(i ,j,k2)=ad_dmvde(i ,j,k2)+cff2*adfac
1258 ad_dmvde(i-1,j,k2)=ad_dmvde(i-1,j,k2)+cff3*adfac
1259 ad_dmvde(i ,j,k1)=ad_dmvde(i ,j,k1)+cff4*adfac
1260!^ tl_cff8=(0.5_r8+SIGN(0.5_r8, dZde_r(i ,j,k1)))* &
1261!^ & tl_dZde_r(i ,j,k1)
1262!^
1263 ad_dzde_r(i ,j,k1)=ad_dzde_r(i ,j,k1)+ &
1264 & (0.5_r8+ &
1265 & sign(0.5_r8, dzde_r(i ,j,k1)))* &
1266 & ad_cff8
1267 ad_cff8=0.0_r8
1268!^ tl_cff7=(0.5_r8+SIGN(0.5_r8, dZde_r(i-1,j,k2)))* &
1269!^ & tl_dZde_r(i-1,j,k2)
1270!^
1271 ad_dzde_r(i-1,j,k2)=ad_dzde_r(i-1,j,k2)+ &
1272 & (0.5_r8+ &
1273 & sign(0.5_r8, dzde_r(i-1,j,k2)))* &
1274 & ad_cff7
1275 ad_cff7=0.0_r8
1276!^ tl_cff6=(0.5_r8+SIGN(0.5_r8,-dZde_r(i ,j,k2)))* &
1277!^ & tl_dZde_r(i ,j,k2)
1278!^
1279 ad_dzde_r(i ,j,k2)=ad_dzde_r(i ,j,k2)+ &
1280 & (0.5_r8+ &
1281 & sign(0.5_r8,-dzde_r(i ,j,k2)))* &
1282 & ad_cff6
1283 ad_cff6=0.0_r8
1284!^ tl_cff5=(0.5_r8+SIGN(0.5_r8,-dZde_r(i-1,j,k1)))* &
1285!^ & tl_dZde_r(i-1,j,k1)
1286!^
1287 ad_dzde_r(i-1,j,k1)=ad_dzde_r(i-1,j,k1)+ &
1288 & (0.5_r8+ &
1289 & sign(0.5_r8,-dzde_r(i-1,j,k1)))* &
1290 & ad_cff5
1291 ad_cff5=0.0_r8
1292!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZdx_r(i ,j,k1)))* &
1293!^ & tl_dZdx_r(i ,j,k1)
1294!^
1295 ad_dzdx_r(i ,j,k1)=ad_dzdx_r(i ,j,k1)+ &
1296 & (0.5_r8+ &
1297 & sign(0.5_r8, dzdx_r(i ,j,k1)))* &
1298 & ad_cff4
1299 ad_cff4=0.0_r8
1300!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZdx_r(i-1,j,k2)))* &
1301!^ & tl_dZdx_r(i-1,j,k2)
1302!^
1303 ad_dzdx_r(i-1,j,k2)=ad_dzdx_r(i-1,j,k2)+ &
1304 & (0.5_r8+ &
1305 & sign(0.5_r8, dzdx_r(i-1,j,k2)))* &
1306 & ad_cff3
1307 ad_cff3=0.0_r8
1308!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZdx_r(i ,j,k2)))* &
1309!^ & tl_dZdx_r(i ,j,k2)
1310!^
1311 ad_dzdx_r(i ,j,k2)=ad_dzdx_r(i ,j,k2)+ &
1312 & (0.5_r8+ &
1313 & sign(0.5_r8,-dzdx_r(i ,j,k2)))* &
1314 & ad_cff2
1315 ad_cff2=0.0_r8
1316!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZdx_r(i-1,j,k1)))* &
1317!^ & tl_dZdx_r(i-1,j,k1)
1318!^
1319 ad_dzdx_r(i-1,j,k1)=ad_dzdx_r(i-1,j,k1)+ &
1320 & (0.5_r8+ &
1321 & sign(0.5_r8,-dzdx_r(i-1,j,k1)))* &
1322 & ad_cff1
1323 ad_cff1=0.0_r8
1324!
1325 cff1=min(dzde_p(i,j ,k1),0.0_r8)
1326 cff2=min(dzde_p(i,j+1,k2),0.0_r8)
1327 cff3=max(dzde_p(i,j ,k2),0.0_r8)
1328 cff4=max(dzde_p(i,j+1,k1),0.0_r8)
1329 cff5=min(dzdx_p(i,j ,k1),0.0_r8)
1330 cff6=min(dzdx_p(i,j+1,k2),0.0_r8)
1331 cff7=max(dzdx_p(i,j ,k2),0.0_r8)
1332 cff8=max(dzdx_p(i,j+1,k1),0.0_r8)
1333#ifdef VISC_3DCOEF
1334!^ tl_UFsx(i,j,k2)=tl_UFsx(i,j,k2)+ &
1335!^ & tl_fac1* &
1336!^ & (cff1*(cff5*dnVdz-dnVdx(i,j ,k1))+ &
1337!^ & cff2*(cff6*dnVdz-dnVdx(i,j+1,k2))+ &
1338!^ & cff3*(cff7*dnVdz-dnVdx(i,j ,k2))+ &
1339!^ & cff4*(cff8*dnVdz-dnVdx(i,j+1,k1)))
1340!^
1341 ad_fac1=ad_fac1+ &
1342 & (cff1*(cff5*dnvdz-dnvdx(i,j ,k1))+ &
1343 & cff2*(cff6*dnvdz-dnvdx(i,j+1,k2))+ &
1344 & cff3*(cff7*dnvdz-dnvdx(i,j ,k2))+ &
1345 & cff4*(cff8*dnvdz-dnvdx(i,j+1,k1)))* &
1346 & ad_ufsx(i,j,k2)
1347#endif
1348!^ tl_UFsx(i,j,k2)=tl_UFsx(i,j,k2)+ &
1349!^ & fac1* &
1350!^ & (tl_cff1*(cff5*dnVdz-dnVdx(i,j ,k1))+ &
1351!^ & tl_cff2*(cff6*dnVdz-dnVdx(i,j+1,k2))+ &
1352!^ & tl_cff3*(cff7*dnVdz-dnVdx(i,j ,k2))+ &
1353!^ & tl_cff4*(cff8*dnVdz-dnVdx(i,j+1,k1))+ &
1354!^ & cff1*(tl_cff5*dnVdz+cff5*tl_dnVdz- &
1355!^ & tl_dnVdx(i,j ,k1))+ &
1356!^ & cff2*(tl_cff6*dnVdz+cff6*tl_dnVdz- &
1357!^ & tl_dnVdx(i,j+1,k2))+ &
1358!^ & cff3*(tl_cff7*dnVdz+cff7*tl_dnVdz- &
1359!^ & tl_dnVdx(i,j ,k2))+ &
1360!^ & cff4*(tl_cff8*dnVdz+cff8*tl_dnVdz- &
1361!^ & tl_dnVdx(i,j+1,k1)))
1362!^
1363 adfac=fac1*ad_ufsx(i,j,k2)
1364 adfac1=adfac*dnvdz
1365 ad_cff1=ad_cff1+(cff5*dnvdz-dnvdx(i,j ,k1))*adfac
1366 ad_cff2=ad_cff2+(cff6*dnvdz-dnvdx(i,j+1,k2))*adfac
1367 ad_cff3=ad_cff3+(cff7*dnvdz-dnvdx(i,j ,k2))*adfac
1368 ad_cff4=ad_cff4+(cff8*dnvdz-dnvdx(i,j+1,k1))*adfac
1369 ad_cff5=ad_cff5+cff1*adfac1
1370 ad_cff6=ad_cff6+cff2*adfac1
1371 ad_cff7=ad_cff7+cff3*adfac1
1372 ad_cff8=ad_cff8+cff4*adfac1
1373 ad_dnvdz=ad_dnvdz+ &
1374 & (cff1*cff5+cff2*cff6+cff3*cff7+cff4*cff8)* &
1375 & adfac
1376 ad_dnvdx(i,j ,k1)=ad_dnvdx(i,j ,k1)-cff1*adfac
1377 ad_dnvdx(i,j+1,k2)=ad_dnvdx(i,j+1,k2)-cff2*adfac
1378 ad_dnvdx(i,j ,k2)=ad_dnvdx(i,j ,k2)-cff3*adfac
1379 ad_dnvdx(i,j+1,k1)=ad_dnvdx(i,j+1,k1)-cff4*adfac
1380!^ tl_cff8=(0.5_r8+SIGN(0.5_r8, dZdx_p(i,j+1,k1)))* &
1381!^ & tl_dZdx_p(i,j+1,k1)
1382!^
1383 ad_dzdx_p(i,j+1,k1)=ad_dzdx_p(i,j+1,k1)+ &
1384 & (0.5_r8+ &
1385 & sign(0.5_r8, dzdx_p(i,j+1,k1)))* &
1386 & ad_cff8
1387 ad_cff8=0.0_r8
1388!^ tl_cff7=(0.5_r8+SIGN(0.5_r8, dZdx_p(i,j ,k2)))* &
1389!^ & tl_dZdx_p(i,j ,k2)
1390!^
1391 ad_dzdx_p(i,j ,k2)=ad_dzdx_p(i,j ,k2)+ &
1392 & (0.5_r8+ &
1393 & sign(0.5_r8, dzdx_p(i,j ,k2)))* &
1394 & ad_cff7
1395 ad_cff7=0.0_r8
1396!^ tl_cff6=(0.5_r8+SIGN(0.5_r8,-dZdx_p(i,j+1,k2)))* &
1397!^ & tl_dZdx_p(i,j+1,k2)
1398!^
1399 ad_dzdx_p(i,j+1,k2)=ad_dzdx_p(i,j+1,k2)+ &
1400 & (0.5_r8+ &
1401 & sign(0.5_r8,-dzdx_p(i,j+1,k2)))* &
1402 & ad_cff6
1403 ad_cff6=0.0_r8
1404!^ tl_cff5=(0.5_r8+SIGN(0.5_r8,-dZdx_p(i,j ,k1)))* &
1405!^ & tl_dZdx_p(i,j ,k1)
1406!^
1407 ad_dzdx_p(i,j ,k1)=ad_dzdx_p(i,j ,k1)+ &
1408 & (0.5_r8+ &
1409 & sign(0.5_r8,-dzdx_p(i,j ,k1)))* &
1410 & ad_cff5
1411 ad_cff5=0.0_r8
1412!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZde_p(i,j+1,k1)))* &
1413!^ & tl_dZde_p(i,j+1,k1)
1414!^
1415 ad_dzde_p(i,j+1,k1)=ad_dzde_p(i,j+1,k1)+ &
1416 & (0.5_r8+ &
1417 & sign(0.5_r8, dzde_p(i,j+1,k1)))* &
1418 & ad_cff4
1419 ad_cff4=0.0_r8
1420!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZde_p(i,j ,k2)))* &
1421!^ & tl_dZde_p(i,j ,k2)
1422!^
1423 ad_dzde_p(i,j ,k2)=ad_dzde_p(i,j ,k2)+ &
1424 & (0.5_r8+ &
1425 & sign(0.5_r8, dzde_p(i,j ,k2)))* &
1426 & ad_cff3
1427 ad_cff3=0.0_r8
1428!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZde_p(i,j+1,k2)))* &
1429!^ & tl_dZde_p(i,j+1,k2)
1430!^
1431 ad_dzde_p(i,j+1,k2)=ad_dzde_p(i,j+1,k2)+ &
1432 & (0.5_r8+ &
1433 & sign(0.5_r8,-dzde_p(i,j+1,k2)))* &
1434 & ad_cff2
1435 ad_cff2=0.0_r8
1436!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZde_p(i,j ,k1)))* &
1437!^ & tl_dZde_p(i,j ,k1)
1438!^
1439 ad_dzde_p(i,j ,k1)=ad_dzde_p(i,j ,k1)+ &
1440 & (0.5_r8+ &
1441 & sign(0.5_r8,-dzde_p(i,j ,k1)))* &
1442 & ad_cff1
1443 ad_cff1=0.0_r8
1444!
1445 cff1=min(dzde_p(i,j ,k1),0.0_r8)
1446 cff2=min(dzde_p(i,j+1,k2),0.0_r8)
1447 cff3=max(dzde_p(i,j ,k2),0.0_r8)
1448 cff4=max(dzde_p(i,j+1,k1),0.0_r8)
1449#ifdef VISC_3DCOEF
1450!^ tl_UFse(i,j,k2)=tl_UFse(i,j,k2)+
1451!^ & tl_fac2* &
1452!^ & (cff1*(cff1*dmUdz-dmUde(i,j ,k1))+ &
1453!^ & cff2*(cff2*dmUdz-dmUde(i,j+1,k2))+ &
1454!^ & cff3*(cff3*dmUdz-dmUde(i,j ,k2))+ &
1455!^ & cff4*(cff4*dmUdz-dmUde(i,j+1,k1)))
1456!^
1457 ad_fac2=ad_fac2+ &
1458 & (cff1*(cff1*dmudz-dmude(i,j ,k1))+ &
1459 & cff2*(cff2*dmudz-dmude(i,j+1,k2))+ &
1460 & cff3*(cff3*dmudz-dmude(i,j ,k2))+ &
1461 & cff4*(cff4*dmudz-dmude(i,j+1,k1)))* &
1462 & ad_ufse(i,j,k2)
1463#endif
1464!^ tl_UFse(i,j,k2)=fac2* &
1465!^ & (tl_cff1*(cff1*dmUdz-dmUde(i,j ,k1))+ &
1466!^ & tl_cff2*(cff2*dmUdz-dmUde(i,j+1,k2))+ &
1467!^ & tl_cff3*(cff3*dmUdz-dmUde(i,j ,k2))+ &
1468!^ & tl_cff4*(cff4*dmUdz-dmUde(i,j+1,k1))+ &
1469!^ & cff1*(tl_cff1*dmUdz+cff1*tl_dmUdz- &
1470!^ & tl_dmUde(i,j ,k1))+ &
1471!^ & cff2*(tl_cff2*dmUdz+cff2*tl_dmUdz- &
1472!^ & tl_dmUde(i,j+1,k2))+ &
1473!^ & cff3*(tl_cff3*dmUdz+cff3*tl_dmUdz- &
1474!^ & tl_dmUde(i,j ,k2))+ &
1475!^ & cff4*(tl_cff4*dmUdz+cff4*tl_dmUdz- &
1476!^ & tl_dmUde(i,j+1,k1)))
1477!^
1478 cff=2.0_r8*dmudz
1479 adfac=fac2*ad_ufse(i,j,k2)
1480 ad_cff1=ad_cff1+(cff1*cff-dmude(i,j ,k1))*adfac
1481 ad_cff2=ad_cff2+(cff2*cff-dmude(i,j+1,k2))*adfac
1482 ad_cff3=ad_cff3+(cff3*cff-dmude(i,j ,k2))*adfac
1483 ad_cff4=ad_cff4+(cff4*cff-dmude(i,j+1,k1))*adfac
1484 ad_dmudz=ad_dmudz+ &
1485 & (cff1*cff1+cff2*cff2+cff3*cff3+cff4*cff4)* &
1486 & adfac
1487 ad_dmude(i,j ,k1)=ad_dmude(i,j ,k1)-cff1*adfac
1488 ad_dmude(i,j+1,k2)=ad_dmude(i,j+1,k2)-cff2*adfac
1489 ad_dmude(i,j ,k2)=ad_dmude(i,j ,k2)-cff3*adfac
1490 ad_dmude(i,j+1,k1)=ad_dmude(i,j+1,k1)-cff4*adfac
1491 ad_ufse(i,j,k2)=0.0_r8
1492!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZde_p(i,j+1,k1)))* &
1493!^ & tl_dZde_p(i,j+1,k1)
1494!^
1495 ad_dzde_p(i,j+1,k1)=ad_dzde_p(i,j+1,k1)+ &
1496 & (0.5_r8+ &
1497 & sign(0.5_r8, dzde_p(i,j+1,k1)))* &
1498 & ad_cff4
1499 ad_cff4=0.0_r8
1500!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZde_p(i,j ,k2)))* &
1501!^ & tl_dZde_p(i,j ,k2)
1502!^
1503 ad_dzde_p(i,j ,k2)=ad_dzde_p(i,j ,k2)+ &
1504 & (0.5_r8+ &
1505 & sign(0.5_r8, dzde_p(i,j ,k2)))* &
1506 & ad_cff3
1507 ad_cff3=0.0_r8
1508!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZde_p(i,j+1,k2)))* &
1509!^ & tl_dZde_p(i,j+1,k2)
1510!^
1511 ad_dzde_p(i,j+1,k2)=ad_dzde_p(i,j+1,k2)+ &
1512 & (0.5_r8+ &
1513 & sign(0.5_r8,-dzde_p(i,j+1,k2)))* &
1514 & ad_cff2
1515 ad_cff2=0.0_r8
1516!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZde_p(i,j ,k1)))* &
1517!^ & tl_dZde_p(i,j ,k1)
1518!^
1519 ad_dzde_p(i,j ,k1)=ad_dzde_p(i,j ,k1)+ &
1520 & (0.5_r8+ &
1521 & sign(0.5_r8,-dzde_p(i,j ,k1)))* &
1522 & ad_cff1
1523 ad_cff1=0.0_r8
1524!
1525 cff1=min(dzdx_r(i-1,j,k1),0.0_r8)
1526 cff2=min(dzdx_r(i ,j,k2),0.0_r8)
1527 cff3=max(dzdx_r(i-1,j,k2),0.0_r8)
1528 cff4=max(dzdx_r(i ,j,k1),0.0_r8)
1529#ifdef VISC_3DCOEF
1530!^ tl_UFsx(i,j,k2)=tl_UFsx(i,j,k2)+ &
1531!^ & tl_fac1* &
1532!^ & (cff1*(cff1*dnUdz-dnUdx(i-1,j,k1))+ &
1533!^ & cff2*(cff2*dnUdz-dnUdx(i ,j,k2))+ &
1534!^ & cff3*(cff3*dnUdz-dnUdx(i-1,j,k2))+ &
1535!^ & cff4*(cff4*dnUdz-dnUdx(i ,j,k1)))
1536!^
1537 ad_fac1=ad_fac1+ &
1538 & (cff1*(cff1*dnudz-dnudx(i-1,j,k1))+ &
1539 & cff2*(cff2*dnudz-dnudx(i ,j,k2))+ &
1540 & cff3*(cff3*dnudz-dnudx(i-1,j,k2))+ &
1541 & cff4*(cff4*dnudz-dnudx(i ,j,k1)))* &
1542 & ad_ufsx(i,j,k2)
1543#endif
1544!^ tl_UFsx(i,j,k2)=fac1* &
1545!^ & (tl_cff1*(cff1*dnUdz-dnUdx(i-1,j,k1))+ &
1546!^ & tl_cff2*(cff2*dnUdz-dnUdx(i ,j,k2))+ &
1547!^ & tl_cff3*(cff3*dnUdz-dnUdx(i-1,j,k2))+ &
1548!^ & tl_cff4*(cff4*dnUdz-dnUdx(i ,j,k1))+ &
1549!^ & cff1*(tl_cff1*dnUdz+cff1*tl_dnUdz- &
1550!^ & tl_dnUdx(i-1,j,k1))+ &
1551!^ & cff2*(tl_cff2*dnUdz+cff2*tl_dnUdz- &
1552!^ & tl_dnUdx(i ,j,k2))+ &
1553!^ & cff3*(tl_cff3*dnUdz+cff3*tl_dnUdz- &
1554!^ & tl_dnUdx(i-1,j,k2))+ &
1555!^ & cff4*(tl_cff4*dnUdz+cff4*tl_dnUdz- &
1556!^ & tl_dnUdx(i ,j,k1)))
1557!^
1558 cff=2.0_r8*dnudz
1559 adfac=fac1*ad_ufsx(i,j,k2)
1560 ad_cff1=ad_cff1+(cff1*cff-dnudx(i-1,j,k1))*adfac
1561 ad_cff2=ad_cff2+(cff2*cff-dnudx(i ,j,k2))*adfac
1562 ad_cff3=ad_cff3+(cff3*cff-dnudx(i-1,j,k2))*adfac
1563 ad_cff4=ad_cff4+(cff4*cff-dnudx(i ,j,k1))*adfac
1564 ad_dnudz=ad_dnudz+ &
1565 & (cff1*cff1+cff2*cff2+cff3*cff3+cff4*cff4)* &
1566 & adfac
1567 ad_dnudx(i-1,j,k1)=ad_dnudx(i-1,j,k1)-cff1*adfac
1568 ad_dnudx(i ,j,k2)=ad_dnudx(i ,j,k2)-cff2*adfac
1569 ad_dnudx(i-1,j,k2)=ad_dnudx(i-1,j,k2)-cff3*adfac
1570 ad_dnudx(i ,j,k1)=ad_dnudx(i ,j,k1)-cff4*adfac
1571 ad_ufsx(i,j,k2)=0.0_r8
1572!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZdx_r(i ,j,k1)))* &
1573!^ & tl_dZdx_r(i ,j,k1)
1574!^
1575 ad_dzdx_r(i ,j,k1)=ad_dzdx_r(i ,j,k1)+ &
1576 & (0.5_r8+ &
1577 & sign(0.5_r8, dzdx_r(i ,j,k1)))* &
1578 & ad_cff4
1579 ad_cff4=0.0_r8
1580!^ tl_cff3=(0.5_r8+SIGN(0.5_r8, dZdx_r(i-1,j,k2)))* &
1581!^ & tl_dZdx_r(i-1,j,k2)
1582!^
1583 ad_dzdx_r(i-1,j,k2)=ad_dzdx_r(i-1,j,k2)+ &
1584 & (0.5_r8+ &
1585 & sign(0.5_r8, dzdx_r(i-1,j,k2)))* &
1586 & ad_cff3
1587 ad_cff3=0.0_r8
1588!^ tl_cff2=(0.5_r8+SIGN(0.5_r8,-dZdx_r(i ,j,k2)))* &
1589!^ & tl_dZdx_r(i ,j,k2)
1590!^
1591 ad_dzdx_r(i ,j,k2)=ad_dzdx_r(i ,j,k2)+ &
1592 & (0.5_r8+ &
1593 & sign(0.5_r8,-dzdx_r(i ,j,k2)))* &
1594 & ad_cff2
1595 ad_cff2=0.0_r8
1596!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZdx_r(i-1,j,k1)))* &
1597!^ & tl_dZdx_r(i-1,j,k1)
1598!^
1599 ad_dzdx_r(i-1,j,k1)=ad_dzdx_r(i-1,j,k1)+ &
1600 & (0.5_r8+ &
1601 & sign(0.5_r8,-dzdx_r(i-1,j,k1)))* &
1602 & ad_cff1
1603 ad_cff1=0.0_r8
1604!
1605 cff=0.5_r8*(pm(i-1,j)+pm(i,j))
1606!^ tl_dmVdz=cff*0.25_r8*(tl_dVdz(i-1,j+1,k2)+ &
1607!^ & tl_dVdz(i ,j+1,k2)+ &
1608!^ & tl_dVdz(i-1,j ,k2)+ &
1609!^ & tl_dVdz(i ,j ,k2))
1610!^
1611 adfac=cff*0.25_r8*ad_dmvdz
1612 ad_dvdz(i-1,j ,k2)=ad_dvdz(i-1,j ,k2)+adfac
1613 ad_dvdz(i ,j ,k2)=ad_dvdz(i ,j ,k2)+adfac
1614 ad_dvdz(i-1,j+1,k2)=ad_dvdz(i-1,j+1,k2)+adfac
1615 ad_dvdz(i ,j+1,k2)=ad_dvdz(i ,j+1,k2)+adfac
1616 ad_dmvdz=0.0_r8
1617!^ tl_dmUdz=cff*tl_dUdz(i,j,k2)
1618!^
1619 ad_dudz(i,j,k2)=ad_dudz(i,j,k2)+cff*ad_dmudz
1620 ad_dmudz=0.0_r8
1621!
1622 cff=0.5_r8*(pn(i-1,j)+pn(i,j))
1623!^ tl_dnVdz=cff*0.25_r8*(tl_dVdz(i-1,j+1,k2)+ &
1624!^ & tl_dVdz(i ,j+1,k2)+ &
1625!^ & tl_dVdz(i-1,j ,k2)+ &
1626!^ & tl_dVdz(i ,j ,k2))
1627!^
1628 adfac=cff*0.25_r8*ad_dnvdz
1629 ad_dvdz(i-1,j ,k2)=ad_dvdz(i-1,j ,k2)+adfac
1630 ad_dvdz(i ,j ,k2)=ad_dvdz(i ,j ,k2)+adfac
1631 ad_dvdz(i-1,j+1,k2)=ad_dvdz(i-1,j+1,k2)+adfac
1632 ad_dvdz(i ,j+1,k2)=ad_dvdz(i ,j+1,k2)+adfac
1633 ad_dnvdz=0.0_r8
1634!^ tl_dnUdz=cff*tl_dUdz(i,j,k2)
1635!^
1636 ad_dudz(i,j,k2)=ad_dudz(i,j,k2)+cff*ad_dnudz
1637 ad_dnudz=0.0_r8
1638#ifdef VISC_3DCOEF
1639!^ tl_fac2=tl_cff*om_u(i,j)
1640!^ tl_fac1=tl_cff*on_u(i,j)
1641!^
1642 ad_cff=ad_cff+ &
1643 & on_u(i,j)*ad_fac1+om_u(i,j)*ad_fac2
1644 ad_fac1=0.0_r8
1645 ad_fac2=0.0_r8
1646!^ tl_cff=0.125_r8* &
1647!^ & (tl_visc3d_r(i-1,j,k )+tl_visc3d_r(i,j,k )+ &
1648!^ & tl_visc3d_r(i-1,j,k+1)+tl_visc3d_r(i,j,k+1))
1649!^
1650 adfac=0.125_r8*ad_cff
1651 ad_visc3d_r(i-1,j,k )=ad_visc3d_r(i-1,j,k )+adfac
1652 ad_visc3d_r(i ,j,k )=ad_visc3d_r(i ,j,k )+adfac
1653 ad_visc3d_r(i-1,j,k+1)=ad_visc3d_r(i-1,j,k+1)+adfac
1654 ad_visc3d_r(i ,j,k+1)=ad_visc3d_r(i ,j,k+1)+adfac
1655 ad_cff=0.0_r8
1656#endif
1657 END DO
1658 END DO
1659 END IF below_surface
1660!
1661 DO j=jstr,jend+1
1662 DO i=istr,iend+1
1663 pm_p=0.25_r8*(pm(i-1,j-1)+pm(i-1,j)+ &
1664 & pm(i ,j-1)+pm(i ,j))
1665 pn_p=0.25_r8*(pn(i-1,j-1)+pn(i-1,j)+ &
1666 & pn(i ,j-1)+pn(i ,j))
1667 cff1=min(dzdx_p(i,j,k1),0.0_r8)
1668 cff2=max(dzdx_p(i,j,k1),0.0_r8)
1669 cff3=min(dzde_p(i,j,k1),0.0_r8)
1670 cff4=max(dzde_p(i,j,k1),0.0_r8)
1671#ifdef VISC_3DCOEF
1672 cff=0.25_r8* &
1673 & (hz(i-1,j ,k)+hz(i,j ,k)+ &
1674 & hz(i-1,j-1,k)+hz(i,j-1,k))* &
1675 & (on_p(i,j)*(dnvdx(i,j,k1)- &
1676 & 0.5_r8*pn_p* &
1677 & (cff1*(dvdz(i-1,j,k1)+ &
1678 & dvdz(i ,j,k2))+ &
1679 & cff2*(dvdz(i-1,j,k2)+ &
1680 & dvdz(i ,j,k1))))+ &
1681 & om_p(i,j)*(dmude(i,j,k1)- &
1682 & 0.5_r8*pm_p* &
1683 & (cff3*(dudz(i,j-1,k1)+ &
1684 & dudz(i,j ,k2))+ &
1685 & cff4*(dudz(i,j-1,k2)+ &
1686 & dudz(i,j ,k1)))))
1687# ifdef MASKING
1688 cff=cff*pmask(i,j)
1689# endif
1690 visc_p=0.25_r8* &
1691 & (visc3d_r(i-1,j-1,k)+visc3d_r(i-1,j,k)+ &
1692 & visc3d_r(i ,j-1,k)+visc3d_r(i ,j,k))
1693!^ tl_VFx(i,j)=on_p(i,j)*on_p(i,j)* &
1694!^ & (tl_visc_p*cff+visc_p*tl_cff)
1695!^
1696 adfac=on_p(i,j)*on_p(i,j)*ad_vfx(i,j)
1697 ad_cff=ad_cff+visc_p*adfac
1698 ad_visc_p=ad_visc_p+cff*adfac
1699 ad_vfx(i,j)=0.0_r8
1700!^ tl_UFe(i,j)=om_p(i,j)*om_p(i,j)* &
1701!^ & (tl_visc_p*cff+visc_p*tl_cff)
1702!^
1703 adfac=om_p(i,j)*om_p(i,j)*ad_ufe(i,j)
1704 ad_cff=ad_cff+visc_p*adfac
1705 ad_visc_p=ad_visc_p+cff*adfac
1706 ad_ufe(i,j)=0.0_r8
1707!^ tl_visc_p=0.25_r8* &
1708!^ & (tl_visc3d_r(i-1,j-1,k)+tl_visc3d_r(i-1,j,k)+ &
1709!^ & tl_visc3d_r(i ,j-1,k)+tl_visc3d_r(i ,j,k))
1710!^
1711 adfac=0.25_r8*ad_visc_p
1712 ad_visc3d_r(i-1,j-1,k)=ad_visc3d_r(i-1,j-1,k)+adfac
1713 ad_visc3d_r(i ,j-1,k)=ad_visc3d_r(i ,j-1,k)+adfac
1714 ad_visc3d_r(i-1,j ,k)=ad_visc3d_r(i-1,j ,k)+adfac
1715 ad_visc3d_r(i ,j ,k)=ad_visc3d_r(i ,j ,k)+adfac
1716 ad_visc_p=0.0_r8
1717#else
1718!^ tl_VFx(i,j)=on_p(i,j)*on_p(i,j)*visc2_p(i,j)*tl_cff
1719!^ tl_UFe(i,j)=om_p(i,j)*om_p(i,j)*visc2_p(i,j)*tl_cff
1720!^
1721 ad_cff=ad_cff+ &
1722 & on_p(i,j)*on_p(i,j)*visc2_p(i,j)*ad_vfx(i,j)+ &
1723 & om_p(i,j)*om_p(i,j)*visc2_p(i,j)*ad_ufe(i,j)
1724 ad_vfx(i,j)=0.0_r8
1725 ad_ufe(i,j)=0.0_r8
1726#endif
1727#ifdef MASKING
1728!^ tl_cff=tl_cff*pmask(i,j)
1729!^
1730 ad_cff=ad_cff*pmask(i,j)
1731#endif
1732!^ tl_cff=0.25_r8* &
1733!^ & ((tl_Hz(i-1,j ,k)+tl_Hz(i,j ,k)+ &
1734!^ & tl_Hz(i-1,j-1,k)+tl_Hz(i,j-1,k))* &
1735!^ & (on_p(i,j)*(dnVdx(i,j,k1)- &
1736!^ & 0.5_r8*pn_p* &
1737!^ & (cff1*(dVdz(i-1,j,k1)+ &
1738!^ & dVdz(i ,j,k2))+ &
1739!^ & cff2*(dVdz(i-1,j,k2)+ &
1740!^ & dVdz(i ,j,k1))))+ &
1741!^ & om_p(i,j)*(dmUde(i,j,k1)- &
1742!^ & 0.5_r8*pm_p* &
1743!^ & (cff3*(dUdz(i,j-1,k1)+ &
1744!^ & dUdz(i,j ,k2))+ &
1745!^ & cff4*(dUdz(i,j-1,k2)+ &
1746!^ & dUdz(i,j ,k1)))))+ &
1747!^ & (Hz(i-1,j ,k)+Hz(i,j ,k)+ &
1748!^ & Hz(i-1,j-1,k)+Hz(i,j-1,k))* &
1749!^ & (on_p(i,j)*(tl_dnVdx(i,j,k1)- &
1750!^ & 0.5_r8*pn_p* &
1751!^ & (tl_cff1*(dVdz(i-1,j,k1)+ &
1752!^ & dVdz(i ,j,k2))+ &
1753!^ & cff1*(tl_dVdz(i-1,j,k1)+ &
1754!^ & tl_dVdz(i ,j,k2))+ &
1755!^ & tl_cff2*(dVdz(i-1,j,k2)+ &
1756!^ & dVdz(i ,j,k1))+ &
1757!^ & cff2*(tl_dVdz(i-1,j,k2)+ &
1758!^ & tl_dVdz(i ,j,k1)))+ &
1759!^ & om_p(i,j)*(tl_dmUde(i,j,k1)- &
1760!^ & 0.5_r8*pm_p* &
1761!^ & (tl_cff3*(dUdz(i,j-1,k1)+ &
1762!^ & dUdz(i,j ,k2))+ &
1763!^ & cff3*(tl_dUdz(i,j-1,k1)+ &
1764!^ & tl_dUdz(i,j ,k2))+ &
1765!^ & tl_cff4*(dUdz(i,j-1,k2)+ &
1766!^ & dUdz(i,j ,k1))+ &
1767!^ & cff4*(tl_dUdz(i,j-1,k2)+ &
1768!^ & tl_dUdz(i,j ,k1)))))))
1769!^
1770 adfac=0.25_r8*ad_cff
1771 ad_cff=0.0_r8
1772 adfac1=adfac*(on_p(i,j)*(dnvdx(i,j,k1)- &
1773 & 0.5_r8*pn_p* &
1774 & (cff1*(dvdz(i-1,j,k1)+ &
1775 & dvdz(i ,j,k2))+ &
1776 & cff2*(dvdz(i-1,j,k2)+ &
1777 & dvdz(i ,j,k1))))+ &
1778 & om_p(i,j)*(dmude(i,j,k1)- &
1779 & 0.5_r8*pm_p* &
1780 & (cff3*(dudz(i,j-1,k1)+ &
1781 & dudz(i,j ,k2))+ &
1782 & cff4*(dudz(i,j-1,k2)+ &
1783 & dudz(i,j ,k1)))))
1784 adfac2=adfac*(hz(i-1,j ,k)+hz(i,j ,k)+ &
1785 & hz(i-1,j-1,k)+hz(i,j-1,k))
1786 adfac3=adfac2*on_p(i,j)
1787 adfac4=adfac3*0.5_r8*pn_p
1788 adfac5=adfac2*om_p(i,j)
1789 adfac6=adfac5*0.5_r8*pm_p
1790 ad_hz(i-1,j-1,k)=ad_hz(i-1,j-1,k)+adfac1
1791 ad_hz(i ,j-1,k)=ad_hz(i ,j-1,k)+adfac1
1792 ad_hz(i-1,j ,k)=ad_hz(i-1,j ,k)+adfac1
1793 ad_hz(i ,j ,k)=ad_hz(i ,j ,k)+adfac1
1794 ad_dnvdx(i,j,k1)=ad_dnvdx(i,j,k1)+adfac3
1795 ad_cff1=ad_cff1- &
1796 & (dvdz(i-1,j,k1)+dvdz(i ,j,k2))*adfac4
1797 ad_cff2=ad_cff2- &
1798 & (dvdz(i-1,j,k2)+dvdz(i ,j,k1))*adfac4
1799 ad_dvdz(i-1,j,k1)=ad_dvdz(i-1,j,k1)-cff1*adfac4
1800 ad_dvdz(i-1,j,k2)=ad_dvdz(i-1,j,k2)-cff2*adfac4
1801 ad_dvdz(i ,j,k1)=ad_dvdz(i ,j,k1)-cff2*adfac4
1802 ad_dvdz(i ,j,k2)=ad_dvdz(i ,j,k2)-cff1*adfac4
1803 ad_dmude(i,j,k1)=ad_dmude(i,j,k1)+adfac5
1804 ad_cff3=ad_cff3- &
1805 & (dudz(i,j-1,k1)+dudz(i,j ,k2))*adfac6
1806 ad_cff4=ad_cff4- &
1807 & (dudz(i,j-1,k2)+dudz(i,j ,k1))*adfac6
1808 ad_dudz(i,j-1,k1)=ad_dudz(i,j-1,k1)-cff3*adfac6
1809 ad_dudz(i,j-1,k2)=ad_dudz(i,j-1,k2)-cff4*adfac6
1810 ad_dudz(i,j ,k1)=ad_dudz(i,j ,k1)-cff4*adfac6
1811 ad_dudz(i,j ,k2)=ad_dudz(i,j ,k2)-cff3*adfac6
1812!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZde_p(i,j,k1)))* &
1813!^ & tl_dZde_p(i,j,k1)
1814!^ tl_cff3=(0.5_r8+SIGN(0.5_r8,-dZde_p(i,j,k1)))* &
1815!^ & tl_dZde_p(i,j,k1)
1816!^
1817 ad_dzde_p(i,j,k1)=ad_dzde_p(i,j,k1)+ &
1818 & (0.5_r8+ &
1819 & sign(0.5_r8, dzde_p(i,j,k1)))* &
1820 & ad_cff4+ &
1821 & (0.5_r8+ &
1822 & sign(0.5_r8,-dzde_p(i,j,k1)))* &
1823 & ad_cff3
1824 ad_cff4=0.0_r8
1825 ad_cff3=0.0_r8
1826!^ tl_cff2=(0.5_r8+SIGN(0.5_r8, dZdx_p(i,j,k1)))* &
1827!^ & tl_dZdx_p(i,j,k1)
1828!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZdx_p(i,j,k1)))* &
1829!^ & tl_dZdx_p(i,j,k1)
1830!^
1831 ad_dzdx_p(i,j,k1)=ad_dzdx_p(i,j,k1)+ &
1832 & (0.5_r8+ &
1833 & sign(0.5_r8, dzdx_p(i,j,k1)))* &
1834 & ad_cff2+ &
1835 & (0.5_r8+ &
1836 & sign(0.5_r8,-dzdx_p(i,j,k1)))* &
1837 & ad_cff1
1838 ad_cff2=0.0_r8
1839 ad_cff1=0.0_r8
1840 END DO
1841 END DO
1842!
1843 DO j=jstrv-1,jend
1844 DO i=istru-1,iend
1845 cff1=min(dzdx_r(i,j,k1),0.0_r8)
1846 cff2=max(dzdx_r(i,j,k1),0.0_r8)
1847 cff3=min(dzde_r(i,j,k1),0.0_r8)
1848 cff4=max(dzde_r(i,j,k1),0.0_r8)
1849#ifdef VISC_3DCOEF
1850 cff=hz(i,j,k)* &
1851 & (on_r(i,j)*(dnudx(i,j,k1)- &
1852 & 0.5_r8*pn(i,j)* &
1853 & (cff1*(dudz(i ,j,k1)+ &
1854 & dudz(i+1,j,k2))+ &
1855 & cff2*(dudz(i ,j,k2)+ &
1856 & dudz(i+1,j,k1))))- &
1857 & om_r(i,j)*(dmvde(i,j,k1)- &
1858 & 0.5_r8*pm(i,j)* &
1859 & (cff3*(dvdz(i,j ,k1)+ &
1860 & dvdz(i,j+1,k2))+ &
1861 & cff4*(dvdz(i,j ,k2)+ &
1862 & dvdz(i,j+1,k1)))))
1863# ifdef MASKING
1864 cff=cff*rmask(i,j)
1865# endif
1866!^ tl_VFe(i,j)=om_r(i,j)*om_r(i,j)* &
1867!^ & (tl_visc3d_r(i,j,k)*cff+ &
1868!^ & visc3d_r(i,j,k)*tl_cff)
1869!^
1870 adfac=om_r(i,j)*om_r(i,j)*ad_vfe(i,j)
1871 ad_cff=ad_cff+visc3d_r(i,j,k)*adfac
1872 ad_visc3d_r(i,j,k)=ad_visc3d_r(i,j,k)+cff*adfac
1873 ad_vfe(i,j)=0.0_r8
1874!^ tl_UFx(i,j)=on_r(i,j)*on_r(i,j)* &
1875!^ & (tl_visc3d_r(i,j,k)*cff+ &
1876!^ & visc3d_r(i,j,k)*tl_cff)
1877!^
1878 adfac=on_r(i,j)*on_r(i,j)*ad_ufx(i,j)
1879 ad_cff=ad_cff+visc3d_r(i,j,k)*adfac
1880 ad_visc3d_r(i,j,k)=ad_visc3d_r(i,j,k)+cff*adfac
1881 ad_ufx(i,j)=0.0_r8
1882#else
1883!^ tl_VFe(i,j)=om_r(i,j)*om_r(i,j)*visc2_r(i,j)*tl_cff
1884!^ tl_UFx(i,j)=on_r(i,j)*on_r(i,j)*visc2_r(i,j)*tl_cff
1885!^
1886 ad_cff=ad_cff+ &
1887 & om_r(i,j)*om_r(i,j)*visc2_r(i,j)*ad_vfe(i,j)+ &
1888 & on_r(i,j)*on_r(i,j)*visc2_r(i,j)*ad_ufx(i,j)
1889 ad_vfe(i,j)=0.0_r8
1890 ad_ufx(i,j)=0.0_r8
1891#endif
1892#ifdef MASKING
1893!^ tl_cff=tl_cff*rmask(i,j)
1894!^
1895 ad_cff=ad_cff*rmask(i,j)
1896#endif
1897!^ tl_cff=tl_Hz(i,j,k)* &
1898!^ & (on_r(i,j)*(dnUdx(i,j,k1)- &
1899!^ & 0.5_r8*pn(i,j)* &
1900!^ & (cff1*(dUdz(i ,j,k1)+ &
1901!^ & dUdz(i+1,j,k2))+ &
1902!^ & cff2*(dUdz(i ,j,k2)+ &
1903!^ & dUdz(i+1,j,k1))))- &
1904!^ & om_r(i,j)*(dmVde(i,j,k1)- &
1905!^ & 0.5_r8*pm(i,j)* &
1906!^ & (cff3*(dVdz(i,j ,k1)+ &
1907!^ & dVdz(i,j+1,k2))+ &
1908!^ & cff4*(dVdz(i,j ,k2)+ &
1909!^ & dVdz(i,j+1,k1)))))+ &
1910!^ & Hz(i,j,k)* &
1911!^ & (on_r(i,j)*(tl_dnUdx(i,j,k1)- &
1912!^ & 0.5_r8*pn(i,j)* &
1913!^ & (tl_cff1*(dUdz(i ,j,k1)+ &
1914!^ & dUdz(i+1,j,k2))+ &
1915!^ & cff1*(tl_dUdz(i ,j,k1)+ &
1916!^ & tl_dUdz(i+1,j,k2))+ &
1917!^ & tl_cff2*(dUdz(i ,j,k2)+ &
1918!^ & dUdz(i+1,j,k1))+ &
1919!^ & cff2*(tl_dUdz(i ,j,k2)+ &
1920!^ & tl_dUdz(i+1,j,k1))))- &
1921!^ & om_r(i,j)*(tl_dmVde(i,j,k1)- &
1922!^ & 0.5_r8*pm(i,j)* &
1923!^ & (tl_cff3*(dVdz(i,j ,k1)+ &
1924!^ & dVdz(i,j+1,k2))+ &
1925!^ & cff3*(tl_dVdz(i,j ,k1)+ &
1926!^ & tl_dVdz(i,j+1,k2))+ &
1927!^ & tl_cff4*(dVdz(i,j ,k2)+ &
1928!^ & dVdz(i,j+1,k1))+ &
1929!^ & cff4*(tl_dVdz(i,j ,k2)+ &
1930!^ & tl_dVdz(i,j+1,k1)))))
1931!^
1932 adfac1=hz(i,j,k)*ad_cff
1933 adfac2=adfac1*on_r(i,j)
1934 adfac3=adfac2*0.5_r8*pn(i,j)
1935 adfac4=adfac1*om_r(i,j)
1936 adfac5=adfac4*0.5_r8*pm(i,j)
1937 ad_hz(i,j,k)=ad_hz(i,j,k)+ &
1938 (on_r(i,j)*(dnudx(i,j,k1)- &
1939 & 0.5_r8*pn(i,j)* &
1940 & (cff1*(dudz(i ,j,k1)+ &
1941 & dudz(i+1,j,k2))+ &
1942 & cff2*(dudz(i ,j,k2)+ &
1943 & dudz(i+1,j,k1))))- &
1944 & om_r(i,j)*(dmvde(i,j,k1)- &
1945 & 0.5_r8*pm(i,j)* &
1946 & (cff3*(dvdz(i,j ,k1)+ &
1947 & dvdz(i,j+1,k2))+ &
1948 & cff4*(dvdz(i,j ,k2)+ &
1949 & dvdz(i,j+1,k1)))))* &
1950 & ad_cff
1951 ad_dnudx(i,j,k1)=ad_dnudx(i,j,k1)+adfac2
1952 ad_cff1=ad_cff1- &
1953 & (dudz(i ,j,k1)+dudz(i+1,j,k2))*adfac3
1954 ad_cff2=ad_cff2- &
1955 (dudz(i ,j,k2)+dudz(i+1,j,k1))*adfac3
1956 ad_dudz(i ,j,k1)=ad_dudz(i ,j,k1)-cff1*adfac3
1957 ad_dudz(i ,j,k2)=ad_dudz(i ,j,k2)-cff2*adfac3
1958 ad_dudz(i+1,j,k1)=ad_dudz(i+1,j,k1)-cff2*adfac3
1959 ad_dudz(i+1,j,k2)=ad_dudz(i+1,j,k2)-cff1*adfac3
1960 ad_dmvde(i,j,k1)=ad_dmvde(i,j,k1)-adfac4
1961 ad_cff3=ad_cff3+ &
1962 & (dvdz(i,j ,k1)+dvdz(i,j+1,k2))*adfac5
1963 ad_cff4=ad_cff4+ &
1964 & (dvdz(i,j ,k2)+dvdz(i,j+1,k1))*adfac5
1965 ad_dvdz(i,j ,k1)=ad_dvdz(i,j ,k1)+cff3*adfac5
1966 ad_dvdz(i,j ,k2)=ad_dvdz(i,j ,k2)+cff4*adfac5
1967 ad_dvdz(i,j+1,k1)=ad_dvdz(i,j+1,k1)+cff4*adfac5
1968 ad_dvdz(i,j+1,k2)=ad_dvdz(i,j+1,k2)+cff3*adfac5
1969 ad_cff=0.0_r8
1970!^ tl_cff4=(0.5_r8+SIGN(0.5_r8, dZde_r(i,j,k1)))* &
1971!^ & tl_dZde_r(i,j,k1)
1972!^ tl_cff3=(0.5_r8+SIGN(0.5_r8,-dZde_r(i,j,k1)))* &
1973!^ & tl_dZde_r(i,j,k1)
1974!^
1975 ad_dzde_r(i,j,k1)=ad_dzde_r(i,j,k1)+ &
1976 & (0.5_r8+ &
1977 & sign(0.5_r8, dzde_r(i,j,k1)))* &
1978 & ad_cff4+ &
1979 & (0.5_r8+ &
1980 & sign(0.5_r8,-dzde_r(i,j,k1)))* &
1981 & ad_cff3
1982 ad_cff4=0.0_r8
1983 ad_cff3=0.0_r8
1984!^ tl_cff2=(0.5_r8+SIGN(0.5_r8, dZdx_r(i,j,k1)))* &
1985!^ & tl_dZdx_r(i,j,k1)
1986!^ tl_cff1=(0.5_r8+SIGN(0.5_r8,-dZdx_r(i,j,k1)))* &
1987!^ & tl_dZdx_r(i,j,k1)
1988!^
1989 ad_dzdx_r(i,j,k1)=ad_dzdx_r(i,j,k1)+ &
1990 & (0.5_r8+ &
1991 & sign(0.5_r8, dzdx_r(i,j,k1)))* &
1992 & ad_cff2+ &
1993 & (0.5_r8+ &
1994 & sign(0.5_r8,-dzdx_r(i,j,k1)))* &
1995 & ad_cff1
1996 ad_cff2=0.0_r8
1997 ad_cff1=0.0_r8
1998 END DO
1999 END DO
2000 END IF above_bottom
2001
2002 IF ((k.eq.0).or.(k.eq.n(ng))) THEN
2003 DO j=jstrv-1,jend+1
2004 DO i=istr-1,iend+1
2005!^ tl_VFse(i,j,k2)=0.0_r8
2006!^
2007 ad_vfse(i,j,k2)=0.0_r8
2008!^ tl_VFsx(i,j,k2)=0.0_r8
2009!^
2010 ad_vfsx(i,j,k2)=0.0_r8
2011 END DO
2012 END DO
2013 DO j=jstr-1,jend+1
2014 DO i=istru-1,iend+1
2015!^ tl_UFse(i,j,k2)=0.0_r8
2016!^
2017 ad_ufse(i,j,k2)=0.0_r8
2018!^ tl_UFsx(i,j,k2)=0.0_r8
2019!^
2020 ad_ufsx(i,j,k2)=0.0_r8
2021 END DO
2022 END DO
2023
2024 DO j=jstrv-1,jend+1
2025 DO i=istr-1,iend+1
2026!^ tl_dVdz(i,j,k2)=0.0_r8
2027!^
2028 ad_dvdz(i,j,k2)=0.0_r8
2029 END DO
2030 END DO
2031 DO j=jstr-1,jend+1
2032 DO i=istru-1,iend+1
2033
2034!^ tl_dUdz(i,j,k2)=0.0_r8
2035!^
2036 ad_dudz(i,j,k2)=0.0_r8
2037 END DO
2038 END DO
2039 ELSE
2040 DO j=jstrv-1,jend+1
2041 DO i=istr-1,iend+1
2042 cff=1.0_r8/(0.5_r8*(z_r(i,j-1,k+1)-z_r(i,j-1,k)+ &
2043 & z_r(i,j ,k+1)-z_r(i,j ,k)))
2044!^ tl_dVdz(i,j,k2)=tl_cff*(v(i,j,k+1,nrhs)- &
2045!^ & v(i,j,k ,nrhs))+ &
2046!^ & cff*(tl_v(i,j,k+1,nrhs)- &
2047!^ & tl_v(i,j,k ,nrhs))
2048!^
2049 adfac=cff*ad_dvdz(i,j,k2)
2050 ad_v(i,j,k ,nrhs)=ad_v(i,j,k ,nrhs)-adfac
2051 ad_v(i,j,k+1,nrhs)=ad_v(i,j,k+1,nrhs)+adfac
2052 ad_cff=ad_cff+(v(i,j,k+1,nrhs)- &
2053 & v(i,j,k ,nrhs))*ad_dvdz(i,j,k2)
2054 ad_dvdz(i,j,k2)=0.0_r8
2055!^ tl_cff=-cff*cff*(0.5_r8*(tl_z_r(i,j-1,k+1)- &
2056!^ & tl_z_r(i,j-1,k )+ &
2057!^ & tl_z_r(i,j ,k+1)- &
2058!^ & tl_z_r(i,j ,k )))
2059!^
2060 adfac=-cff*cff*0.5_r8*ad_cff
2061 ad_z_r(i,j-1,k )=ad_z_r(i,j-1,k )-adfac
2062 ad_z_r(i,j-1,k+1)=ad_z_r(i,j-1,k+1)+adfac
2063 ad_z_r(i,j ,k )=ad_z_r(i,j ,k )-adfac
2064 ad_z_r(i,j ,k+1)=ad_z_r(i,j ,k+1)+adfac
2065 ad_cff=0.0_r8
2066 END DO
2067 END DO
2068
2069 DO j=jstr-1,jend+1
2070 DO i=istru-1,iend+1
2071 cff=1.0_r8/(0.5_r8*(z_r(i-1,j,k+1)-z_r(i-1,j,k)+ &
2072 & z_r(i ,j,k+1)-z_r(i ,j,k)))
2073!^ tl_dUdz(i,j,k2)=tl_cff*(u(i,j,k+1,nrhs)- &
2074!^ & u(i,j,k ,nrhs))+ &
2075!^ & cff*(tl_u(i,j,k+1,nrhs)- &
2076!^ & tl_u(i,j,k ,nrhs))
2077!^
2078 adfac=cff*ad_dudz(i,j,k2)
2079 ad_u(i,j,k ,nrhs)=ad_u(i,j,k ,nrhs)-adfac
2080 ad_u(i,j,k+1,nrhs)=ad_u(i,j,k+1,nrhs)+adfac
2081 ad_cff=ad_cff+(u(i,j,k+1,nrhs)- &
2082 & u(i,j,k ,nrhs))*ad_dudz(i,j,k2)
2083 ad_dudz(i,j,k2)=0.0_r8
2084!^ tl_cff=-cff*cff*(0.5_r8*((tl_z_r(i-1,j,k+1)- &
2085!^ & tl_z_r(i-1,j,k )+ &
2086!^ & tl_z_r(i ,j,k+1)- &
2087!^ & tl_z_r(i ,j,k)))
2088!^
2089 adfac=-cff*cff*0.5_r8*ad_cff
2090 ad_z_r(i-1,j,k )=ad_z_r(i-1,j,k )-adfac
2091 ad_z_r(i-1,j,k+1)=ad_z_r(i-1,j,k+1)+adfac
2092 ad_z_r(i ,j,k )=ad_z_r(i ,j,k )-adfac
2093 ad_z_r(i ,j,k+1)=ad_z_r(i ,j,k+1)+adfac
2094 ad_cff=0.0_r8
2095 END DO
2096 END DO
2097 END IF
2098
2099 IF (k.lt.n(ng)) THEN
2100 DO j=jstrv-1,jend
2101 DO i=istru-1,iend
2102 cff=0.5_r8*pn(i,j)
2103#ifdef MASKING
2104 cff=cff*rmask(i,j)
2105#endif
2106!^ tl_dmVde(i,j,k2)=cff*((pm(i,j )+pm(i,j+1))* &
2107!^ & tl_v(i,j+1,k+1,nrhs)- &
2108!^ & (pm(i,j-1)+pm(i,j ))* &
2109!^ & tl_v(i,j ,k+1,nrhs))
2110!^
2111 adfac=cff*ad_dmvde(i,j,k2)
2112 ad_v(i,j ,k+1,nrhs)=ad_v(i,j ,k+1,nrhs)- &
2113 & (pm(i,j-1)+pm(i,j ))*adfac
2114 ad_v(i,j+1,k+1,nrhs)=ad_v(i,j+1,k+1,nrhs)+ &
2115 & (pm(i,j )+pm(i,j+1))*adfac
2116 ad_dmvde(i,j,k2)=0.0_r8
2117 END DO
2118 END DO
2119
2120 DO j=jstr,jend+1
2121 DO i=istru-1,iend+1
2122 cff=0.125_r8*(pm(i-1,j )+pm(i,j )+ &
2123 & pm(i-1,j-1)+pm(i,j-1))
2124#ifdef MASKING
2125 cff=cff*pmask(i,j)
2126#endif
2127!^ tl_dnVdx(i,j,k2)=cff*((pn(i ,j-1)+pn(i ,j))* &
2128!^ & tl_v(i ,j,k+1,nrhs)- &
2129!^ & (pn(i-1,j-1)+pn(i-1,j))* &
2130!^ & tl_v(i-1,j,k+1,nrhs))
2131!^
2132 adfac=cff*ad_dnvdx(i,j,k2)
2133 ad_v(i-1,j,k+1,nrhs)=ad_v(i-1,j,k+1,nrhs)- &
2134 & (pn(i-1,j-1)+pn(i-1,j))*adfac
2135 ad_v(i ,j,k+1,nrhs)=ad_v(i ,j,k+1,nrhs)+ &
2136 & (pn(i ,j-1)+pn(i ,j))*adfac
2137 ad_dnvdx(i,j,k2)=0.0_r8
2138 END DO
2139 END DO
2140
2141 DO j=jstr,jend+1
2142 DO i=istr,iend+1
2143 cff=0.125_r8*(pn(i-1,j )+pn(i,j )+ &
2144 & pn(i-1,j-1)+pn(i,j-1))
2145#ifdef MASKING
2146 cff=cff*pmask(i,j)
2147#endif
2148!^ tl_dmUde(i,j,k2)=cff*((pm(i-1,j )+pm(i,j ))* &
2149!^ & tl_u(i,j ,k+1,nrhs)- &
2150!^ & (pm(i-1,j-1)+pm(i,j-1))* &
2151!^ & tl_u(i,j-1,k+1,nrhs))
2152!^
2153 adfac=cff*ad_dmude(i,j,k2)
2154 ad_u(i,j-1,k+1,nrhs)=ad_u(i,j-1,k+1,nrhs)- &
2155 & (pm(i-1,j-1)+pm(i,j-1))*adfac
2156 ad_u(i,j ,k+1,nrhs)=ad_u(i,j ,k+1,nrhs)+ &
2157 & (pm(i-1,j )+pm(i,j ))*adfac
2158 ad_dmude(i,j,k2)=0.0_r8
2159 END DO
2160 END DO
2161
2162 DO j=jstrv-1,jend
2163 DO i=istru-1,iend
2164 cff=0.5_r8*pm(i,j)
2165#ifdef MASKING
2166 cff=cff*rmask(i,j)
2167#endif
2168!^ tl_dnUdx(i,j,k2)=cff*((pn(i ,j)+pn(i+1,j))* &
2169!^ & tl_u(i+1,j,k+1,nrhs)- &
2170!^ & (pn(i-1,j)+pn(i ,j))* &
2171!^ & tl_u(i ,j,k+1,nrhs))
2172!^
2173 adfac=cff*ad_dnudx(i,j,k2)
2174 ad_u(i ,j,k+1,nrhs)=ad_u(i ,j,k+1,nrhs)- &
2175 & (pn(i-1,j)+pn(i ,j))*adfac
2176 ad_u(i+1,j,k+1,nrhs)=ad_u(i+1,j,k+1,nrhs)+ &
2177 & (pn(i ,j)+pn(i+1,j))*adfac
2178 ad_dnudx(i,j,k2)=0.0_r8
2179 END DO
2180 END DO
2181!
2182! Compute slopes (nondimensional) at RHO- and PSI-points.
2183!
2184 DO j=jstrv-1,jend
2185 DO i=istru-1,iend
2186!^ tl_dZde_r(i,j,k2)=0.5_r8*(tl_VFe(i,j )+ &
2187!^ & tl_VFe(i,j+1))
2188!^
2189 adfac=0.5_r8*ad_dzde_r(i,j,k2)
2190 ad_vfe(i,j )=ad_vfe(i,j )+adfac
2191 ad_vfe(i,j+1)=ad_vfe(i,j+1)+adfac
2192 ad_dzde_r(i,j,k2)=0.0_r8
2193!^ tl_dZdx_r(i,j,k2)=0.5_r8*(tl_UFx(i ,j)+ &
2194!^ & tl_UFx(i+1,j))
2195!^
2196 adfac=0.5_r8*ad_dzdx_r(i,j,k2)
2197 ad_ufx(i ,j)=ad_ufx(i ,j)+adfac
2198 ad_ufx(i+1,j)=ad_ufx(i+1,j)+adfac
2199 ad_dzdx_r(i,j,k2)=0.0_r8
2200 END DO
2201 END DO
2202
2203 DO j=jstr,jend+1
2204 DO i=istr,iend+1
2205!^ tl_dZde_p(i,j,k2)=0.5_r8*(tl_VFe(i-1,j)+ &
2206!^ & tl_VFe(i ,j))
2207!^
2208 adfac=0.5_r8*ad_dzde_p(i,j,k2)
2209 ad_vfe(i-1,j)=ad_vfe(i-1,j)+adfac
2210 ad_vfe(i ,j)=ad_vfe(i ,j)+adfac
2211 ad_dzde_p(i,j,k2)=0.0_r8
2212!^ tl_dZdx_p(i,j,k2)=0.5_r8*(tl_UFx(i,j-1)+ &
2213!^ & tl_UFx(i,j ))
2214!^
2215 adfac=0.5_r8*ad_dzdx_p(i,j,k2)
2216 ad_ufx(i,j-1)=ad_ufx(i,j-1)+adfac
2217 ad_ufx(i,j )=ad_ufx(i,j )+adfac
2218 ad_dzdx_p(i,j,k2)=0.0_r8
2219 END DO
2220 END DO
2221!
2222 DO j=jstrv-1,jend+1
2223 DO i=istr-1,iend+1
2224 cff=0.5_r8*(pn(i,j-1)+pn(i,j))
2225#ifdef MASKING
2226 cff=cff*vmask(i,j)
2227#endif
2228!^ tl_VFe(i,j)=cff*(tl_z_r(i,j ,k+1)- &
2229!^ & tl_z_r(i,j-1,k+1))
2230!^
2231 adfac=cff*ad_vfe(i,j)
2232 ad_z_r(i,j-1,k+1)=ad_z_r(i,j-1,k+1)-adfac
2233 ad_z_r(i,j ,k+1)=ad_z_r(i,j ,k+1)+adfac
2234 ad_vfe(i,j)=0.0_r8
2235 END DO
2236 END DO
2237
2238 DO j=jstr-1,jend+1
2239 DO i=istru-1,iend+1
2240 cff=0.5_r8*(pm(i-1,j)+pm(i,j))
2241#ifdef MASKING
2242 cff=cff*umask(i,j)
2243#endif
2244!^ tl_UFx(i,j)=cff*(tl_z_r(i ,j,k+1)- &
2245!^ & tl_z_r(i-1,j,k+1))
2246!^
2247 adfac=cff*ad_ufx(i,j)
2248 ad_z_r(i-1,j,k+1)=ad_z_r(i-1,j,k+1)-adfac
2249 ad_z_r(i ,j,k+1)=ad_z_r(i ,j,k+1)+adfac
2250 ad_ufx(i,j)=0.0_r8
2251 END DO
2252 END DO
2253 END IF
2254!
2255! Compute new recursive storage indices.
2256!
2257 kt=k2
2258 k2=k1
2259 k1=kt
2260 END DO k_loop
2261!
2262 RETURN
integer, dimension(:), allocatable n
Definition mod_param.F:479
real(dp), dimension(:), allocatable dt

References mod_scalars::dt.

Referenced by ad_uv3dmix2().

Here is the caller graph for this function:

◆ ad_uv3dmix2_s_tile()

subroutine ad_uv3dmix2_mod::ad_uv3dmix2_s_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) nnew,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmask,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(in) hz,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng)), intent(inout) ad_hz,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) om_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) on_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pm,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pmon_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pn,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) pnom_r,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_p,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(in) visc2_r,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(in) u,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(in) v,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) ad_rufrc,
real(r8), dimension(lbi:ubi,lbj:ubj), intent(inout) ad_rvfrc,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(inout) ad_u,
real(r8), dimension(lbi:ubi,lbj:ubj,n(ng),2), intent(inout) ad_v )
private

Definition at line 113 of file ad_uv3dmix2_s.h.

132!***********************************************************************
133!
134 USE mod_param
135 USE mod_scalars
136!
137! Imported variable declarations.
138!
139 integer, intent(in) :: ng, tile
140 integer, intent(in) :: LBi, UBi, LBj, UBj
141 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
142 integer, intent(in) :: nrhs, nnew
143
144#ifdef ASSUMED_SHAPE
145#ifdef MASKING
146 real(r8), intent(in) :: pmask(LBi:,LBj:)
147#endif
148 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
149 real(r8), intent(in) :: om_p(LBi:,LBj:)
150 real(r8), intent(in) :: om_r(LBi:,LBj:)
151 real(r8), intent(in) :: on_p(LBi:,LBj:)
152 real(r8), intent(in) :: on_r(LBi:,LBj:)
153 real(r8), intent(in) :: pm(LBi:,LBj:)
154 real(r8), intent(in) :: pmon_p(LBi:,LBj:)
155 real(r8), intent(in) :: pmon_r(LBi:,LBj:)
156 real(r8), intent(in) :: pn(LBi:,LBj:)
157 real(r8), intent(in) :: pnom_p(LBi:,LBj:)
158 real(r8), intent(in) :: pnom_r(LBi:,LBj:)
159 real(r8), intent(in) :: visc2_p(LBi:,LBj:)
160 real(r8), intent(in) :: visc2_r(LBi:,LBj:)
161
162 real(r8), intent(in) :: u(LBi:,LBj:,:,:)
163 real(r8), intent(in) :: v(LBi:,LBj:,:,:)
164
165 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
166 real(r8), intent(inout) :: ad_rufrc(LBi:,LBj:)
167 real(r8), intent(inout) :: ad_rvfrc(LBi:,LBj:)
168 real(r8), intent(inout) :: ad_u(LBi:,LBj:,:,:)
169 real(r8), intent(inout) :: ad_v(LBi:,LBj:,:,:)
170#else
171#ifdef MASKING
172 real(r8), intent(in) :: pmask(LBi:UBi,LBj:UBj)
173#endif
174 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
175 real(r8), intent(in) :: om_p(LBi:UBi,LBj:UBj)
176 real(r8), intent(in) :: om_r(LBi:UBi,LBj:UBj)
177 real(r8), intent(in) :: on_p(LBi:UBi,LBj:UBj)
178 real(r8), intent(in) :: on_r(LBi:UBi,LBj:UBj)
179 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
180 real(r8), intent(in) :: pmon_p(LBi:UBi,LBj:UBj)
181 real(r8), intent(in) :: pmon_r(LBi:UBi,LBj:UBj)
182 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
183 real(r8), intent(in) :: pnom_p(LBi:UBi,LBj:UBj)
184 real(r8), intent(in) :: pnom_r(LBi:UBi,LBj:UBj)
185 real(r8), intent(in) :: visc2_p(LBi:UBi,LBj:UBj)
186 real(r8), intent(in) :: visc2_r(LBi:UBi,LBj:UBj)
187
188 real(r8), intent(in) :: u(LBi:UBi,LBj:UBj,N(ng),2)
189 real(r8), intent(in) :: v(LBi:UBi,LBj:UBj,N(ng),2)
190
191 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
192 real(r8), intent(inout) :: ad_rufrc(LBi:UBi,LBj:UBj)
193 real(r8), intent(inout) :: ad_rvfrc(LBi:UBi,LBj:UBj)
194 real(r8), intent(inout) :: ad_u(LBi:UBi,LBj:UBj,N(ng),2)
195 real(r8), intent(inout) :: ad_v(LBi:UBi,LBj:UBj,N(ng),2)
196#endif
197!
198! Local variable declarations.
199!
200 integer :: i, j, k
201
202 real(r8) :: cff, ad_cff, ad_cff1, ad_cff2
203 real(r8) :: adfac, adfac1, adfac2, adfac3, adfac4
204
205 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFe
206 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFe
207 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_UFx
208 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_VFx
209
210#include "set_bounds.h"
211!
212!-----------------------------------------------------------------------
213! Initialize adjoint private variables.
214!-----------------------------------------------------------------------
215!
216 ad_cff=0.0_r8
217 ad_cff1=0.0_r8
218 ad_cff2=0.0_r8
219
220 ad_ufe(imins:imaxs,jmins:jmaxs)=0.0_r8
221 ad_vfe(imins:imaxs,jmins:jmaxs)=0.0_r8
222 ad_ufx(imins:imaxs,jmins:jmaxs)=0.0_r8
223 ad_vfx(imins:imaxs,jmins:jmaxs)=0.0_r8
224!
225!-----------------------------------------------------------------------
226! Compute horizontal harmonic viscosity along constant S-surfaces.
227!-----------------------------------------------------------------------
228!
229 k_loop : DO k=1,n(ng)
230!
231! Time-step harmonic, S-surfaces viscosity term. Notice that momentum
232! at this stage is HzU and HzV and has m2/s units. Add contribution for
233! barotropic forcing terms.
234!
235 DO j=jstrv,jend
236 DO i=istr,iend
237 cff=0.25_r8*(pm(i,j)+pm(i,j-1))*(pn(i,j)+pn(i,j-1))
238!^ tl_v(i,j,k,nnew)=tl_v(i,j,k,nnew)+tl_cff2
239!^
240 ad_cff2=ad_cff2+ad_v(i,j,k,nnew)
241!^ tl_rvfrc(i,j)=tl_rvfrc(i,j)+tl_cff1
242!^
243 ad_cff1=ad_cff1+ad_rvfrc(i,j)
244!^ tl_cff2=dt(ng)*cff*tl_cff1
245!^
246 ad_cff1=ad_cff1+dt(ng)*cff*ad_cff2
247 ad_cff2=0.0_r8
248!^ tl_cff1=0.5_r8*((pn(i,j-1)+pn(i,j))* &
249!^ & (tl_VFx(i+1,j)-tl_VFx(i,j ))- &
250!^ & (pm(i,j-1)+pm(i,j))* &
251!^ & (tl_VFe(i ,j)-tl_VFe(i,j-1)))
252!^
253 adfac=0.5_r8*ad_cff1
254 adfac1=adfac*(pn(i,j-1)+pn(i,j))
255 adfac2=adfac*(pm(i,j-1)+pm(i,j))
256 ad_vfx(i ,j)=ad_vfx(i ,j)-adfac1
257 ad_vfx(i+1,j)=ad_vfx(i+1,j)+adfac1
258 ad_vfe(i,j-1)=ad_vfe(i,j-1)+adfac2
259 ad_vfe(i,j )=ad_vfe(i,j )-adfac2
260 ad_cff1=0.0_r8
261 END DO
262 END DO
263 DO j=jstr,jend
264 DO i=istru,iend
265 cff=0.25_r8*(pm(i-1,j)+pm(i,j))*(pn(i-1,j)+pn(i,j))
266!^ tl_u(i,j,k,nnew)=tl_u(i,j,k,nnew)+tl_cff2
267!^
268 ad_cff2=ad_cff2+ad_u(i,j,k,nnew)
269!^ tl_rufrc(i,j)=tl_rufrc(i,j)+tl_cff1
270!^
271 ad_cff1=ad_cff1+ad_rufrc(i,j)
272!^ tl_cff2=dt(ng)*cff*tl_cff1
273!^
274 ad_cff1=ad_cff1+dt(ng)*cff*ad_cff2
275 ad_cff2=0.0_r8
276!^ cff1=0.5_r8*((pn(i-1,j)+pn(i,j))* &
277!^ & (UFx(i,j )-UFx(i-1,j))+ &
278!^ & (pm(i-1,j)+pm(i,j))* &
279!^ & (UFe(i,j+1)-UFe(i ,j)))
280!^
281!^
282 adfac=0.5_r8*ad_cff1
283 adfac1=adfac*(pn(i-1,j)+pn(i,j))
284 adfac2=adfac*(pm(i-1,j)+pm(i,j))
285 ad_ufx(i-1,j)=ad_ufx(i-1,j)-adfac1
286 ad_ufx(i ,j)=ad_ufx(i ,j)+adfac1
287 ad_ufe(i,j )=ad_ufe(i,j )-adfac2
288 ad_ufe(i,j+1)=ad_ufe(i,j+1)+adfac2
289 ad_cff1=0.0_r8
290 END DO
291 END DO
292!
293! Compute flux-components of the horizontal divergence of the stress
294! tensor (m5/s2) in XI- and ETA-directions.
295!
296 DO j=jstr,jend+1
297 DO i=istr,iend+1
298!^ tl_VFx(i,j)=on_p(i,j)*on_p(i,j)*tl_cff
299!^ tl_UFe(i,j)=om_p(i,j)*om_p(i,j)*tl_cff
300!^
301 ad_cff=ad_cff+ &
302 & on_p(i,j)*on_p(i,j)*ad_vfx(i,j)+ &
303 & om_p(i,j)*om_p(i,j)*ad_ufe(i,j)
304 ad_vfx(i,j)=0.0_r8
305 ad_ufe(i,j)=0.0_r8
306#ifdef MASKING
307!^ tl_cff=tl_cff*pmask(i,j)
308!^
309 ad_cff=ad_cff*pmask(i,j)
310#endif
311!^ tl_cff=visc2_p(i,j)*0.125_r8* &
312!^ & ((tl_Hz(i-1,j ,k)+tl_Hz(i,j ,k)+ &
313!^ & tl_Hz(i-1,j-1,k)+tl_Hz(i,j-1,k))* &
314!^ & (pmon_p(i,j)* &
315!^ & ((pn(i ,j-1)+pn(i ,j))*v(i ,j,k,nrhs)- &
316!^ & (pn(i-1,j-1)+pn(i-1,j))*v(i-1,j,k,nrhs))+ &
317!^ & pnom_p(i,j)* &
318!^ & ((pm(i-1,j )+pm(i,j ))*u(i,j ,k,nrhs)- &
319!^ & (pm(i-1,j-1)+pm(i,j-1))*u(i,j-1,k,nrhs)))+ &
320!^ & (Hz(i-1,j ,k)+Hz(i,j ,k)+ &
321!^ & Hz(i-1,j-1,k)+Hz(i,j-1,k))* &
322!^ & (pmon_p(i,j)* &
323!^ & ((pn(i ,j-1)+pn(i ,j))*tl_v(i ,j,k,nrhs)- &
324!^ & (pn(i-1,j-1)+pn(i-1,j))*tl_v(i-1,j,k,nrhs))+ &
325!^ & pnom_p(i,j)* &
326!^ & ((pm(i-1,j )+pm(i,j ))*tl_u(i,j ,k,nrhs)- &
327!^ & (pm(i-1,j-1)+pm(i,j-1))*tl_u(i,j-1,k,nrhs))))
328!^
329
330 adfac=visc2_p(i,j)*0.125_r8*ad_cff
331 adfac1=adfac* &
332 & (pmon_p(i,j)* &
333 & ((pn(i ,j-1)+pn(i ,j))*v(i ,j,k,nrhs)- &
334 & (pn(i-1,j-1)+pn(i-1,j))*v(i-1,j,k,nrhs))+ &
335 & pnom_p(i,j)* &
336 & ((pm(i-1,j )+pm(i,j ))*u(i,j ,k,nrhs)- &
337 & (pm(i-1,j-1)+pm(i,j-1))*u(i,j-1,k,nrhs)))
338 adfac2=adfac*(hz(i-1,j ,k)+hz(i,j ,k)+ &
339 & hz(i-1,j-1,k)+hz(i,j-1,k))
340 adfac3=adfac2*pmon_p(i,j)
341 adfac4=adfac2*pnom_p(i,j)
342 ad_hz(i-1,j-1,k)=ad_hz(i-1,j-1,k)+adfac1
343 ad_hz(i ,j-1,k)=ad_hz(i ,j-1,k)+adfac1
344 ad_hz(i-1,j ,k)=ad_hz(i-1,j ,k)+adfac1
345 ad_hz(i ,j ,k)=ad_hz(i ,j ,k)+adfac1
346 ad_v(i-1,j,k,nrhs)=ad_v(i-1,j,k,nrhs)- &
347 & (pn(i-1,j-1)+pn(i-1,j))*adfac3
348 ad_v(i ,j,k,nrhs)=ad_v(i ,j,k,nrhs)+ &
349 & (pn(i ,j-1)+pn(i ,j))*adfac3
350 ad_u(i,j-1,k,nrhs)=ad_u(i,j-1,k,nrhs)- &
351 & (pm(i-1,j-1)+pm(i,j-1))*adfac4
352 ad_u(i,j ,k,nrhs)=ad_u(i,j ,k,nrhs)+ &
353 & (pm(i-1,j )+pm(i,j ))*adfac4
354 ad_cff=0.0_r8
355 END DO
356 END DO
357 DO j=jstrv-1,jend
358 DO i=istru-1,iend
359!^ tl_VFe(i,j)=om_r(i,j)*om_r(i,j)*tl_cff
360!^ tl_UFx(i,j)=on_r(i,j)*on_r(i,j)*tl_cff
361!^
362 ad_cff=ad_cff+ &
363 & om_r(i,j)*om_r(i,j)*ad_vfe(i,j)+ &
364 & on_r(i,j)*on_r(i,j)*ad_ufx(i,j)
365 ad_vfe(i,j)=0.0_r8
366 ad_ufx(i,j)=0.0_r8
367!^ tl_cff=visc2_r(i,j)*0.5_r8* &
368!^ & (tl_Hz(i,j,k)* &
369!^ & (pmon_r(i,j)* &
370!^ & ((pn(i ,j)+pn(i+1,j))*u(i+1,j,k,nrhs)- &
371!^ & (pn(i-1,j)+pn(i ,j))*u(i ,j,k,nrhs))- &
372!^ & pnom_r(i,j)* &
373!^ & ((pm(i,j )+pm(i,j+1))*v(i,j+1,k,nrhs)- &
374!^ & (pm(i,j-1)+pm(i,j ))*v(i,j ,k,nrhs)))+ &
375!^ & Hz(i,j,k)* &
376!^ & (pmon_r(i,j)* &
377!^ & ((pn(i ,j)+pn(i+1,j))*tl_u(i+1,j,k,nrhs)- &
378!^ & (pn(i-1,j)+pn(i ,j))*tl_u(i ,j,k,nrhs))- &
379!^ & pnom_r(i,j)* &
380!^ & ((pm(i,j )+pm(i,j+1))*tl_v(i,j+1,k,nrhs)- &
381!^ & (pm(i,j-1)+pm(i,j ))*tl_v(i,j ,k,nrhs))))
382!^
383 adfac=visc2_r(i,j)*0.5_r8*ad_cff
384 adfac1=adfac*hz(i,j,k)
385 adfac2=adfac1*pmon_r(i,j)
386 adfac3=adfac1*pnom_r(i,j)
387 ad_hz(i,j,k)=ad_hz(i,j,k)+ &
388 & (pmon_r(i,j)* &
389 & ((pn(i ,j)+pn(i+1,j))*u(i+1,j,k,nrhs)- &
390 & (pn(i-1,j)+pn(i ,j))*u(i ,j,k,nrhs))- &
391 & pnom_r(i,j)* &
392 & ((pm(i,j )+pm(i,j+1))*v(i,j+1,k,nrhs)- &
393 & (pm(i,j-1)+pm(i,j ))*v(i,j ,k,nrhs)))* &
394 & adfac
395 ad_u(i ,j,k,nrhs)=ad_u(i ,j,k,nrhs)- &
396 & (pn(i-1,j)+pn(i ,j))*adfac2
397 ad_u(i+1,j,k,nrhs)=ad_u(i+1,j,k,nrhs)+ &
398 & (pn(i ,j)+pn(i+1,j))*adfac2
399 ad_v(i,j ,k,nrhs)=ad_v(i,j ,k,nrhs)+ &
400 & (pm(i,j-1)+pm(i,j ))*adfac3
401 ad_v(i,j+1,k,nrhs)=ad_v(i,j+1,k,nrhs)- &
402 & (pm(i,j )+pm(i,j+1))*adfac3
403 ad_cff=0.0_r8
404 END DO
405 END DO
406 END DO k_loop
407!
408 RETURN

References mod_scalars::dt.