ROMS
Loading...
Searching...
No Matches
ad_prsgrd32.h
Go to the documentation of this file.
1 MODULE ad_prsgrd_mod
2!
3!git $Id$
4!================================================== Hernan G. Arango ===
5! Copyright (c) 2002-2025 The ROMS Group Andrew M. Moore !
6! Licensed under a MIT/X style license !
7! See License_ROMS.md !
8!=======================================================================
9! !
10! This sub routine evaluates the adjoint baroclinic, hydrostatic !
11! pressure gradient term using a nonconservative Density-Jacobian !
12! scheme, based on cubic polynomial fits for "rho" and "z_r" as !
13! functions of nondimensional coordinates (XI,ETA,s), that is, its !
14! respective array indices. The cubic polynomials are monotonized !
15! by using harmonic mean instead of linear averages to interpolate !
16! slopes. This scheme retains exact anti-symmetry: !
17! !
18! J(rho,z_r)=-J(z_r,rho). !
19! !
20! If parameter OneFifth (below) is set to zero, the scheme becomes !
21! identical to standard Jacobian. !
22! !
23! Reference: !
24! !
25! Shchepetkin A.F and J.C. McWilliams, 2003: A method for !
26! computing horizontal pressure gradient force in an ocean !
27! model with non-aligned vertical coordinate, JGR, 108, !
28! 1-34. !
29! !
30!=======================================================================
31!
32 implicit none
33!
34 PRIVATE
35 PUBLIC :: ad_prsgrd
36!
37 CONTAINS
38!
39!***********************************************************************
40 SUBROUTINE ad_prsgrd (ng, tile)
41!***********************************************************************
42!
43 USE mod_param
44#ifdef DIAGNOSTICS
45!! USE mod_diags
46#endif
47#ifdef ATM_PRESS
48 USE mod_forces
49#endif
50 USE mod_grid
51 USE mod_ocean
52 USE mod_stepping
53!
54! Imported variable declarations.
55!
56 integer, intent(in) :: ng, tile
57!
58! Local variable declarations.
59!
60 character (len=*), parameter :: MyFile = &
61 & __FILE__
62!
63#include "tile.h"
64!
65#ifdef PROFILE
66 CALL wclock_on (ng, iadm, 23, __line__, myfile)
67#endif
68 CALL ad_prsgrd32_tile (ng, tile, &
69 & lbi, ubi, lbj, ubj, &
70 & imins, imaxs, jmins, jmaxs, &
71 & nrhs(ng), &
72#ifdef MASKING
73 & grid(ng) % umask, &
74 & grid(ng) % vmask, &
75#endif
76 & grid(ng) % om_v, &
77 & grid(ng) % on_u, &
78 & grid(ng) % Hz, &
79 & grid(ng) % ad_Hz, &
80 & grid(ng) % z_r, &
81 & grid(ng) % ad_z_r, &
82 & grid(ng) % z_w, &
83 & grid(ng) % ad_z_w, &
84 & ocean(ng) % rho, &
85 & ocean(ng) % ad_rho, &
86#ifdef TIDE_GENERATING_FORCES
87 & ocean(ng) % eq_tide, &
88 & ocean(ng) % ad_eq_tide, &
89#endif
90#ifdef ATM_PRESS
91 & forces(ng) % Pair, &
92#endif
93#ifdef DIAGNOSTICS_UV
94!! & DIAGS(ng) % DiaRU, &
95!! & DIAGS(ng) % DiaRV, &
96#endif
97 & ocean(ng) % ad_ru, &
98 & ocean(ng) % ad_rv)
99#ifdef PROFILE
100 CALL wclock_off (ng, iadm, 23, __line__, myfile)
101#endif
102!
103 RETURN
104 END SUBROUTINE ad_prsgrd
105!
106!***********************************************************************
107 SUBROUTINE ad_prsgrd32_tile (ng, tile, &
108 & LBi, UBi, LBj, UBj, &
109 & IminS, ImaxS, JminS, JmaxS, &
110 & nrhs, &
111#ifdef MASKING
112 & umask, vmask, &
113#endif
114 & om_v, on_u, &
115 & Hz, ad_Hz, &
116 & z_r, ad_z_r, &
117 & z_w, ad_z_w, &
118 & rho, ad_rho, &
119#ifdef TIDE_GENERATING_FORCES
120 & eq_tide, ad_eq_tide, &
121#endif
122#ifdef ATM_PRESS
123 & Pair, &
124#endif
125#ifdef DIAGNOSTICS_UV
126!! & DiaRU, DiaRV, &
127#endif
128 & ad_ru, ad_rv)
129!***********************************************************************
130!
131 USE mod_param
132 USE mod_scalars
133!
134! Imported variable declarations.
135!
136 integer, intent(in) :: ng, tile
137 integer, intent(in) :: LBi, UBi, LBj, UBj
138 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
139 integer, intent(in) :: nrhs
140
141#ifdef ASSUMED_SHAPE
142# ifdef MASKING
143 real(r8), intent(in) :: umask(LBi:,LBj:)
144 real(r8), intent(in) :: vmask(LBi:,LBj:)
145# endif
146 real(r8), intent(in) :: om_v(LBi:,LBj:)
147 real(r8), intent(in) :: on_u(LBi:,LBj:)
148 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
149 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
150 real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
151 real(r8), intent(in) :: rho(LBi:,LBj:,:)
152# ifdef ATM_PRESS
153 real(r8), intent(in) :: Pair(LBi:,LBj:)
154# endif
155# ifdef TIDE_GENERATING_FORCES
156 real(r8), intent(in) :: eq_tide(LBi:,LBj:)
157 real(r8), intent(inout) :: ad_eq_tide(LBi:,LBj:)
158# endif
159# ifdef DIAGNOSTICS_UV
160!! real(r8), intent(inout) :: DiaRU(LBi:,LBj:,:,:,:)
161!! real(r8), intent(inout) :: DiaRV(LBi:,LBj:,:,:,:)
162# endif
163 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
164 real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:)
165 real(r8), intent(inout) :: ad_z_w(LBi:,LBj:,0:)
166 real(r8), intent(inout) :: ad_rho(LBi:,LBj:,:)
167 real(r8), intent(inout) :: ad_ru(LBi:,LBj:,0:,:)
168 real(r8), intent(inout) :: ad_rv(LBi:,LBj:,0:,:)
169#else
170# ifdef MASKING
171 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
172 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
173# endif
174 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
175 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
176 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
177 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
178 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
179 real(r8), intent(in) :: rho(LBi:UBi,LBj:UBj,N(ng))
180# ifdef ATM_PRESS
181 real(r8), intent(in) :: Pair(LBi:UBi,LBj:UBj)
182# endif
183# ifdef TIDE_GENERATING_FORCES
184 real(r8), intent(in) :: eq_tide(LBi:UBi,LBj:UBj)
185 real(r8), intent(inout) :: ad_eq_tide(LBi:UBi,LBj:UBj)
186# endif
187# ifdef DIAGNOSTICS_UV
188!! real(r8), intent(inout) :: DiaRU(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
189!! real(r8), intent(inout) :: DiaRV(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
190# endif
191 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
192 real(r8), intent(inout) :: ad_z_r(LBi:UBi,LBj:UBj,N(ng))
193 real(r8), intent(inout) :: ad_z_w(LBi:UBi,LBj:UBj,0:N(ng))
194 real(r8), intent(inout) :: ad_rho(LBi:UBi,LBj:UBj,N(ng))
195 real(r8), intent(inout) :: ad_ru(LBi:UBi,LBj:UBj,0:N(ng),2)
196 real(r8), intent(inout) :: ad_rv(LBi:UBi,LBj:UBj,0:N(ng),2)
197#endif
198!
199! Local variable declarations.
200!
201 integer :: i, j, k
202
203 real(r8), parameter :: OneFifth = 0.2_r8
204 real(r8), parameter :: OneTwelfth = 1.0_r8/12.0_r8
205 real(r8), parameter :: eps = 1.0e-10_r8
206
207 real(r8) :: GRho, GRho0, HalfGRho
208 real(r8) :: cff, cff1, cff2
209#ifdef ATM_PRESS
210 real(r8) :: OneAtm, fac
211#endif
212 real(r8) :: ad_cff, ad_cff1, ad_cff2, adfac
213 real(r8) :: adfac1, adfac2, adfac3, adfac4, adfac5, adfac6
214 real(r8) :: adfac7, adfac8, adfac9, adfac10, adfac11, adfac12
215
216 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: P
217
218 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: ad_P
219
220 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dR
221 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dR1
222 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dZ
223 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dZ1
224
225 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_dR
226 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_dZ
227
228 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FC
229 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: aux
230 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: dRx
231 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: dZx
232
233 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_FC
234 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_aux
235 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_dRx
236 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: ad_dZx
237
238#include "set_bounds.h"
239!
240!-----------------------------------------------------------------------
241! Initialize adjoint private variables.
242!-----------------------------------------------------------------------
243!
244 ad_cff=0.0_r8
245 ad_cff1=0.0_r8
246 ad_cff2=0.0_r8
247 DO j=jmins,jmaxs
248 DO i=imins,imaxs
249 ad_fc(i,j)=0.0_r8
250 ad_aux(i,j)=0.0_r8
251 ad_drx(i,j)=0.0_r8
252 ad_dzx(i,j)=0.0_r8
253 END DO
254 DO k=1,n(ng)
255 DO i=imins,imaxs
256 ad_p(i,j,k)=0.0_r8
257 END DO
258 END DO
259 END DO
260 DO k=0,n(ng)
261 DO i=imins,imaxs
262 ad_dr(i,k)=0.0_r8
263 ad_dz(i,k)=0.0_r8
264 END DO
265 END DO
266!
267!-----------------------------------------------------------------------
268! Preliminary step (same for XI- and ETA-components):
269!-----------------------------------------------------------------------
270!
271! Compute BASIC STATE dynamic pressure, P.
272!
273 grho=g/rho0
274 grho0=1000.0_r8*grho
275 halfgrho=0.5_r8*grho
276#ifdef ATM_PRESS
277 oneatm=1013.25_r8 ! 1 atm = 1013.25 mb
278 fac=100.0_r8/rho0
279#endif
280!
281 DO j=jstrv-1,jend
282 DO k=1,n(ng)-1
283 DO i=istru-1,iend
284 dr(i,k)=rho(i,j,k+1)-rho(i,j,k)
285 dz(i,k)=z_r(i,j,k+1)-z_r(i,j,k)
286 END DO
287 END DO
288 DO i=istru-1,iend
289 dr(i,n(ng))=dr(i,n(ng)-1)
290 dz(i,n(ng))=dz(i,n(ng)-1)
291 dr(i,0)=dr(i,1)
292 dz(i,0)=dz(i,1)
293 END DO
294 DO k=n(ng),1,-1
295 DO i=istru-1,iend
296 cff=2.0_r8*dr(i,k)*dr(i,k-1)
297 IF (cff.gt.eps) THEN
298 dr(i,k)=cff/(dr(i,k)+dr(i,k-1))
299 ELSE
300 dr(i,k)=0.0_r8
301 END IF
302 dz(i,k)=2.0_r8*dz(i,k)*dz(i,k-1)/(dz(i,k)+dz(i,k-1))
303 END DO
304 END DO
305 DO i=istru-1,iend
306 cff1=1.0_r8/(z_r(i,j,n(ng))-z_r(i,j,n(ng)-1))
307 cff2=0.5_r8*(rho(i,j,n(ng))-rho(i,j,n(ng)-1))* &
308 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))*cff1
309 p(i,j,n(ng))=g*z_w(i,j,n(ng))+ &
310#ifdef ATM_PRESS
311 & fac*(pair(i,j)-oneatm)+ &
312#endif
313 & grho*(rho(i,j,n(ng))+cff2)* &
314 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))
315#ifdef TIDE_GENERATING_FORCES
316 p(i,j,n(ng))=p(i,j,n(ng))-g*eq_tide(i,j)
317#endif
318 END DO
319 DO k=n(ng)-1,1,-1
320 DO i=istru-1,iend
321 p(i,j,k)=p(i,j,k+1)+ &
322 & halfgrho*((rho(i,j,k+1)+rho(i,j,k))* &
323 & (z_r(i,j,k+1)-z_r(i,j,k))- &
324 & onefifth* &
325 & ((dr(i,k+1)-dr(i,k))* &
326 & (z_r(i,j,k+1)-z_r(i,j,k)- &
327 & onetwelfth* &
328 & (dz(i,k+1)+dz(i,k)))- &
329 & (dz(i,k+1)-dz(i,k))* &
330 & (rho(i,j,k+1)-rho(i,j,k)- &
331 & onetwelfth* &
332 & (dr(i,k+1)+dr(i,k)))))
333 END DO
334 END DO
335 END DO
336!
337!-----------------------------------------------------------------------
338! Adjoint ETA-component pressure gradient term.
339!-----------------------------------------------------------------------
340!
341 DO k=1,n(ng)
342!
343! Compute BASIC STATE variables.
344!
345 DO j=jstrv-1,jend+1
346 DO i=istr,iend
347 aux(i,j)=z_r(i,j,k)-z_r(i,j-1,k)
348# ifdef MASKING
349 aux(i,j)=aux(i,j)*vmask(i,j)
350# endif
351 fc(i,j)=rho(i,j,k)-rho(i,j-1,k)
352# ifdef MASKING
353 fc(i,j)=fc(i,j)*vmask(i,j)
354# endif
355 END DO
356 END DO
357!
358 DO j=jstrv-1,jend
359 DO i=istr,iend
360 cff=2.0_r8*aux(i,j)*aux(i,j+1)
361 IF (cff.gt.eps) THEN
362 cff1=1.0_r8/(aux(i,j)+aux(i,j+1))
363 dzx(i,j)=cff*cff1
364 ELSE
365 dzx(i,j)=0.0_r8
366 END IF
367 cff1=2.0_r8*fc(i,j)*fc(i,j+1)
368 IF (cff1.gt.eps) THEN
369 cff2=1.0_r8/(fc(i,j)+fc(i,j+1))
370 drx(i,j)=cff1*cff2
371 ELSE
372 drx(i,j)=0.0_r8
373 END IF
374 END DO
375 END DO
376!
377 DO j=jstrv,jend
378 DO i=istr,iend
379#ifdef DIAGNOSTICS_UV
380!! DiaRV(i,j,k,nrhs,M3pgrd)=rv(i,j,k,nrhs)
381#endif
382!^ tl_rv(i,j,k,nrhs)=om_v(i,j)*0.5_r8* &
383!^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* &
384!^ & (P(i,j-1,k)-P(i,j,k)- &
385!^ & HalfGRho* &
386!^ & ((rho(i,j,k)+rho(i,j-1,k))* &
387!^ & (z_r(i,j,k)-z_r(i,j-1,k))- &
388!^ & OneFifth* &
389!^ & ((dRx(i,j)-dRx(i,j-1))* &
390!^ & (z_r(i,j,k)-z_r(i,j-1,k)- &
391!^ & OneTwelfth* &
392!^ & (dZx(i,j)+dZx(i,j-1)))- &
393!^ & (dZx(i,j)-dZx(i,j-1))* &
394!^ & (rho(i,j,k)-rho(i,j-1,k)- &
395!^ & OneTwelfth* &
396!^ & (dRx(i,j)+dRx(i,j-1))))))+ &
397!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
398!^ & (tl_P(i,j-1,k)-tl_P(i,j,k)- &
399!^ & HalfGRho* &
400!^ & ((tl_rho(i,j,k)+tl_rho(i,j-1,k))* &
401!^ & (z_r(i,j,k)-z_r(i,j-1,k))+ &
402!^ & (rho(i,j,k)+rho(i,j-1,k))* &
403!^ & (tl_z_r(i,j,k)-tl_z_r(i,j-1,k))- &
404!^ & OneFifth* &
405!^ & ((tl_dRx(i,j)-tl_dRx(i,j-1))* &
406!^ & (z_r(i,j,k)-z_r(i,j-1,k)- &
407!^ & OneTwelfth* &
408!^ & (dZx(i,j)+dZx(i,j-1)))+ &
409!^ & (dRx(i,j)-dRx(i,j-1))* &
410!^ & (tl_z_r(i,j,k)-tl_z_r(i,j-1,k)- &
411!^ & OneTwelfth* &
412!^ & (tl_dZx(i,j)+tl_dZx(i,j-1)))- &
413!^ & (tl_dZx(i,j)-tl_dZx(i,j-1))* &
414!^ & (rho(i,j,k)-rho(i,j-1,k)- &
415!^ & OneTwelfth* &
416!^ & (dRx(i,j)+dRx(i,j-1)))- &
417!^ & (dZx(i,j)-dZx(i,j-1))* &
418!^ & (tl_rho(i,j,k)-tl_rho(i,j-1,k)- &
419!^ & OneTwelfth* &
420!^ & (tl_dRx(i,j)+tl_dRx(i,j-1)))))))
421!^
422 adfac=om_v(i,j)*0.5_r8*ad_rv(i,j,k,nrhs)
423 adfac1=adfac*(p(i,j-1,k)-p(i,j,k)- &
424 & halfgrho* &
425 & ((rho(i,j,k)+rho(i,j-1,k))* &
426 & (z_r(i,j,k)-z_r(i,j-1,k))- &
427 & onefifth* &
428 & ((drx(i,j)-drx(i,j-1))* &
429 & (z_r(i,j,k)-z_r(i,j-1,k)- &
430 & onetwelfth* &
431 & (dzx(i,j)+dzx(i,j-1)))- &
432 & (dzx(i,j)-dzx(i,j-1))* &
433 & (rho(i,j,k)-rho(i,j-1,k)- &
434 & onetwelfth* &
435 & (drx(i,j)+drx(i,j-1))))))
436 adfac2=adfac*(hz(i,j,k)+hz(i,j-1,k))
437 adfac3=adfac2*halfgrho
438 adfac4=adfac3*(z_r(i,j,k)-z_r(i,j-1,k))
439 adfac5=adfac3*(rho(i,j,k)+rho(i,j-1,k))
440 adfac6=adfac3*onefifth
441 adfac7=adfac6*(z_r(i,j,k)-z_r(i,j-1,k)- &
442 & onetwelfth*(dzx(i,j)+dzx(i,j-1)))
443 adfac8=adfac6*(drx(i,j)-drx(i,j-1))
444 adfac9=adfac8*onetwelfth
445 adfac10=adfac6*(rho(i,j,k)-rho(i,j-1,k)- &
446 & onetwelfth*(drx(i,j)+drx(i,j-1)))
447 adfac11=adfac6*(dzx(i,j)-dzx(i,j-1))
448 adfac12=adfac11*onetwelfth
449 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac1
450 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac1
451 ad_p(i,j-1,k)=ad_p(i,j-1,k)+adfac2
452 ad_p(i,j ,k)=ad_p(i,j ,k)-adfac2
453 ad_rho(i,j-1,k)=ad_rho(i,j-1,k)-adfac4+adfac11
454 ad_rho(i,j ,k)=ad_rho(i,j ,k)-adfac4-adfac11
455 ad_z_r(i,j-1,k)=ad_z_r(i,j-1,k)+adfac5-adfac8
456 ad_z_r(i,j ,k)=ad_z_r(i,j ,k)-adfac5+adfac8
457 ad_drx(i,j-1)=ad_drx(i,j-1)-adfac7+adfac12
458 ad_drx(i,j )=ad_drx(i,j )+adfac7+adfac12
459 ad_dzx(i,j-1)=ad_dzx(i,j-1)-adfac9+adfac10
460 ad_dzx(i,j )=ad_dzx(i,j )-adfac9-adfac10
461 ad_rv(i,j,k,nrhs)=0.0_r8
462 END DO
463 END DO
464!
465 DO j=jstrv-1,jend
466 DO i=istr,iend
467 cff1=2.0_r8*fc(i,j)*fc(i,j+1)
468 IF (cff1.gt.eps) THEN
469 cff2=1.0_r8/(fc(i,j)+fc(i,j+1))
470!^ tl_dRx(i,j)=tl_cff1*cff2+cff1*tl_cff2
471!^
472 ad_cff1=ad_cff1+cff2*ad_drx(i,j)
473 ad_cff2=ad_cff2+cff1*ad_drx(i,j)
474 ad_drx(i,j)=0.0_r8
475!^ tl_cff2=-cff2*cff2*(tl_FC(i,j)+tl_FC(i,j+1))
476!^
477 adfac=-cff2*cff2*ad_cff2
478 ad_fc(i,j )=ad_fc(i,j )+adfac
479 ad_fc(i,j+1)=ad_fc(i,j+1)+adfac
480 ad_cff2=0.0_r8
481 ELSE
482!^ tl_dRx(i,j)=0.0_r8
483!^
484 ad_drx(i,j)=0.0_r8
485 END IF
486!^ tl_cff1=2.0_r8*(tl_FC(i,j)*FC(i,j+1)+ &
487!^ & FC(i,j)*tl_FC(i,j+1))
488!^
489 adfac=2.0_r8*ad_cff1
490 ad_fc(i,j )=ad_fc(i,j )+fc(i,j+1)*adfac
491 ad_fc(i,j+1)=ad_fc(i,j+1)+fc(i,j )*adfac
492 ad_cff1=0.0_r8
493
494 cff=2.0_r8*aux(i,j)*aux(i,j+1)
495 IF (cff.gt.eps) THEN
496 cff1=1.0_r8/(aux(i,j)+aux(i,j+1))
497!^ tl_dZx(i,j)=tl_cff*cff1+cff*tl_cff1
498!^
499 ad_cff=ad_cff+cff1*ad_dzx(i,j)
500 ad_cff1=ad_cff1+cff*ad_dzx(i,j)
501 ad_dzx(i,j)=0.0_r8
502!^ tl_cff1=-cff1*cff1*(tl_aux(i,j)+tl_aux(i,j+1))
503!^
504 adfac=-cff1*cff1*ad_cff1
505 ad_aux(i,j )=ad_aux(i,j )+adfac
506 ad_aux(i,j+1)=ad_aux(i,j+1)+adfac
507 ad_cff1=0.0_r8
508 ELSE
509!^ tl_dZx(i,j)=0.0_r8
510!^
511 ad_dzx(i,j)=0.0_r8
512 END IF
513!^ tl_cff=2.0_r8*(tl_aux(i,j)*aux(i,j+1)+ &
514!^ & aux(i,j)*tl_aux(i,j+1))
515!^
516 adfac=2.0_r8*ad_cff
517 ad_aux(i,j )=ad_aux(i,j )+aux(i,j+1)*adfac
518 ad_aux(i,j+1)=ad_aux(i,j+1)+aux(i,j )*adfac
519 ad_cff=0.0_r8
520 END DO
521 END DO
522 DO j=jstrv-1,jend+1
523 DO i=istr,iend
524#ifdef MASKING
525!^ tl_FC(i,j)=tl_FC(i,j)*vmask(i,j)
526!^
527 ad_fc(i,j)=ad_fc(i,j)*vmask(i,j)
528#endif
529!^ tl_FC(i,j)=tl_rho(i,j,k)-tl_rho(i,j-1,k)
530!^
531 ad_rho(i,j-1,k)=ad_rho(i,j-1,k)-ad_fc(i,j)
532 ad_rho(i,j, k)=ad_rho(i,j ,k)+ad_fc(i,j)
533 ad_fc(i,j)=0.0_r8
534#ifdef MASKING
535!^ tl_aux(i,j)=tl_aux(i,j)*vmask(i,j)
536!^
537 ad_aux(i,j)=ad_aux(i,j)*vmask(i,j)
538#endif
539!^ tl_aux(i,j)=tl_z_r(i,j,k)-tl_z_r(i,j-1,k)
540!^
541 ad_z_r(i,j-1,k)=ad_z_r(i,j-1,k)-ad_aux(i,j)
542 ad_z_r(i,j ,k)=ad_z_r(i,j ,k)+ad_aux(i,j)
543 ad_aux(i,j)=0.0_r8
544 END DO
545 END DO
546 END DO
547!
548!-----------------------------------------------------------------------
549! Compute adjoint XI-component pressure gradient term.
550!-----------------------------------------------------------------------
551!
552 DO k=1,n(ng)
553!
554! Compute BASIC STATE variables.
555!
556 DO j=jstr,jend
557 DO i=istru-1,iend+1
558 aux(i,j)=z_r(i,j,k)-z_r(i-1,j,k)
559# ifdef MASKING
560 aux(i,j)=aux(i,j)*umask(i,j)
561# endif
562 fc(i,j)=rho(i,j,k)-rho(i-1,j,k)
563# ifdef MASKING
564 fc(i,j)=fc(i,j)*umask(i,j)
565# endif
566 END DO
567 END DO
568!
569 DO j=jstr,jend
570 DO i=istru-1,iend
571 cff=2.0_r8*aux(i,j)*aux(i+1,j)
572 IF (cff.gt.eps) THEN
573 cff1=1.0_r8/(aux(i,j)+aux(i+1,j))
574 dzx(i,j)=cff*cff1
575 ELSE
576 dzx(i,j)=0.0_r8
577 END IF
578 cff1=2.0_r8*fc(i,j)*fc(i+1,j)
579 IF (cff1.gt.eps) THEN
580 cff2=1.0_r8/(fc(i,j)+fc(i+1,j))
581 drx(i,j)=cff1*cff2
582 ELSE
583 drx(i,j)=0.0_r8
584 END IF
585 END DO
586 END DO
587!
588 DO j=jstr,jend
589 DO i=istru,iend
590#ifdef DIAGNOSTICS_UV
591!! DiaRU(i,j,k,nrhs,M3pgrd)=ru(i,j,k,nrhs)
592#endif
593!^ tl_ru(i,j,k,nrhs)=on_u(i,j)*0.5_r8* &
594!^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* &
595!^ & (P(i-1,j,k)-P(i,j,k)- &
596!^ & HalfGRho* &
597!^ & ((rho(i,j,k)+rho(i-1,j,k))* &
598!^ & (z_r(i,j,k)-z_r(i-1,j,k))- &
599!^ & OneFifth* &
600!^ & ((dRx(i,j)-dRx(i-1,j))* &
601!^ & (z_r(i,j,k)-z_r(i-1,j,k)- &
602!^ & OneTwelfth* &
603!^ & (dZx(i,j)+dZx(i-1,j)))- &
604!^ & (dZx(i,j)-dZx(i-1,j))* &
605!^ & (rho(i,j,k)-rho(i-1,j,k)- &
606!^ & OneTwelfth* &
607!^ & (dRx(i,j)+dRx(i-1,j))))))+ &
608!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
609!^ & (tl_P(i-1,j,k)-tl_P(i,j,k)- &
610!^ & HalfGRho* &
611!^ & ((tl_rho(i,j,k)+tl_rho(i-1,j,k))* &
612!^ & (z_r(i,j,k)-z_r(i-1,j,k))+ &
613!^ & (rho(i,j,k)+rho(i-1,j,k))* &
614!^ & (tl_z_r(i,j,k)-tl_z_r(i-1,j,k))- &
615!^ & OneFifth* &
616!^ & ((tl_dRx(i,j)-tl_dRx(i-1,j))* &
617!^ & (z_r(i,j,k)-z_r(i-1,j,k)- &
618!^ & OneTwelfth* &
619!^ & (dZx(i,j)+dZx(i-1,j)))+ &
620!^ & (dRx(i,j)-dRx(i-1,j))* &
621!^ & (tl_z_r(i,j,k)-tl_z_r(i-1,j,k)- &
622!^ & OneTwelfth* &
623!^ & (tl_dZx(i,j)+tl_dZx(i-1,j)))- &
624!^ & (tl_dZx(i,j)-tl_dZx(i-1,j))* &
625!^ & (rho(i,j,k)-rho(i-1,j,k)- &
626!^ & OneTwelfth* &
627!^ & (dRx(i,j)+dRx(i-1,j)))- &
628!^ & (dZx(i,j)-dZx(i-1,j))* &
629!^ & (tl_rho(i,j,k)-tl_rho(i-1,j,k)- &
630!^ & OneTwelfth* &
631!^ & (tl_dRx(i,j)+tl_dRx(i-1,j)))))))
632!^
633 adfac=on_u(i,j)*0.5_r8*ad_ru(i,j,k,nrhs)
634 adfac1=adfac*(p(i-1,j,k)-p(i,j,k)- &
635 & halfgrho* &
636 & ((rho(i,j,k)+rho(i-1,j,k))* &
637 & (z_r(i,j,k)-z_r(i-1,j,k))- &
638 & onefifth* &
639 & ((drx(i,j)-drx(i-1,j))* &
640 & (z_r(i,j,k)-z_r(i-1,j,k)- &
641 & onetwelfth* &
642 & (dzx(i,j)+dzx(i-1,j)))- &
643 & (dzx(i,j)-dzx(i-1,j))* &
644 & (rho(i,j,k)-rho(i-1,j,k)- &
645 & onetwelfth* &
646 & (drx(i,j)+drx(i-1,j))))))
647 adfac2=adfac*(hz(i,j,k)+hz(i-1,j,k))
648 adfac3=adfac2*halfgrho
649 adfac4=adfac3*(z_r(i,j,k)-z_r(i-1,j,k))
650 adfac5=adfac3*(rho(i,j,k)+rho(i-1,j,k))
651 adfac6=adfac3*onefifth
652 adfac7=adfac6*(z_r(i,j,k)-z_r(i-1,j,k)- &
653 & onetwelfth*(dzx(i,j)+dzx(i-1,j)))
654 adfac8=adfac6*(drx(i,j)-drx(i-1,j))
655 adfac9=adfac8*onetwelfth
656 adfac10=adfac6*(rho(i,j,k)-rho(i-1,j,k)- &
657 & onetwelfth*(drx(i,j)+drx(i-1,j)))
658 adfac11=adfac6*(dzx(i,j)-dzx(i-1,j))
659 adfac12=adfac11*onetwelfth
660 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac1
661 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac1
662 ad_p(i-1,j,k)=ad_p(i-1,j,k)+adfac2
663 ad_p(i ,j,k)=ad_p(i ,j,k)-adfac2
664 ad_rho(i-1,j,k)=ad_rho(i-1,j,k)-adfac4+adfac11
665 ad_rho(i ,j,k)=ad_rho(i ,j,k)-adfac4-adfac11
666 ad_z_r(i-1,j,k)=ad_z_r(i-1,j,k)+adfac5-adfac8
667 ad_z_r(i ,j,k)=ad_z_r(i ,j,k)-adfac5+adfac8
668 ad_drx(i-1,j)=ad_drx(i-1,j)-adfac7+adfac12
669 ad_drx(i ,j)=ad_drx(i ,j)+adfac7+adfac12
670 ad_dzx(i-1,j)=ad_dzx(i-1,j)-adfac9+adfac10
671 ad_dzx(i ,j)=ad_dzx(i ,j)-adfac9-adfac10
672 ad_ru(i,j,k,nrhs)=0.0_r8
673 END DO
674 END DO
675!
676 DO j=jstr,jend
677 DO i=istru-1,iend
678 cff1=2.0_r8*fc(i,j)*fc(i+1,j)
679 IF (cff1.gt.eps) THEN
680 cff2=1.0_r8/(fc(i,j)+fc(i+1,j))
681!^ tl_dRx(i,j)=tl_cff1*cff2+cff1*tl_cff2
682!^
683 ad_cff1=ad_cff1+cff2*ad_drx(i,j)
684 ad_cff2=ad_cff2+cff1*ad_drx(i,j)
685 ad_drx(i,j)=0.0_r8
686!^ tl_cff2=-cff2*cff2*(tl_FC(i,j)+tl_FC(i+1,j))
687!^
688 adfac=-cff2*cff2*ad_cff2
689 ad_fc(i ,j)=ad_fc(i ,j)+adfac
690 ad_fc(i+1,j)=ad_fc(i+1,j)+adfac
691 ad_cff2=0.0_r8
692 ELSE
693!^ tl_dRx(i,j)=0.0_r8
694!^
695 ad_drx(i,j)=0.0_r8
696 END IF
697!^ tl_cff1=2.0_r8*(tl_FC(i,j)*FC(i+1,j)+ &
698!^ & FC(i,j)*tl_FC(i+1,j)
699!^
700 adfac=2.0_r8*ad_cff1
701 ad_fc(i ,j)=ad_fc(i ,j)+fc(i+1,j)*adfac
702 ad_fc(i+1,j)=ad_fc(i+1,j)+fc(i ,j)*adfac
703 ad_cff1=0.0_r8
704
705 cff=2.0_r8*aux(i,j)*aux(i+1,j)
706 IF (cff.gt.eps) THEN
707 cff1=1.0_r8/(aux(i,j)+aux(i+1,j))
708!^ tl_dZx(i,j)=tl_cff*cff1+cff*tl_cff1
709!^
710 ad_cff=ad_cff+cff1*ad_dzx(i,j)
711 ad_cff1=ad_cff1+cff*ad_dzx(i,j)
712 ad_dzx(i,j)=0.0_r8
713!^ tl_cff1=-cff1*cff1*(tl_aux(i,j)+tl_aux(i+1,j))
714!^
715 adfac=-cff1*cff1*ad_cff1
716 ad_aux(i ,j)=ad_aux(i ,j)+adfac
717 ad_aux(i+1,j)=ad_aux(i+1,j)+adfac
718 ad_cff1=0.0_r8
719 ELSE
720!^ tl_dZx(i,j)=0.0_r8
721!^
722 ad_dzx(i,j)=0.0_r8
723 END IF
724!^ tl_cff=2.0_r8*(tl_aux(i,j)*aux(i+1,j)+ &
725!^ & aux(i,j)*tl_aux(i+1,j)
726!^
727 adfac=2.0_r8*ad_cff
728 ad_aux(i ,j)=ad_aux(i ,j)+aux(i+1,j)*adfac
729 ad_aux(i+1,j)=ad_aux(i+1,j)+aux(i ,j)*adfac
730 ad_cff=0.0_r8
731 END DO
732 END DO
733 DO j=jstr,jend
734 DO i=istru-1,iend+1
735#ifdef MASKING
736!^ tl_FC(i,j)=tl_FC(i,j)*umask(i,j)
737!^
738 ad_fc(i,j)=ad_fc(i,j)*umask(i,j)
739#endif
740!^ tl_FC(i,j)=tl_rho(i,j,k)-tl_rho(i-1,j,k)
741!^
742 ad_rho(i-1,j,k)=ad_rho(i-1,j,k)-ad_fc(i,j)
743 ad_rho(i ,j,k)=ad_rho(i ,j,k)+ad_fc(i,j)
744 ad_fc(i,j)=0.0_r8
745#ifdef MASKING
746!^ tl_aux(i,j)=tl_aux(i,j)*umask(i,j)
747!^
748 ad_aux(i,j)=ad_aux(i,j)*umask(i,j)
749#endif
750!^ tl_aux(i,j)=tl_z_r(i,j,k)-tl_z_r(i-1,j,k)
751!^
752 ad_z_r(i-1,j,k)=ad_z_r(i-1,j,k)-ad_aux(i,j)
753 ad_z_r(i ,j,k)=ad_z_r(i ,j,k)+ad_aux(i,j)
754 ad_aux(i,j)=0.0_r8
755 END DO
756 END DO
757 END DO
758!
759!-----------------------------------------------------------------------
760! Adjoint of Preliminary step (same for XI- and ETA-components):
761!-----------------------------------------------------------------------
762!
763 DO j=jstrv-1,jend
764 DO k=1,n(ng)-1
765 DO i=istru-1,iend
766 dr(i,k)=rho(i,j,k+1)-rho(i,j,k)
767 dz(i,k)=z_r(i,j,k+1)-z_r(i,j,k)
768 dr1(i,k)=dr(i,k)
769 dz1(i,k)=dz(i,k)
770 END DO
771 END DO
772 DO i=istru-1,iend
773 dr(i,n(ng))=dr(i,n(ng)-1)
774 dr1(i,n(ng))=dr(i,n(ng))
775 dz(i,n(ng))=dz(i,n(ng)-1)
776 dz1(i,n(ng))=dz(i,n(ng))
777 dr(i,0)=dr(i,1)
778 dr1(i,0)=dr(i,0)
779 dz(i,0)=dz(i,1)
780 dz1(i,0)=dz(i,0)
781 END DO
782 DO k=n(ng),1,-1
783 DO i=istru-1,iend
784 cff=2.0_r8*dr(i,k)*dr(i,k-1)
785 IF (cff.gt.eps) THEN
786 dr(i,k)=cff/(dr(i,k)+dr(i,k-1))
787 ELSE
788 dr(i,k)=0.0_r8
789 END IF
790 dz(i,k)=2.0_r8*dz(i,k)*dz(i,k-1)/(dz(i,k)+dz(i,k-1))
791 END DO
792 END DO
793!
794 DO k=1,n(ng)-1
795 DO i=istru-1,iend
796!^ tl_P(i,j,k)=tl_P(i,j,k+1)+tl_cff
797!^
798 ad_cff=ad_cff+ad_p(i,j,k)
799!^ tl_cff=HalfGRho*((tl_rho(i,j,k+1)+tl_rho(i,j,k))* &
800!^ & (z_r(i,j,k+1)-z_r(i,j,k))+ &
801!^ & (rho(i,j,k+1)+rho(i,j,k))* &
802!^ & (tl_z_r(i,j,k+1)-tl_z_r(i,j,k))- &
803!^ & OneFifth* &
804!^ & ((tl_dR(i,k+1)-tl_dR(i,k))* &
805!^ & (z_r(i,j,k+1)-z_r(i,j,k)- &
806!^ & OneTwelfth* &
807!^ & (dZ(i,k+1)+dZ(i,k)))+ &
808!^ & (dR(i,k+1)-dR(i,k))* &
809!^ & (tl_z_r(i,j,k+1)-tl_z_r(i,j,k)- &
810!^ & OneTwelfth* &
811!^ & (tl_dZ(i,k+1)+tl_dZ(i,k)))- &
812!^ & (tl_dZ(i,k+1)-tl_dZ(i,k))* &
813!^ & (rho(i,j,k+1)-rho(i,j,k)- &
814!^ & OneTwelfth* &
815!^ & (dR(i,k+1)+dR(i,k)))- &
816!^ & (dZ(i,k+1)-dZ(i,k))* &
817!^ & (tl_rho(i,j,k+1)-tl_rho(i,j,k)- &
818!^ & OneTwelfth* &
819!^ & (tl_dR(i,k+1)+tl_dR(i,k)))))
820!^
821 adfac=halfgrho*ad_cff
822 adfac1=adfac*(z_r(i,j,k+1)-z_r(i,j,k))
823 adfac2=adfac*(rho(i,j,k+1)+rho(i,j,k))
824 adfac3=adfac*onefifth
825 adfac4=adfac3*(z_r(i,j,k+1)-z_r(i,j,k)- &
826 & onetwelfth*(dz(i,k+1)+dz(i,k)))
827 adfac5=adfac3*(dr(i,k+1)-dr(i,k))
828 adfac6=adfac5*onetwelfth
829 adfac7=adfac3*(rho(i,j,k+1)-rho(i,j,k)- &
830 & onetwelfth*(dr(i,k+1)+dr(i,k)))
831 adfac8=adfac3*(dz(i,k+1)-dz(i,k))
832 adfac9=adfac8*onetwelfth
833 ad_p(i,j,k+1)=ad_p(i,j,k+1)+ad_p(i,j,k)
834 ad_rho(i,j,k )=ad_rho(i,j,k )+adfac1-adfac8
835 ad_rho(i,j,k+1)=ad_rho(i,j,k+1)+adfac1+adfac8
836 ad_z_r(i,j,k )=ad_z_r(i,j,k )-adfac2+adfac5
837 ad_z_r(i,j,k+1)=ad_z_r(i,j,k+1)+adfac2-adfac5
838 ad_dr(i,k )=ad_dr(i,k )+adfac4-adfac9
839 ad_dr(i,k+1)=ad_dr(i,k+1)-adfac4-adfac9
840 ad_dz(i,k )=ad_dz(i,k )+adfac6-adfac7
841 ad_dz(i,k+1)=ad_dz(i,k+1)+adfac6+adfac7
842 ad_cff=0.0_r8
843 END DO
844 END DO
845 DO i=istru-1,iend
846 cff1=1.0_r8/(z_r(i,j,n(ng))-z_r(i,j,n(ng)-1))
847 cff2=0.5_r8*(rho(i,j,n(ng))-rho(i,j,n(ng)-1))* &
848 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))*cff1
849#ifdef TIDE_GENERATING_FORCES
850!^ tl_P(i,j,N(ng))=tl_P(i,j,N(ng))-g*tl_eq_tide(i,j)
851!^
852 ad_eq_tide(i,j)=ad_eq_tide(i,j)-g*ad_p(i,j,n(ng))
853#endif
854!^ tl_P(i,j,N(ng))=g*tl_z_w(i,j,N(ng))+ &
855!^ & GRho*((tl_rho(i,j,N(ng))+tl_cff2)* &
856!^ & (z_w(i,j,N(ng))-z_r(i,j,N(ng)))+ &
857!^ & (rho(i,j,N(ng))+cff2)* &
858!^ & (tl_z_w(i,j,N(ng))-tl_z_r(i,j,N(ng))))
859!^
860 adfac=grho*ad_p(i,j,n(ng))
861 adfac1=adfac*(z_w(i,j,n(ng))-z_r(i,j,n(ng)))
862 adfac2=adfac*(rho(i,j,n(ng))+cff2)
863 ad_z_r(i,j,n(ng))=ad_z_r(i,j,n(ng))-adfac2
864 ad_z_w(i,j,n(ng))=ad_z_w(i,j,n(ng))+adfac2+ &
865 & g*ad_p(i,j,n(ng))
866 ad_rho(i,j,n(ng))=ad_rho(i,j,n(ng))+adfac1
867 ad_cff2=ad_cff2+adfac1
868 ad_p(i,j,n(ng))=0.0_r8
869!^ tl_cff2=0.5_r8*((tl_rho(i,j,N(ng))-tl_rho(i,j,N(ng)-1))* &
870!^ & (z_w(i,j,N(ng))-z_r(i,j,N(ng)))*cff1+ &
871!^ & (rho(i,j,N(ng))-rho(i,j,N(ng)-1))* &
872!^ & ((tl_z_w(i,j,N(ng))-tl_z_r(i,j,N(ng)))*cff1+ &
873!^ & (z_w(i,j,N(ng))-z_r(i,j,N(ng)))*tl_cff1))
874!^
875 adfac=0.5_r8*ad_cff2
876 adfac1=adfac*(z_w(i,j,n(ng))-z_r(i,j,n(ng)))*cff1
877 adfac2=adfac*(rho(i,j,n(ng))-rho(i,j,n(ng)-1))
878 adfac3=adfac2*cff1
879 ad_rho(i,j,n(ng)-1)=ad_rho(i,j,n(ng)-1)-adfac1
880 ad_rho(i,j,n(ng) )=ad_rho(i,j,n(ng) )+adfac1
881 ad_z_r(i,j,n(ng))=ad_z_r(i,j,n(ng))-adfac3
882 ad_z_w(i,j,n(ng))=ad_z_w(i,j,n(ng))+adfac3
883 ad_cff1=ad_cff1+(z_w(i,j,n(ng))-z_r(i,j,n(ng)))*adfac2
884 ad_cff2=0.0_r8
885!^ tl_cff1=-cff1*cff1*(tl_z_r(i,j,N(ng))-tl_z_r(i,j,N(ng)-1))
886!^
887 adfac=-cff1*cff1*ad_cff1
888 ad_z_r(i,j,n(ng)-1)=ad_z_r(i,j,n(ng)-1)-adfac
889 ad_z_r(i,j,n(ng) )=ad_z_r(i,j,n(ng) )+adfac
890 ad_cff1=0.0_r8
891 END DO
892!
893! The forward code has recursive statements, so we need to dR1(ik) and
894! dZ1(i,k) for BOTH terms in denominator of adjoint code.
895!
896 DO k=1,n(ng)
897 DO i=istru-1,iend
898!^ tl_dZ(i,k)=(2.0_r8*(tl_dZ(i,k)*dZ1(i,k-1)+ &
899!^ & dZ1(i,k)*tl_dZ(i,k-1))-
900!^ & dZ(i,k)*(tl_dZ(i,k)+tl_dZ(i,k-1)))/
901!^ & (dZ1(i,k)+dZ1(i,k-1))
902!^ recursive
903 adfac=ad_dz(i,k)/(dz1(i,k)+dz1(i,k-1))
904 adfac1=adfac*2.0_r8
905 adfac2=adfac*dz(i,k)
906 ad_dz(i,k-1)=ad_dz(i,k-1)+dz1(i,k)*adfac1-adfac2
907 ad_dz(i,k )=dz1(i,k-1)*adfac1-adfac2
908
909 cff=2.0_r8*dr1(i,k)*dr1(i,k-1)
910 IF (cff.gt.eps) THEN
911!^ tl_dR(i,k)=(tl_cff-dR(i,k)*(tl_dR(i,k)+tl_dR(i,k-1)))/ &
912!^ & (dR1(i,k)+dR1(i,k-1))
913!^
914 adfac=ad_dr(i,k)/(dr1(i,k)+dr1(i,k-1))
915 adfac1=adfac*dr(i,k)
916 ad_dr(i,k-1)=ad_dr(i,k-1)-adfac1
917 ad_dr(i,k )=-adfac1
918 ad_cff=ad_cff+adfac
919 ELSE
920!^ tl_dR(i,k)=0.0_r8
921!^
922 ad_dr(i,k)=0.0_r8
923 END IF
924!^ tl_cff=2.0_r8*(tl_dR(i,k)*dR1(i,k-1)+ &
925!^ & dR1(i,k)*tl_dR(i,k-1))
926!^
927 adfac=2.0_r8*ad_cff
928 ad_dr(i,k-1)=ad_dr(i,k-1)+dr1(i,k )*adfac
929 ad_dr(i,k )=ad_dr(i,k )+dr1(i,k-1)*adfac
930 ad_cff=0.0_r8
931 END DO
932 END DO
933 DO i=istru-1,iend
934!^ tl_dZ(i,0)=tl_dZ(i,1)
935!^
936 ad_dz(i,1)=ad_dz(i,1)+ad_dz(i,0)
937 ad_dz(i,0)=0.0_r8
938!^ tl_dR(i,0)=tl_dR(i,1)
939!^
940 ad_dr(i,1)=ad_dr(i,1)+ad_dr(i,0)
941 ad_dr(i,0)=0.0_r8
942!^ tl_dZ(i,N(ng))=tl_dZ(i,N(ng)-1)
943!^
944 ad_dz(i,n(ng)-1)=ad_dz(i,n(ng)-1)+ad_dz(i,n(ng))
945 ad_dz(i,n(ng))=0.0_r8
946
947!^ tl_dR(i,N(ng))=tl_dR(i,N(ng)-1)
948!^
949 ad_dr(i,n(ng)-1)=ad_dr(i,n(ng)-1)+ad_dr(i,n(ng))
950 ad_dr(i,n(ng))=0.0_r8
951 END DO
952 DO k=n(ng)-1,1,-1
953 DO i=istru-1,iend
954!^ tl_dZ(i,k)=tl_z_r(i,j,k+1)-tl_z_r(i,j,k)
955!^
956 ad_z_r(i,j,k )=ad_z_r(i,j,k )-ad_dz(i,k)
957 ad_z_r(i,j,k+1)=ad_z_r(i,j,k+1)+ad_dz(i,k)
958 ad_dz(i,k)=0.0_r8
959
960!^ tl_dR(i,k)=tl_rho(i,j,k+1)-tl_rho(i,j,k)
961!^
962 ad_rho(i,j,k )=ad_rho(i,j,k )-ad_dr(i,k)
963 ad_rho(i,j,k+1)=ad_rho(i,j,k+1)+ad_dr(i,k)
964 ad_dr(i,k)=0.0_r8
965 END DO
966 END DO
967 END DO
968!
969 RETURN
970 END SUBROUTINE ad_prsgrd32_tile
971
972 END MODULE ad_prsgrd_mod
subroutine, public ad_prsgrd(ng, tile)
Definition ad_prsgrd31.h:37
subroutine ad_prsgrd32_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, umask, vmask, om_v, on_u, hz, ad_hz, z_r, ad_z_r, z_w, ad_z_w, rho, ad_rho, eq_tide, ad_eq_tide, pair, ad_ru, ad_rv)
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter iadm
Definition mod_param.F:665
real(dp) g
real(dp) rho0
integer, dimension(:), allocatable nrhs
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