ROMS
Loading...
Searching...
No Matches
ad_prsgrd31.h
Go to the documentation of this file.
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 routine evalutes the adjoint baroclinic hydrostatic pressure !
11! gradient term using the standard and weighted Jacobians scheme of !
12! Song and Wright (1997). Notice that horizontal gradients are !
13! computed before of the vertical integration. !
14! !
15! The pressure gradient terms (m4/s2) are loaded into right-hand- !
16! side arrays "ad_ru" and "ad_rv". !
17! !
18! Reference: !
19! !
20! Song, Y.T. and D.G. Wright, 1997: A general pressure gradient !
21! formutlation for numerical ocean models. Part I: Scheme !
22! design and diagnostic analysis. DRAFT. !
23! !
24! BASIC STATE variables needed: Hz, rho, z_r, z_w !
25! !
26!=======================================================================
27!
28 implicit none
29!
30 PRIVATE
31 PUBLIC :: ad_prsgrd
32!
33 CONTAINS
34!
35!***********************************************************************
36 SUBROUTINE ad_prsgrd (ng, tile)
37!***********************************************************************
38!
39 USE mod_param
40#ifdef DIAGNOSTICS
41!! USE mod_diags
42#endif
43#ifdef ATM_PRESS
44 USE mod_forces
45#endif
46 USE mod_grid
47 USE mod_ocean
48 USE mod_stepping
49!
50! Imported variable declarations.
51!
52 integer, intent(in) :: ng, tile
53!
54! Local variable declarations.
55!
56 character (len=*), parameter :: myfile = &
57 & __FILE__
58!
59#include "tile.h"
60!
61#ifdef PROFILE
62 CALL wclock_on (ng, iadm, 23, __line__, myfile)
63#endif
64 CALL ad_prsgrd31_tile (ng, tile, &
65 & lbi, ubi, lbj, ubj, &
66 & imins, imaxs, jmins, jmaxs, &
67 & nrhs(ng), &
68 & grid(ng) % om_v, &
69 & grid(ng) % on_u, &
70 & grid(ng) % Hz, &
71 & grid(ng) % ad_Hz, &
72 & grid(ng) % z_r, &
73 & grid(ng) % ad_z_r, &
74 & grid(ng) % z_w, &
75 & grid(ng) % ad_z_w, &
76 & ocean(ng) % rho, &
77 & ocean(ng) % ad_rho, &
78#ifdef TIDE_GENERATING_FORCES
79 & ocean(ng) % eq_tide, &
80 & ocean(ng) % ad_eq_tide, &
81#endif
82#ifdef ATM_PRESS
83 & forces(ng) % Pair, &
84#endif
85#ifdef DIAGNOSTICS_UV
86!! & DIAGS(ng) % DiaRU, &
87!! & DIAGS(ng) % DiaRV, &
88#endif
89 & ocean(ng) % ad_ru, &
90 & ocean(ng) % ad_rv)
91#ifdef PROFILE
92 CALL wclock_off (ng, iadm, 23, __line__, myfile)
93#endif
94!
95 RETURN
96 END SUBROUTINE ad_prsgrd
97!
98!***********************************************************************
99 SUBROUTINE ad_prsgrd31_tile (ng, tile, &
100 & LBi, UBi, LBj, UBj, &
101 & IminS, ImaxS, JminS, JmaxS, &
102 & nrhs, &
103 & om_v, on_u, &
104 & Hz, ad_Hz, &
105 & z_r, ad_z_r, &
106 & z_w, ad_z_w, &
107 & rho, ad_rho, &
108#ifdef TIDE_GENERATING_FORCES
109 & eq_tide, ad_eq_tide, &
110#endif
111#ifdef ATM_PRESS
112 & Pair, &
113#endif
114#ifdef DIAGNOSTICS_UV
115!! & DiaRU, DiaRV, &
116#endif
117 & ad_ru, ad_rv)
118!***********************************************************************
119!
120 USE mod_param
121 USE mod_scalars
122!
123! Imported variable declarations.
124!
125 integer, intent(in) :: ng, tile
126 integer, intent(in) :: LBi, UBi, LBj, UBj
127 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
128 integer, intent(in) :: nrhs
129
130#ifdef ASSUMED_SHAPE
131 real(r8), intent(in) :: om_v(LBi:,LBj:)
132 real(r8), intent(in) :: on_u(LBi:,LBj:)
133 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
134 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
135 real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
136 real(r8), intent(in) :: rho(LBi:,LBj:,:)
137# ifdef ATM_PRESS
138 real(r8), intent(in) :: Pair(LBi:,LBj:)
139# endif
140# ifdef TIDE_GENERATING_FORCES
141 real(r8), intent(in) :: eq_tide(LBi:,LBj:)
142 real(r8), intent(inout) :: ad_eq_tide(LBi:,LBj:)
143# endif
144# ifdef DIAGNOSTICS_UV
145!! real(r8), intent(inout) :: DiaRU(LBi:,LBj:,:,:,:)
146!! real(r8), intent(inout) :: DiaRV(LBi:,LBj:,:,:,:)
147# endif
148 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
149 real(r8), intent(inout) :: ad_z_r(LBi:,LBj:,:)
150 real(r8), intent(inout) :: ad_z_w(LBi:,LBj:,0:)
151 real(r8), intent(inout) :: ad_rho(LBi:,LBj:,:)
152 real(r8), intent(inout) :: ad_ru(LBi:,LBj:,0:,:)
153 real(r8), intent(inout) :: ad_rv(LBi:,LBj:,0:,:)
154#else
155 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
156 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
157 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
158 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
159 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
160 real(r8), intent(in) :: rho(LBi:UBi,LBj:UBj,N(ng))
161# ifdef ATM_PRESS
162 real(r8), intent(in) :: Pair(LBi:UBi,LBj:UBj)
163# endif
164# ifdef TIDE_GENERATING_FORCES
165 real(r8), intent(in) :: eq_tide(LBi:,LBj:)
166 real(r8), intent(out) :: ad_eq_tide(LBi:,LBj:)
167# endif
168# ifdef DIAGNOSTICS_UV
169!! real(r8), intent(inout) :: DiaRU(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
170!! real(r8), intent(inout) :: DiaRV(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
171# endif
172 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
173 real(r8), intent(inout) :: ad_z_r(LBi:UBi,LBj:UBj,N(ng))
174 real(r8), intent(inout) :: ad_z_w(LBi:UBi,LBj:UBj,0:N(ng))
175 real(r8), intent(inout) :: ad_rho(LBi:UBi,LBj:UBj,N(ng))
176 real(r8), intent(inout) :: ad_ru(LBi:UBi,LBj:UBj,0:N(ng),2)
177 real(r8), intent(inout) :: ad_rv(LBi:UBi,LBj:UBj,0:N(ng),2)
178#endif
179!
180! Local variable declarations.
181!
182 integer :: i, j, k, kk
183
184 real(r8) :: fac, fac1, fac2, fac3
185 real(r8) :: cff1, cff2, cff3, cff4
186 real(r8) :: adfac, adfac1, adfac2
187 real(r8) :: ad_cff1, ad_cff2, ad_cff3, ad_cff4
188#ifdef WJ_GRADP
189 real(r8) :: gamma, ad_gamma
190#endif
191
192 real(r8), dimension(IminS:ImaxS) :: phie
193 real(r8), dimension(IminS:ImaxS) :: phix
194
195 real(r8), dimension(IminS:ImaxS) :: ad_phie
196 real(r8), dimension(IminS:ImaxS) :: ad_phix
197
198#include "set_bounds.h"
199!
200!-----------------------------------------------------------------------
201! Initialize adjoint private variables.
202!-----------------------------------------------------------------------
203!
204 ad_cff1=0.0_r8
205 ad_cff2=0.0_r8
206 ad_cff3=0.0_r8
207 ad_cff4=0.0_r8
208
209#ifdef WJ_GRADP
210 ad_gamma=0.0_r8
211#endif
212 DO i=imins,imaxs
213 ad_phie(i)=0.0_r8
214 ad_phix(i)=0.0_r8
215 END DO
216!
217!-----------------------------------------------------------------------
218! Calculate adjoint pressure gradient in the ETA-direction (m4/s2).
219!-----------------------------------------------------------------------
220!
221#ifdef ATM_PRESS
222 fac=100.0_r8/rho0
223#endif
224 fac1=0.5_r8*g/rho0
225 fac2=1000.0_r8*g/rho0
226 fac3=0.25_r8*g/rho0
227
228 j_loop : DO j=jstr,jend
229 IF (j.ge.jstrv) THEN
230!
231! Compute appropriate BASIC STATE "phie". Notice that a reverse
232! vertical integration of "phie" is carried out over kk-index.
233!
234!^ DO k=N(ng)-1,1,-1
235!^
236 DO k=1,n(ng)-1
237 DO i=istr,iend
238#ifdef TIDE_GENERATING_FORCES
239 cff1=z_w(i,j ,n(ng))-eq_tide(i,j )-z_r(i,j ,n(ng))+ &
240 & z_w(i,j-1,n(ng))-eq_tide(i,j-1)-z_r(i,j-1,n(ng))
241#else
242 cff1=z_w(i,j ,n(ng))-z_r(i,j ,n(ng))+ &
243 & z_w(i,j-1,n(ng))-z_r(i,j-1,n(ng))
244#endif
245 phie(i)=fac1*(rho(i,j,n(ng))-rho(i,j-1,n(ng)))*cff1
246#ifdef ATM_PRESS
247 phie(i)=phie(i)+fac*(pair(i,j)-pair(i,j-1))
248#endif
249#ifdef RHO_SURF
250 phie(i)=phie(i)+ &
251 & (fac2+fac1*(rho(i,j,n(ng))+rho(i,j-1,n(ng))))* &
252 & (z_w(i,j,n(ng))-z_w(i,j-1,n(ng)))
253#endif
254 END DO
255 DO kk=n(ng)-1,k,-1
256 DO i=istr,iend
257#ifdef WJ_GRADP
258 cff1=1.0_r8/((z_r(i,j ,kk+1)-z_r(i,j ,kk))* &
259 & (z_r(i,j-1,kk+1)-z_r(i,j-1,kk)))
260 cff2=z_r(i,j ,kk )-z_r(i,j-1,kk )+ &
261 & z_r(i,j ,kk+1)-z_r(i,j-1,kk+1)
262 cff3=z_r(i,j ,kk+1)-z_r(i,j ,kk )- &
263 & z_r(i,j-1,kk+1)+z_r(i,j-1,kk )
264 gamma=0.125_r8*cff1*cff2*cff3
265
266 cff1=(1.0_r8+gamma)*(rho(i,j,kk+1)-rho(i,j-1,kk+1))+ &
267 & (1.0_r8-gamma)*(rho(i,j,kk )-rho(i,j-1,kk ))
268 cff2=rho(i,j,kk+1)+rho(i,j-1,kk+1)- &
269 & rho(i,j,kk )-rho(i,j-1,kk )
270 cff3=z_r(i,j,kk+1)+z_r(i,j-1,kk+1)- &
271 & z_r(i,j,kk )-z_r(i,j-1,kk )
272 cff4=(1.0_r8+gamma)*(z_r(i,j,kk+1)-z_r(i,j-1,kk+1))+ &
273 & (1.0_r8-gamma)*(z_r(i,j,kk )-z_r(i,j-1,kk ))
274 phie(i)=phie(i)+ &
275 & fac3*(cff1*cff3-cff2*cff4)
276#else
277 cff1=rho(i,j,kk+1)-rho(i,j-1,kk+1)+ &
278 & rho(i,j,kk )-rho(i,j-1,kk )
279 cff2=rho(i,j,kk+1)+rho(i,j-1,kk+1)- &
280 & rho(i,j,kk )-rho(i,j-1,kk )
281 cff3=z_r(i,j,kk+1)+z_r(i,j-1,kk+1)- &
282 & z_r(i,j,kk )-z_r(i,j-1,kk )
283 cff4=z_r(i,j,kk+1)-z_r(i,j-1,kk+1)+ &
284 & z_r(i,j,kk )-z_r(i,j-1,kk )
285 phie(i)=phie(i)+ &
286 & fac3*(cff1*cff3-cff2*cff4)
287#endif
288 END DO
289 END DO
290!
291! Compute interior adjoint baroclinic pressure gradient. Differentiate
292! and then vertically integrate.
293!
294 DO i=istr,iend
295#ifdef DIAGNOSTICS_UV
296!! DiaRV(i,j,k,nrhs,M3pgrd)=rv(i,j,k,nrhs)
297#endif
298!^ tl_rv(i,j,k,nrhs)=-0.5_r8*om_v(i,j)* &
299!^ & ((tl_Hz(i,j,k)+tl_Hz(i,j-1,k))* &
300!^ & phie(i)+ &
301!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
302!^ & tl_phie(i))
303!^
304 adfac=-0.5_r8*om_v(i,j)*ad_rv(i,j,k,nrhs)
305 adfac1=adfac*phie(i)
306 ad_phie(i)=ad_phie(i)+ &
307 & (hz(i,j,k)+hz(i,j-1,k))*adfac
308 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac1
309 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac1
310 ad_rv(i,j,k,nrhs)=0.0_r8
311#ifdef WJ_GRADP
312 cff1=1.0_r8/((z_r(i,j ,k+1)-z_r(i,j ,k))* &
313 & (z_r(i,j-1,k+1)-z_r(i,j-1,k)))
314 cff2=z_r(i,j ,k )-z_r(i,j-1,k )+ &
315 & z_r(i,j ,k+1)-z_r(i,j-1,k+1)
316 cff3=z_r(i,j ,k+1)-z_r(i,j ,k )- &
317 & z_r(i,j-1,k+1)+z_r(i,j-1,k )
318 gamma=0.125_r8*cff1*cff2*cff3
319
320 cff1=(1.0_r8+gamma)*(rho(i,j,k+1)-rho(i,j-1,k+1))+ &
321 & (1.0_r8-gamma)*(rho(i,j,k )-rho(i,j-1,k ))
322 cff2=rho(i,j,k+1)+rho(i,j-1,k+1)- &
323 & rho(i,j,k )-rho(i,j-1,k )
324 cff3=z_r(i,j,k+1)+z_r(i,j-1,k+1)- &
325 & z_r(i,j,k )-z_r(i,j-1,k )
326 cff4=(1.0_r8+gamma)*(z_r(i,j,k+1)-z_r(i,j-1,k+1))+ &
327 & (1.0_r8-gamma)*(z_r(i,j,k )-z_r(i,j-1,k ))
328!^ tl_phie(i)=tl_phie(i)+ &
329!^ & fac3*(tl_cff1*cff3+ &
330!^ & cff1*tl_cff3- &
331!^ & tl_cff2*cff4- &
332!^ & cff2*tl_cff4)
333!^
334 adfac=fac3*ad_phie(i)
335 ad_cff1=ad_cff1+cff3*adfac
336 ad_cff2=ad_cff2-cff4*adfac
337 ad_cff3=ad_cff3+cff1*adfac
338 ad_cff4=ad_cff4-cff2*adfac
339!^ tl_cff4=tl_gamma*(z_r(i,j,k+1)-z_r(i,j-1,k+1)- &
340!^ & z_r(i,j,k )+z_r(i,j-1,k ))+ &
341!^ & (1.0_r8+gamma)*(tl_z_r(i,j ,k+1)- &
342!^ & tl_z_r(i,j-1,k+1))+ &
343!^ & (1.0_r8-gamma)*(tl_z_r(i,j ,k )- &
344!^ & tl_z_r(i,j-1,k ))
345!^ tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i,j-1,k+1)- &
346!^ & tl_z_r(i,j,k )-tl_z_r(i,j-1,k )
347!^
348 adfac1=(1.0_r8+gamma)*ad_cff4
349 adfac2=(1.0_r8-gamma)*ad_cff4
350 ad_z_r(i,j-1,k )=ad_z_r(i,j-1,k )-adfac2-ad_cff3
351 ad_z_r(i,j ,k )=ad_z_r(i,j ,k )+adfac2-ad_cff3
352 ad_z_r(i,j-1,k+1)=ad_z_r(i,j-1,k+1)-adfac1+ad_cff3
353 ad_z_r(i,j ,k+1)=ad_z_r(i,j ,k+1)+adfac1+ad_cff3
354 ad_gamma=ad_gamma+ &
355 & (z_r(i,j,k+1)-z_r(i,j-1,k+1)- &
356 & z_r(i,j,k )+z_r(i,j-1,k ))*ad_cff4
357 ad_cff4=0.0_r8
358 ad_cff3=0.0_r8
359!^ tl_cff2=tl_rho(i,j,k+1)+tl_rho(i,j-1,k+1)- &
360!^ & tl_rho(i,j,k )-tl_rho(i,j-1,k )
361!^ tl_cff1=tl_gamma*(rho(i,j,k+1)-rho(i,j-1,k+1)- &
362!^ & rho(i,j,k )+rho(i,j-1,k ))+ &
363!^ & (1.0_r8+gamma)*(tl_rho(i,j ,k+1)- &
364!^ & tl_rho(i,j-1,k+1))+ &
365!^ & (1.0_r8-gamma)*(tl_rho(i,j ,k )- &
366!^ & tl_rho(i,j-1,k ))
367!^
368 adfac1=(1.0_r8+gamma)*ad_cff1
369 adfac2=(1.0_r8-gamma)*ad_cff1
370 ad_rho(i,j-1,k )=ad_rho(i,j-1,k )-adfac2-ad_cff2
371 ad_rho(i,j ,k )=ad_rho(i,j ,k )+adfac2-ad_cff2
372 ad_rho(i,j-1,k+1)=ad_rho(i,j-1,k+1)-adfac1+ad_cff2
373 ad_rho(i,j ,k+1)=ad_rho(i,j ,k+1)+adfac1+ad_cff2
374 ad_gamma=ad_gamma+ &
375 & (rho(i,j,k+1)-rho(i,j-1,k+1)- &
376 & rho(i,j,k )+rho(i,j-1,k ))*ad_cff1
377 ad_cff2=0.0_r8
378 ad_cff1=0.0_r8
379!
380 cff1=1.0_r8/((z_r(i,j ,k+1)-z_r(i,j ,k))* &
381 & (z_r(i,j-1,k+1)-z_r(i,j-1,k)))
382 cff2=z_r(i,j ,k )-z_r(i,j-1,k )+ &
383 & z_r(i,j ,k+1)-z_r(i,j-1,k+1)
384 cff3=z_r(i,j ,k+1)-z_r(i,j ,k )- &
385 & z_r(i,j-1,k+1)+z_r(i,j-1,k )
386
387!^ tl_gamma=0.125_r8*(tl_cff1*cff2*cff3+ &
388!^ & cff1*(tl_cff2*cff3+ &
389!^ & cff2*tl_cff3))
390!^
391 adfac=0.125_r8*ad_gamma
392 adfac1=adfac*cff1
393 ad_cff3=ad_cff3+cff2*adfac1
394 ad_cff2=ad_cff2+cff3*adfac1
395 ad_cff1=ad_cff1+cff2*cff3*adfac
396 ad_gamma=0.0_r8
397!^ tl_cff3=tl_z_r(i,j ,k+1)-tl_z_r(i,j ,k )- &
398!^ & tl_z_r(i,j-1,k+1)+tl_z_r(i,j-1,k )
399!^ tl_cff2=tl_z_r(i,j ,k )-tl_z_r(i,j-1,k )+ &
400!^ & tl_z_r(i,j ,k+1)-tl_z_r(i,j-1,k+1)
401!^
402 ad_z_r(i,j-1,k )=ad_z_r(i,j-1,k )-ad_cff2+ad_cff3
403 ad_z_r(i,j ,k )=ad_z_r(i,j ,k )+ad_cff2-ad_cff3
404 ad_z_r(i,j-1,k+1)=ad_z_r(i,j-1,k+1)-ad_cff2-ad_cff3
405 ad_z_r(i,j ,k+1)=ad_z_r(i,j ,k+1)+ad_cff2+ad_cff3
406 ad_cff3=0.0_r8
407 ad_cff2=0.0_r8
408!^ tl_cff1=-cff1*cff1*((tl_z_r(i,j ,k+1)-tl_z_r(i,j ,k))* &
409!^ & (z_r(i,j-1,k+1)-z_r(i,j-1,k))+ &
410!^ & (z_r(i,j ,k+1)-z_r(i,j ,k))* &
411!^ & (tl_z_r(i,j-1,k+1)-tl_z_r(i,j-1,k)))
412!^
413 adfac=-cff1*cff1*ad_cff1
414 adfac1=adfac*(z_r(i,j-1,k+1)-z_r(i,j-1,k))
415 adfac2=adfac*(z_r(i,j ,k+1)-z_r(i,j ,k))
416 ad_z_r(i,j-1,k )=ad_z_r(i,j-1,k )-adfac2
417 ad_z_r(i,j ,k )=ad_z_r(i,j ,k )-adfac1
418 ad_z_r(i,j-1,k+1)=ad_z_r(i,j-1,k+1)+adfac2
419 ad_z_r(i,j ,k+1)=ad_z_r(i,j ,k+1)+adfac1
420 ad_cff1=0.0_r8
421#else
422!
423 cff1=rho(i,j,k+1)-rho(i,j-1,k+1)+ &
424 & rho(i,j,k )-rho(i,j-1,k )
425 cff2=rho(i,j,k+1)+rho(i,j-1,k+1)- &
426 & rho(i,j,k )-rho(i,j-1,k )
427 cff3=z_r(i,j,k+1)+z_r(i,j-1,k+1)- &
428 & z_r(i,j,k )-z_r(i,j-1,k )
429 cff4=z_r(i,j,k+1)-z_r(i,j-1,k+1)+ &
430 & z_r(i,j,k )-z_r(i,j-1,k )
431
432!^ tl_phie(i)=tl_phie(i)+ &
433!^ & fac3*(tl_cff1*cff3+ &
434!^ & cff1*tl_cff3- &
435!^ & tl_cff2*cff4- &
436!^ & cff2*tl_cff4)
437!^
438 adfac=fac3*ad_phie(i)
439 ad_cff1=ad_cff1+cff3*adfac
440 ad_cff2=ad_cff2-cff4*adfac
441 ad_cff3=ad_cff3+cff1*adfac
442 ad_cff4=ad_cff4-cff2*adfac
443!^ tl_cff4=tl_z_r(i,j,k+1)-tl_z_r(i,j-1,k+1)+ &
444!^ & tl_z_r(i,j,k )-tl_z_r(i,j-1,k )
445!^ tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i,j-1,k+1)- &
446!^ & tl_z_r(i,j,k )-tl_z_r(i,j-1,k )
447!^
448 ad_z_r(i,j-1,k )=ad_z_r(i,j-1,k )-ad_cff3-ad_cff4
449 ad_z_r(i,j ,k )=ad_z_r(i,j ,k )-ad_cff3+ad_cff4
450 ad_z_r(i,j-1,k+1)=ad_z_r(i,j-1,k+1)+ad_cff3-ad_cff4
451 ad_z_r(i,j ,k+1)=ad_z_r(i,j ,k+1)+ad_cff3+ad_cff4
452 ad_cff4=0.0_r8
453 ad_cff3=0.0_r8
454!^ tl_cff2=tl_rho(i,j,k+1)+tl_rho(i,j-1,k+1)- &
455!^ & tl_rho(i,j,k )-tl_rho(i,j-1,k )
456!^ tl_cff1=tl_rho(i,j,k+1)-tl_rho(i,j-1,k+1)+ &
457!^ & tl_rho(i,j,k )-tl_rho(i,j-1,k )
458!^
459 ad_rho(i,j-1,k )=ad_rho(i,j-1,k )-ad_cff1-ad_cff2
460 ad_rho(i,j ,k )=ad_rho(i,j ,k )+ad_cff1-ad_cff2
461 ad_rho(i,j-1,k+1)=ad_rho(i,j-1,k+1)-ad_cff1+ad_cff2
462 ad_rho(i,j ,k+1)=ad_rho(i,j ,k+1)+ad_cff1+ad_cff2
463 ad_cff2=0.0_r8
464 ad_cff1=0.0_r8
465#endif
466 END DO
467 END DO
468!
469! Compute surface adjoint baroclinic pressure gradient.
470!
471 DO i=istr,iend
472#ifdef TIDE_GENERATING_FORCES
473 cff1=z_w(i,j ,n(ng))-eq_tide(i,j )-z_r(i,j ,n(ng))+ &
474 & z_w(i,j-1,n(ng))-eq_tide(i,j-1)-z_r(i,j-1,n(ng))
475#else
476 cff1=z_w(i,j ,n(ng))-z_r(i,j ,n(ng))+ &
477 & z_w(i,j-1,n(ng))-z_r(i,j-1,n(ng))
478#endif
479 phie(i)=fac1*(rho(i,j,n(ng))-rho(i,j-1,n(ng)))*cff1
480#ifdef ATM_PRESS
481 phie(i)=phie(i)+fac*(pair(i,j)-pair(i,j-1))
482#endif
483#ifdef RHO_SURF
484 phie(i)=phie(i)+ &
485 & (fac2+fac1*(rho(i,j,n(ng))+rho(i,j-1,n(ng))))* &
486 & (z_w(i,j,n(ng))-z_w(i,j-1,n(ng)))
487#endif
488# ifdef DIAGNOSTICS_UV
489!! DiaRV(i,j,N(ng),nrhs,M3pgrd)=rv(i,j,N(ng),nrhs)
490# endif
491!^ tl_rv(i,j,N(ng),nrhs)=-0.5_r8*om_v(i,j)* &
492!^ & ((tl_Hz(i,j ,N(ng))+ &
493!^ & tl_Hz(i,j-1,N(ng)))*phie(i)+ &
494!^ & (Hz(i,j ,N(ng))+ &
495!^ & Hz(i,j-1,N(ng)))*tl_phie(i))
496!^
497 adfac=-0.5_r8*om_v(i,j)*ad_rv(i,j,n(ng),nrhs)
498 adfac1=adfac*phie(i)
499 ad_phie(i)=ad_phie(i)+(hz(i,j ,n(ng))+ &
500 & hz(i,j-1,n(ng)))*adfac
501 ad_hz(i,j-1,n(ng))=ad_hz(i,j-1,n(ng))+adfac1
502 ad_hz(i,j ,n(ng))=ad_hz(i,j ,n(ng))+adfac1
503 ad_rv(i,j,n(ng),nrhs)=0.0
504#ifdef RHO_SURF
505!^ tl_phie(i)=tl_phie(i)+ &
506!^ & (fac1*(tl_rho(i,j,N(ng))+tl_rho(i,j-1,N(ng))))* &
507!^ & (z_w(i,j,N(ng))-z_w(i,j-1,N(ng)))+ &
508!^ & (fac2+fac1*(rho(i,j,N(ng))+rho(i,j-1,N(ng))))* &
509!^ & (tl_z_w(i,j,N(ng))-tl_z_w(i,j-1,N(ng)))
510!^
511 adfac1=fac1*(z_w(i,j,n(ng))-z_w(i,j-1,n(ng)))* &
512 & ad_phie(i)
513 adfac2=(fac2+fac1*(rho(i,j,n(ng))+rho(i,j-1,n(ng))))* &
514 & ad_phie(i)
515 ad_rho(i,j-1,n(ng))=ad_rho(i,j-1,n(ng))+adfac1
516 ad_rho(i,j ,n(ng))=ad_rho(i,j ,n(ng))+adfac1
517 ad_z_w(i,j-1,n(ng))=ad_z_w(i,j-1,n(ng))-adfac2
518 ad_z_w(i,j ,n(ng))=ad_z_w(i,j ,n(ng))+adfac2
519#endif
520!^ tl_phie(i)=fac1* &
521!^ & ((tl_rho(i,j,N(ng))-tl_rho(i,j-1,N(ng)))*cff1+ &
522!^ & (rho(i,j,N(ng))-rho(i,j-1,N(ng)))*tl_cff1)
523!^
524 adfac=fac1*ad_phie(i)
525 adfac1=adfac*cff1
526 ad_rho(i,j-1,n(ng))=ad_rho(i,j-1,n(ng))-adfac1
527 ad_rho(i,j ,n(ng))=ad_rho(i,j ,n(ng))+adfac1
528 ad_cff1=ad_cff1+ &
529 & (rho(i,j,n(ng))-rho(i,j-1,n(ng)))*adfac
530 ad_phie(i)=0.0_r8
531#ifdef TIDE_GENERATING_FORCES
532!^ tl_cff1=tl_z_w(i,j ,N(ng))-tl_eq_tide(i,j )- &
533!^ & tl_z_r(i,j ,N(ng))+ &
534!^ & tl_z_w(i,j-1,N(ng))-tl_eq_tide(i,j-1)- &
535!^ & tl_z_r(i,j-1,N(ng))
536!^
537 ad_eq_tide(i,j-1)=ad_eq_tide(i,j-1)-ad_cff1
538 ad_eq_tide(i,j )=ad_eq_tide(i,j )-ad_cff1
539 ad_z_r(i,j-1,n(ng))=ad_z_r(i,j-1,n(ng))-ad_cff1
540 ad_z_r(i,j ,n(ng))=ad_z_r(i,j ,n(ng))-ad_cff1
541 ad_z_w(i,j-1,n(ng))=ad_z_w(i,j-1,n(ng))+ad_cff1
542 ad_z_w(i,j ,n(ng))=ad_z_w(i,j ,n(ng))+ad_cff1
543 ad_cff1=0.0_r8
544#else
545!^ tl_cff1=tl_z_w(i,j ,N(ng))-tl_z_r(i,j ,N(ng))+ &
546!^ & tl_z_w(i,j-1,N(ng))-tl_z_r(i,j-1,N(ng))
547!^
548 ad_z_r(i,j-1,n(ng))=ad_z_r(i,j-1,n(ng))-ad_cff1
549 ad_z_r(i,j ,n(ng))=ad_z_r(i,j ,n(ng))-ad_cff1
550 ad_z_w(i,j-1,n(ng))=ad_z_w(i,j-1,n(ng))+ad_cff1
551 ad_z_w(i,j ,n(ng))=ad_z_w(i,j ,n(ng))+ad_cff1
552 ad_cff1=0.0_r8
553#endif
554 END DO
555 END IF
556!
557!-----------------------------------------------------------------------
558! Calculate adjoint pressure gradient in the XI-direction (m4/s2).
559!-----------------------------------------------------------------------
560!
561! Compute appropriate BASIC STATE "phix". Notice that a reverse
562! vertical integration of "phix" is carried out over kk-index.
563!
564!^ DO k=N(ng)-1,1,-1
565!^
566 DO k=1,n(ng)-1
567 DO i=istru,iend
568#ifdef TIDE_GENERATING_FORCES
569 cff1=z_w(i ,j,n(ng))-eq_tide(i ,j)-z_r(i ,j,n(ng))+ &
570 & z_w(i-1,j,n(ng))-eq_tide(i-1,j)-z_r(i-1,j,n(ng))
571#else
572 cff1=z_w(i ,j,n(ng))-z_r(i ,j,n(ng))+ &
573 & z_w(i-1,j,n(ng))-z_r(i-1,j,n(ng))
574#endif
575 phix(i)=fac1*(rho(i,j,n(ng))-rho(i-1,j,n(ng)))*cff1
576#ifdef ATM_PRESS
577 phix(i)=phix(i)+fac*(pair(i,j)-pair(i-1,j))
578#endif
579#ifdef RHO_SURF
580 phix(i)=phix(i)+ &
581 & (fac2+fac1*(rho(i,j,n(ng))+rho(i-1,j,n(ng))))* &
582 & (z_w(i,j,n(ng))-z_w(i-1,j,n(ng)))
583#endif
584 END DO
585 DO kk=n(ng)-1,k,-1
586 DO i=istru,iend
587#ifdef WJ_GRADP
588 cff1=1.0_r8/((z_r(i ,j,kk+1)-z_r(i ,j,kk))* &
589 & (z_r(i-1,j,kk+1)-z_r(i-1,j,kk)))
590 cff2=z_r(i ,j,kk )-z_r(i-1,j,kk )+ &
591 & z_r(i ,j,kk+1)-z_r(i-1,j,k+1)
592 cff3=z_r(i ,j,kk+1)-z_r(i ,j,kk )- &
593 & z_r(i-1,j,kk+1)+z_r(i-1,j,kk )
594 gamma=0.125_r8*cff1*cff2*cff3
595
596 cff1=(1.0_r8+gamma)*(rho(i,j,kk+1)-rho(i-1,j,kk+1))+ &
597 & (1.0_r8-gamma)*(rho(i,j,kk )-rho(i-1,j,kk ))
598 cff2=rho(i,j,kk+1)+rho(i-1,j,kk+1)- &
599 & rho(i,j,kk )-rho(i-1,j,kk )
600 cff3=z_r(i,j,kk+1)+z_r(i-1,j,kk+1)- &
601 & z_r(i,j,kk )-z_r(i-1,j,kk )
602 cff4=(1.0_r8+gamma)*(z_r(i,j,kk+1)-z_r(i-1,j,kk+1))+ &
603 & (1.0_r8-gamma)*(z_r(i,j,kk )-z_r(i-1,j,kk ))
604 phix(i)=phix(i)+ &
605 & fac3*(cff1*cff3-cff2*cff4)
606#else
607 cff1=rho(i,j,kk+1)-rho(i-1,j,kk+1)+ &
608 & rho(i,j,kk )-rho(i-1,j,kk )
609 cff2=rho(i,j,kk+1)+rho(i-1,j,kk+1)- &
610 & rho(i,j,kk )-rho(i-1,j,kk )
611 cff3=z_r(i,j,kk+1)+z_r(i-1,j,kk+1)- &
612 & z_r(i,j,kk )-z_r(i-1,j,kk )
613 cff4=z_r(i,j,kk+1)-z_r(i-1,j,kk+1)+ &
614 & z_r(i,j,kk )-z_r(i-1,j,kk )
615 phix(i)=phix(i)+ &
616 & fac3*(cff1*cff3-cff2*cff4)
617#endif
618 END DO
619 END DO
620!
621! Compute interior adjoint baroclinic pressure gradient. Differentiate
622! and then vertically integrate.
623!
624 DO i=istru,iend
625# ifdef DIAGNOSTICS_UV
626!! DiaRU(i,j,k,nrhs,M3pgrd)=ru(i,j,k,nrhs)
627# endif
628!^ tl_ru(i,j,k,nrhs)=-0.5_r8*on_u(i,j)* &
629!^ & ((tl_Hz(i,j,k)+tl_Hz(i-1,j,k))* &
630!^ & phix(i)+ &
631!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
632!^ & tl_phix(i))
633!^
634 adfac=-0.5_r8*on_u(i,j)*ad_ru(i,j,k,nrhs)
635 adfac1=adfac*phix(i)
636 ad_phix(i)=ad_phix(i)+ &
637 & (hz(i,j,k)+hz(i-1,j,k))*adfac
638 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac1
639 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac1
640 ad_ru(i,j,k,nrhs)=0.0
641#ifdef WJ_GRADP
642 cff1=1.0_r8/((z_r(i ,j,k+1)-z_r(i ,j,k))* &
643 & (z_r(i-1,j,k+1)-z_r(i-1,j,k)))
644 cff2=z_r(i ,j,k )-z_r(i-1,j,k )+ &
645 & z_r(i ,j,k+1)-z_r(i-1,j,k+1)
646 cff3=z_r(i ,j,k+1)-z_r(i ,j,k )- &
647 & z_r(i-1,j,k+1)+z_r(i-1,j,k )
648 gamma=0.125_r8*cff1*cff2*cff3
649
650 cff1=(1.0_r8+gamma)*(rho(i,j,k+1)-rho(i-1,j,k+1))+ &
651 & (1.0_r8-gamma)*(rho(i,j,k )-rho(i-1,j,k ))
652 cff2=rho(i,j,k+1)+rho(i-1,j,k+1)- &
653 & rho(i,j,k )-rho(i-1,j,k )
654 cff3=z_r(i,j,k+1)+z_r(i-1,j,k+1)- &
655 & z_r(i,j,k )-z_r(i-1,j,k )
656 cff4=(1.0_r8+gamma)*(z_r(i,j,k+1)-z_r(i-1,j,k+1))+ &
657 & (1.0_r8-gamma)*(z_r(i,j,k )-z_r(i-1,j,k ))
658!^ tl_phix(i)=tl_phix(i)+ &
659!^ & fac3*(tl_cff1*cff3+ &
660!^ & cff1*tl_cff3- &
661!^ & tl_cff2*cff4- &
662!^ & cff2*tl_cff4)
663!^
664 adfac=fac3*ad_phix(i)
665 ad_cff1=ad_cff1+cff3*adfac
666 ad_cff2=ad_cff2-cff4*adfac
667 ad_cff3=ad_cff3+cff1*adfac
668 ad_cff4=ad_cff4-cff2*adfac
669!^ tl_cff4=tl_gamma*(z_r(i,j,k+1)-z_r(i-1,j,k+1)- &
670!^ & z_r(i,j,k )+z_r(i-1,j,k ))+ &
671!^ & (1.0_r8+gamma)*(tl_z_r(i ,j,k+1)- &
672!^ & tl_z_r(i-1,j,k+1))+ &
673!^ & (1.0_r8-gamma)*(tl_z_r(i ,j,k )- &
674!^ & tl_z_r(i-1,j,k ))
675!^ tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i-1,j,k+1)- &
676!^ & tl_z_r(i,j,k )-tl_z_r(i-1,j,k )
677!^
678 adfac1=(1.0_r8+gamma)*ad_cff4
679 adfac2=(1.0_r8-gamma)*ad_cff4
680 ad_z_r(i-1,j,k )=ad_z_r(i-1,j,k )-adfac2-ad_cff3
681 ad_z_r(i ,j,k )=ad_z_r(i ,j,k )+adfac2-ad_cff3
682 ad_z_r(i-1,j,k+1)=ad_z_r(i-1,j,k+1)-adfac1+ad_cff3
683 ad_z_r(i ,j,k+1)=ad_z_r(i ,j,k+1)+adfac1+ad_cff3
684 ad_gamma=ad_gamma+ &
685 & (z_r(i,j,k+1)-z_r(i-1,j,k+1)- &
686 & z_r(i,j,k )+z_r(i-1,j,k ))*ad_cff4
687 ad_cff4=0.0_r8
688 ad_cff3=0.0_r8
689!^ & tl_cff2=tl_rho(i,j,k+1)+tl_rho(i-1,j,k+1)- &
690!^ & tl_rho(i,j,k )-tl_rho(i-1,j,k )
691!^ tl_cff1=tl_gamma*(rho(i,j,k+1)-rho(i-1,j,k+1)- &
692!^ & rho(i,j,k )+rho(i-1,j,k ))+ &
693!^ & (1.0_r8+gamma)*(tl_rho(i ,j,k+1)- &
694!^ & tl_rho(i-1,j,k+1))+ &
695!^ & (1.0_r8-gamma)*(tl_rho(i ,j,k )- &
696!^ & tl_rho(i-1,j,k ))
697!^
698 adfac1=(1.0_r8+gamma)*ad_cff1
699 adfac2=(1.0_r8-gamma)*ad_cff1
700 ad_rho(i-1,j,k )=ad_rho(i-1,j,k )-adfac2-ad_cff2
701 ad_rho(i ,j,k )=ad_rho(i ,j,k )+adfac2-ad_cff2
702 ad_rho(i-1,j,k+1)=ad_rho(i-1,j,k+1)-adfac1+ad_cff2
703 ad_rho(i ,j,k+1)=ad_rho(i ,j,k+1)+adfac1+ad_cff2
704 ad_gamma=ad_gamma+ &
705 & (rho(i,j,k+1)-rho(i-1,j,k+1)- &
706 & rho(i,j,k )+rho(i-1,j,k ))*ad_cff1
707 ad_cff2=0.0_r8
708 ad_cff1=0.0_r8
709!
710 cff1=1.0_r8/((z_r(i ,j,k+1)-z_r(i ,j,k))* &
711 & (z_r(i-1,j,k+1)-z_r(i-1,j,k)))
712 cff2=z_r(i ,j,k )-z_r(i-1,j,k )+ &
713 & z_r(i ,j,k+1)-z_r(i-1,j,k+1)
714 cff3=z_r(i ,j,k+1)-z_r(i ,j,k )- &
715 & z_r(i-1,j,k+1)+z_r(i-1,j,k )
716
717!^ tl_gamma=0.125_r8*(tl_cff1*cff2*cff3+ &
718!^ & cff1*(tl_cff2*cff3+ &
719!^ & cff2*tl_cff3))
720!^
721 adfac=0.125_r8*ad_gamma
722 adfac1=adfac*cff1
723 ad_cff3=ad_cff3+cff2*adfac1
724 ad_cff2=ad_cff2+cff3*adfac1
725 ad_cff1=ad_cff1+cff2*cff3*adfac
726 ad_gamma=0.0_r8
727!^ & tl_cff3=tl_z_r(i ,j,k+1)-tl_z_r(i ,j,k )- &
728!^ & tl_z_r(i-1,j,k+1)+tl_z_r(i-1,j,k )
729!^ tl_cff2=tl_z_r(i ,j,k )-tl_z_r(i-1,j,k )+ &
730!^ & tl_z_r(i ,j,k+1)-tl_z_r(i-1,j,k+1)
731!^
732 ad_z_r(i-1,j,k )=ad_z_r(i-1,j,k )-ad_cff2+ad_cff3
733 ad_z_r(i ,j,k )=ad_z_r(i ,j,k )+ad_cff2-ad_cff3
734 ad_z_r(i-1,j,k+1)=ad_z_r(i-1,j,k+1)-ad_cff2-ad_cff3
735 ad_z_r(i ,j,k+1)=ad_z_r(i ,j,k+1)+ad_cff2+ad_cff3
736 ad_cff3=0.0
737 ad_cff2=0.0
738!^ tl_cff1=-cff1*cff1*((tl_z_r(i ,j,k+1)-tl_z_r(i ,j,k))* &
739!^ & (z_r(i-1,j,k+1)-z_r(i-1,j,k))+ &
740!^ & (z_r(i ,j,k+1)-z_r(i ,j,k))* &
741!^ & (tl_z_r(i-1,j,k+1)-tl_z_r(i-1,j,k)))
742!^
743 adfac=-cff1*cff1*ad_cff1
744 adfac1=adfac*(z_r(i-1,j,k+1)-z_r(i-1,j,k))
745 adfac2=adfac*(z_r(i ,j,k+1)-z_r(i ,j,k))
746 ad_z_r(i-1,j,k )=ad_z_r(i-1,j,k )-adfac2
747 ad_z_r(i ,j,k )=ad_z_r(i ,j,k )-adfac1
748 ad_z_r(i-1,j,k+1)=ad_z_r(i-1,j,k+1)+adfac2
749 ad_z_r(i ,j,k+1)=ad_z_r(i ,j,k+1)+adfac1
750 ad_cff1=0.0_r8
751#else
752 cff1=rho(i,j,k+1)-rho(i-1,j,k+1)+ &
753 & rho(i,j,k )-rho(i-1,j,k )
754 cff2=rho(i,j,k+1)+rho(i-1,j,k+1)- &
755 & rho(i,j,k )-rho(i-1,j,k )
756 cff3=z_r(i,j,k+1)+z_r(i-1,j,k+1)- &
757 & z_r(i,j,k )-z_r(i-1,j,k )
758 cff4=z_r(i,j,k+1)-z_r(i-1,j,k+1)+ &
759 & z_r(i,j,k )-z_r(i-1,j,k )
760!^ tl_phix(i)=tl_phix(i)+ &
761!^ & fac3*(tl_cff1*cff3+ &
762!^ & cff1*tl_cff3- &
763!^ & tl_cff2*cff4- &
764!^ & cff2*tl_cff4)
765!^
766 adfac=fac3*ad_phix(i)
767 ad_cff1=ad_cff1+cff3*adfac
768 ad_cff2=ad_cff2-cff4*adfac
769 ad_cff3=ad_cff3+cff1*adfac
770 ad_cff4=ad_cff4-cff2*adfac
771!^ tl_cff4=tl_z_r(i,j,k+1)-tl_z_r(i-1,j,k+1)+ &
772!^ & tl_z_r(i,j,k )-tl_z_r(i-1,j,k )
773!^ tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i-1,j,k+1)- &
774!^ & tl_z_r(i,j,k )-tl_z_r(i-1,j,k )
775!^
776 ad_z_r(i-1,j,k )=ad_z_r(i-1,j,k )-ad_cff3-ad_cff4
777 ad_z_r(i ,j,k )=ad_z_r(i ,j,k )-ad_cff3+ad_cff4
778 ad_z_r(i-1,j,k+1)=ad_z_r(i-1,j,k+1)+ad_cff3-ad_cff4
779 ad_z_r(i ,j,k+1)=ad_z_r(i ,j,k+1)+ad_cff3+ad_cff4
780 ad_cff4=0.0_r8
781 ad_cff3=0.0_r8
782!^ tl_cff1=tl_rho(i,j,k+1)-tl_rho(i-1,j,k+1)+ &
783!^ & tl_rho(i,j,k )-tl_rho(i-1,j,k )
784!^ tl_cff2=tl_rho(i,j,k+1)+tl_rho(i-1,j,k+1)- &
785!^ & tl_rho(i,j,k )-tl_rho(i-1,j,k )
786!^
787 ad_rho(i-1,j,k )=ad_rho(i-1,j,k )-ad_cff2-ad_cff1
788 ad_rho(i ,j,k )=ad_rho(i ,j,k )-ad_cff2+ad_cff1
789 ad_rho(i-1,j,k+1)=ad_rho(i-1,j,k+1)+ad_cff2-ad_cff1
790 ad_rho(i ,j,k+1)=ad_rho(i ,j,k+1)+ad_cff2+ad_cff1
791 ad_cff2=0.0_r8
792 ad_cff1=0.0_r8
793#endif
794 END DO
795 END DO
796!
797! Compute surface adjoint baroclinic pressure gradient.
798!
799 DO i=istru,iend
800#ifdef TIDE_GENERATING_FORCES
801 cff1=z_w(i ,j,n(ng))-eq_tide(i ,j)-z_r(i ,j,n(ng))+ &
802 & z_w(i-1,j,n(ng))-eq_tide(i-1,j)-z_r(i-1,j,n(ng))
803#else
804 cff1=z_w(i ,j,n(ng))-z_r(i ,j,n(ng))+ &
805 & z_w(i-1,j,n(ng))-z_r(i-1,j,n(ng))
806#endif
807 phix(i)=fac1*(rho(i,j,n(ng))-rho(i-1,j,n(ng)))*cff1
808#ifdef ATM_PRESS
809 phix(i)=phix(i)+fac*(pair(i,j)-pair(i-1,j))
810#endif
811#ifdef RHO_SURF
812 phix(i)=phix(i)+ &
813 & (fac2+fac1*(rho(i,j,n(ng))+rho(i-1,j,n(ng))))* &
814 & (z_w(i,j,n(ng))-z_w(i-1,j,n(ng)))
815#endif
816#ifdef DIAGNOSTICS_UV
817!! DiaRU(i,j,N(ng),nrhs,M3pgrd)=ru(i,j,N(ng),nrhs)
818#endif
819!^ tl_ru(i,j,N(ng),nrhs)=-0.5_r8*on_u(i,j)* &
820!^ & ((tl_Hz(i ,j,N(ng))+ &
821!^ & tl_Hz(i-1,j,N(ng)))*phix(i)+ &
822!^ & (Hz(i ,j,N(ng))+ &
823!^ & Hz(i-1,j,N(ng)))*tl_phix(i))
824!^
825 adfac=-0.5_r8*on_u(i,j)*ad_ru(i,j,n(ng),nrhs)
826 adfac1=adfac*phix(i)
827 ad_phix(i)=ad_phix(i)+(hz(i ,j,n(ng))+ &
828 & hz(i-1,j,n(ng)))*adfac
829 ad_hz(i-1,j,n(ng))=ad_hz(i-1,j,n(ng))+adfac1
830 ad_hz(i ,j,n(ng))=ad_hz(i ,j,n(ng))+adfac1
831 ad_ru(i,j,n(ng),nrhs)=0.0_r8
832#ifdef RHO_SURF
833!^ tl_phix(i)=tl_phix(i)+ &
834!^ & (fac1*(tl_rho(i,j,N(ng))+tl_rho(i-1,j,N(ng))))* &
835!^ & (z_w(i,j,N(ng))-z_w(i-1,j,N(ng)))+ &
836!^ & (fac2+fac1*(rho(i,j,N(ng))+rho(i-1,j,N(ng))))*
837!^ & (tl_z_w(i,j,N(ng))-tl_z_w(i-1,j,N(ng)))
838!^
839 adfac1=fac1*(z_w(i,j,n(ng))-z_w(i-1,j,n(ng)))* &
840 & ad_phix(i)
841 adfac2=(fac2+fac1*(rho(i,j,n(ng))+rho(i-1,j,n(ng))))* &
842 & ad_phix(i)
843 ad_rho(i-1,j,n(ng))=ad_rho(i-1,j,n(ng))+adfac1
844 ad_rho(i ,j,n(ng))=ad_rho(i ,j,n(ng))+adfac1
845 ad_z_w(i-1,j,n(ng))=ad_z_w(i-1,j,n(ng))-adfac2
846 ad_z_w(i ,j,n(ng))=ad_z_w(i ,j,n(ng))+adfac2
847#endif
848!^ tl_phix(i)=fac1* &
849!^ & ((tl_rho(i,j,N(ng))-tl_rho(i-1,j,N(ng)))*cff1+ &
850!^ & (rho(i,j,N(ng))-rho(i-1,j,N(ng)))*tl_cff1)
851!^
852 adfac=fac1*ad_phix(i)
853 adfac1=adfac*cff1
854 ad_rho(i-1,j,n(ng))=ad_rho(i-1,j,n(ng))-adfac1
855 ad_rho(i ,j,n(ng))=ad_rho(i ,j,n(ng))+adfac1
856 ad_cff1=ad_cff1+ &
857 & (rho(i,j,n(ng))-rho(i-1,j,n(ng)))*adfac
858 ad_phix(i)=0.0
859#ifdef TIDE_GENERATING_FORCES
860!^ tl_cff1=tl_z_w(i ,j,N(ng))-tl_eq_tide(i ,j)- &
861!^ & tl_z_r(i ,j,N(ng))+ &
862!^ & tl_z_w(i-1,j,N(ng))-tl_eq_tide(i-1,j)- &
863!^ & tl_z_r(i-1,j,N(ng))
864!^
865 ad_eq_tide(i-1,j)=ad_eq_tide(i-1,j)-ad_cff1
866 ad_eq_tide(i ,j)=ad_eq_tide(i ,j)-ad_cff1
867 ad_z_r(i-1,j,n(ng))=ad_z_r(i-1,j,n(ng))-ad_cff1
868 ad_z_r(i ,j,n(ng))=ad_z_r(i ,j,n(ng))-ad_cff1
869 ad_z_w(i-1,j,n(ng))=ad_z_w(i-1,j,n(ng))+ad_cff1
870 ad_z_w(i ,j,n(ng))=ad_z_w(i ,j,n(ng))+ad_cff1
871 ad_cff1=0.0_r8
872#else
873!^ tl_cff1=tl_z_w(i ,j,N(ng))-tl_z_r(i ,j,N(ng))+ &
874!^ & tl_z_w(i-1,j,N(ng))-tl_z_r(i-1,j,N(ng))
875!^
876 ad_z_r(i-1,j,n(ng))=ad_z_r(i-1,j,n(ng))-ad_cff1
877 ad_z_r(i ,j,n(ng))=ad_z_r(i ,j,n(ng))-ad_cff1
878 ad_z_w(i-1,j,n(ng))=ad_z_w(i-1,j,n(ng))+ad_cff1
879 ad_z_w(i ,j,n(ng))=ad_z_w(i ,j,n(ng))+ad_cff1
880 ad_cff1=0.0_r8
881#endif
882 END DO
883 END DO j_loop
884!
885 RETURN
886 END SUBROUTINE ad_prsgrd31_tile
887
888 END MODULE ad_prsgrd_mod
subroutine, public ad_prsgrd(ng, tile)
Definition ad_prsgrd31.h:37
subroutine ad_prsgrd31_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, 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