ROMS
Loading...
Searching...
No Matches
ad_prsgrd40.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 subroutine evaluates the adjoint baroclinic, hydrostatic !
11! pressure gradient term using the finite-volume pressure Jacobian !
12! scheme of Lin (1997). !
13! !
14! The pressure gradient terms (m4/s2) are loaded into right-hand- !
15! side arrays "tl_ru" and "tl_rv". !
16! !
17! Reference: !
18! !
19! Lin, Shian-Jiann, 1997: A finite volume integration method !
20! for computing pressure gradient force in general vertical !
21! coordinates, Q. J. R. Meteorol. Soc., 123, 1749-1762. !
22! !
23!=======================================================================
24!
25 implicit none
26!
27 PRIVATE
28 PUBLIC :: ad_prsgrd
29!
30 CONTAINS
31!
32!***********************************************************************
33 SUBROUTINE ad_prsgrd (ng, tile)
34!***********************************************************************
35!
36 USE mod_param
37#ifdef DIAGNOSTICS
38!! USE mod_diags
39#endif
40#ifdef ATM_PRESS
41 USE mod_forces
42#endif
43 USE mod_grid
44 USE mod_ocean
45 USE mod_stepping
46!
47! Imported variable declarations.
48!
49 integer, intent(in) :: ng, tile
50!
51! Local variable declarations.
52!
53 character (len=*), parameter :: MyFile = &
54 & __FILE__
55!
56#include "tile.h"
57!
58#ifdef PROFILE
59 CALL wclock_on (ng, iadm, 23, __line__, myfile)
60#endif
61 CALL ad_prsgrd40_tile (ng, tile, &
62 & lbi, ubi, lbj, ubj, &
63 & imins, imaxs, jmins, jmaxs, &
64 & nrhs(ng), &
65 & grid(ng) % om_v, &
66 & grid(ng) % on_u, &
67 & grid(ng) % Hz, &
68 & grid(ng) % ad_Hz, &
69 & grid(ng) % z_w, &
70 & grid(ng) % ad_z_w, &
71 & ocean(ng) % rho, &
72 & ocean(ng) % ad_rho, &
73#ifdef TIDE_GENERATING_FORCES
74 & ocean(ng) % eq_tide, &
75 & ocean(ng) % ad_eq_tide, &
76#endif
77#ifdef ATM_PRESS
78 & forces(ng) % Pair, &
79#endif
80#ifdef DIAGNOSTICS_UV
81!! & DIAGS(ng) % DiaRU, &
82!! & DIAGS(ng) % DiaRV, &
83#endif
84 & ocean(ng) % ad_ru, &
85 & ocean(ng) % ad_rv)
86#ifdef PROFILE
87 CALL wclock_off (ng, iadm, 23, __line__, myfile)
88#endif
89!
90 RETURN
91 END SUBROUTINE ad_prsgrd
92!
93!***********************************************************************
94 SUBROUTINE ad_prsgrd40_tile (ng, tile, &
95 & LBi, UBi, LBj, UBj, &
96 & IminS, ImaxS, JminS, JmaxS, &
97 & nrhs, &
98 & om_v, on_u, &
99 & Hz, ad_Hz, &
100 & z_w, ad_z_w, &
101 & rho, ad_rho, &
102#ifdef TIDE_GENERATING_FORCES
103 & eq_tide, ad_eq_tide, &
104#endif
105#ifdef ATM_PRESS
106 & Pair, &
107#endif
108#ifdef DIAGNOSTICS_UV
109!! & DiaRU, DiaRV, &
110#endif
111 & ad_ru, ad_rv)
112!***********************************************************************
113!
114 USE mod_param
115 USE mod_scalars
116!
117! Imported variable declarations.
118!
119 integer, intent(in) :: ng, tile
120 integer, intent(in) :: LBi, UBi, LBj, UBj
121 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
122 integer, intent(in) :: nrhs
123
124#ifdef ASSUMED_SHAPE
125 real(r8), intent(in) :: om_v(LBi:,LBj:)
126 real(r8), intent(in) :: on_u(LBi:,LBj:)
127 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
128 real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
129 real(r8), intent(in) :: rho(LBi:,LBj:,:)
130# ifdef ATM_PRESS
131 real(r8), intent(in) :: Pair(LBi:,LBj:)
132# endif
133# ifdef TIDE_GENERATING_FORCES
134 real(r8), intent(in) :: eq_tide(LBi:,LBj:)
135 real(r8), intent(inout) :: ad_eq_tide(LBi:,LBj:)
136# endif
137# ifdef DIAGNOSTICS_UV
138!! real(r8), intent(inout) :: DiaRU(LBi:,LBj:,:,:,:)
139!! real(r8), intent(inout) :: DiaRV(LBi:,LBj:,:,:,:)
140# endif
141 real(r8), intent(inout) :: ad_Hz(LBi:,LBj:,:)
142 real(r8), intent(inout) :: ad_z_w(LBi:,LBj:,0:)
143 real(r8), intent(inout) :: ad_rho(LBi:,LBj:,:)
144 real(r8), intent(inout) :: ad_ru(LBi:,LBj:,0:,:)
145 real(r8), intent(inout) :: ad_rv(LBi:,LBj:,0:,:)
146#else
147 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
148 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
149 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
150 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
151 real(r8), intent(in) :: rho(LBi:UBi,LBj:UBj,N(ng))
152# ifdef ATM_PRESS
153 real(r8), intent(in) :: Pair(LBi:UBi,LBj:UBj)
154# endif
155# ifdef TIDE_GENERATING_FORCES
156 real(r8), intent(in) :: eq_tide(LBi:UBi,LBj:UBj)
157 real(r8), intent(inout) :: ad_eq_tide(LBi:UBi,LBj:UBj)
158# endif
159# ifdef DIAGNOSTICS_UV
160!! real(r8), intent(inout) :: DiaRU(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
161!! real(r8), intent(inout) :: DiaRV(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
162# endif
163 real(r8), intent(inout) :: ad_Hz(LBi:UBi,LBj:UBj,N(ng))
164 real(r8), intent(inout) :: ad_z_w(LBi:UBi,LBj:UBj,0:N(ng))
165 real(r8), intent(inout) :: ad_rho(LBi:UBi,LBj:UBj,N(ng))
166 real(r8), intent(inout) :: ad_ru(LBi:UBi,LBj:UBj,0:N(ng),2)
167 real(r8), intent(inout) :: ad_rv(LBi:UBi,LBj:UBj,0:N(ng),2)
168#endif
169!
170! Local variable declarations.
171!
172 integer :: i, j, k
173
174 real(r8) :: cff, cff1, dh
175 real(r8) :: ad_dh
176#ifdef ATM_PRESS
177 real(r8) :: OneAtm, fac
178#endif
179 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: ad_FC
180
181 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: ad_FX
182
183 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,0:N(ng)) :: P
184 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,0:N(ng)) :: ad_P
185
186#include "set_bounds.h"
187!
188!-----------------------------------------------------------------------
189! Initialize adjoint private variables.
190!-----------------------------------------------------------------------
191!
192 ad_dh=0.0_r8
193 DO j=jmins,jmaxs
194 DO i=imins,imaxs
195 ad_fx=(i,j)=0.0_r8
196 END DO
197 DO k=0,n(ng)
198 DO i=imins,imaxs
199 ad_p(i,j,k)=0.0_r8
200 END DO
201 END DO
202 END DO
203 DO k=0,n(ng)
204 DO i=imins,imaxs
205 ad_fc=(i,k)=0.0_r8
206 END DO
207 END DO
208!
209!-----------------------------------------------------------------------
210! Finite Volume pressure gradient algorithm (Lin, 1997).
211!-----------------------------------------------------------------------
212!
213! Compute BASIC STATE dynamic pressure, P.
214!
215#ifdef ATM_PRESS
216 oneatm=1013.25_r8 ! 1 atm = 1013.25 mb
217 fac=100.0_r8/g
218#endif
219 DO j=jstrv-1,jend
220 DO i=istru-1,iend
221 p(i,j,n(ng))=0.0_r8
222#ifdef ATM_PRESS
223 p(i,j,n(ng))=p(i,j,n(ng))+fac*(pair(i,j)-oneatm)
224#endif
225#ifdef TIDE_GENERATING_FORCES
226 p(i,j,n(ng))=p(i,j,n(ng))-g*eq_tide(i,j)
227#endif
228 END DO
229 DO k=n(ng),1,-1
230 DO i=istru-1,iend
231 p(i,j,k-1)=p(i,j,k)+hz(i,j,k)*rho(i,j,k)
232 END DO
233 END DO
234 END DO
235!
236! Calculate adjoint of pressure gradient in the ETA-direction (m4/s2).
237!
238 j_loop: DO j=jend,jstrv-1,-1
239 IF (j.ge.jstrv) THEN
240 cff=0.5_r8*g
241 cff1=g/rho0
242 DO k=1,n(ng)
243 DO i=istr,iend
244 dh=z_w(i,j,k-1)-z_w(i,j-1,k-1)
245#ifdef DIAGNOSTICS_UV
246!! DiaRV(i,j,k,nrhs,M3pgrd)=rv(i,j,k,nrhs)
247#endif
248!^ tl_rv(i,j,k,nrhs)=(cff*((tl_Hz(i,j-1,k)+ &
249!^ & tl_Hz(i,j ,k))* &
250!^ & (z_w(i,j-1,N(ng))- &
251!^ & z_w(i,j ,N(ng)))+ &
252!^ & (Hz(i,j-1,k)+ &
253!^ & Hz(i,j ,k))* &
254!^ & (tl_z_w(i,j-1,N(ng))- &
255!^ & tl_z_w(i,j ,N(ng))))+ &
256!^ & cff1*(tl_FX(i,j-1,k)- &
257!^ & tl_FX(i,j ,k)+ &
258!^ & tl_FC(i,k )- &
259!^ & tl_FC(i,k-1)))*om_v(i,j)
260!^
261 adfac=om_v(i,j)*ad_rv(i,j,k,nrhs)
262 adfac1=adfac*cff
263 adfac2=adfac1*(z_w(i,j-1,n(ng))- &
264 & z_w(i,j ,n(ng)))
265 adfac3=adfac1*(hz(i,j-1,k)+ &
266 & hz(i,j ,k))
267 adfac4=adfac*cff1
268 ad_hz(i,j-1,k)=ad_hz(i,j-1,k)+adfac2
269 ad_hz(i,j ,k)=ad_hz(i,j ,k)+adfac2
270 ad_z_w(i,j-1,n(ng))=ad_z_w(i,j-1,n(ng))+adfac3
271 ad_z_w(i,j ,n(ng))=ad_z_w(i,j ,n(ng))-adfac3
272 ad_fx(i,j-1,k)=ad_fx(i,j-1,k)+adfac4
273 ad_fx(i,j ,k)=ad_fx(i,j ,k)-adfac4
274 ad_fc(i,k-1)=ad_fc(i,k-1)-adfac4
275 ad_fc(i,k )=ad_fc(i,k )+adfac4
276 ad_rv(i,j,k,nrhs)=0.0_r8
277!^ tl_FC(i,k-1)=0.5_r8* &
278!^ & (tl_dh*(P(i,j,k-1)+P(i,j-1,k-1))+ &
279!^ & dh*(tl_P(i,j,k-1)+tl_P(i,j-1,k-1))
280!^
281 adfac=0.5_r8*ad_fc(i,k-1)
282 adfac1=adfac*dh
283 ad_dh=ad_dh+(p(i,j,k-1)+p(i,j-1,k-1))*adfac
284 ad_p(i,j-1,k-1)=ad_p(i,j-1,k-1)+adfac1
285 ad_p(i,j ,k-1)=ad_p(i,j ,k-1)+adfac1
286 ad_fc(i,k-1)=0.0_r8
287!^ tl_dh=tl_z_w(i,j,k-1)-tl_z_w(i,j-1,k-1)
288!^
289 ad_z_w(i,j-1,k-1)=ad_z_w(i,j-1,k-1)-ad_dh
290 ad_z_w(i,j ,k-1)=ad_z_w(i,j ,k-1)+ad_dh
291 ad_dh=0.0_r8
292 END DO
293 END DO
294 DO i=istr,iend
295!^ tl_FC(i,N(ng))=0.0_r8
296!^
297 ad_fc(i,n(ng))=0.0_r8
298 END DO
299 END IF
300!
301! Calculate the adjoint of pressure gradient in the XI-direction (m4/s2).
302!
303 IF (j.ge.jstr) THEN
304 cff=0.5_r8*g
305 cff1=g/rho0
306 DO k=1,n(ng)
307 DO i=iend,istru,-1
308 dh=z_w(i,j,k-1)-z_w(i-1,j,k-1)
309#ifdef DIAGNOSTICS_UV
310!! DiaRU(i,j,k,nrhs,M3pgrd)=ru(i,j,k,nrhs)
311#endif
312!^ tl_ru(i,j,k,nrhs)=(cff*((tl_Hz(i-1,j,k)+ &
313!^ & tl_Hz(i ,j,k))* &
314!^ & (z_w(i-1,j,N(ng))- &
315!^ & z_w(i ,j,N(ng)))+ &
316!^ & (Hz(i-1,j,k)+ &
317!^ & Hz(i ,j,k))* &
318!^ & (tl_z_w(i-1,j,N(ng))- &
319!^ & tl_z_w(i ,j,N(ng))))+ &
320!^ & cff1*(tl_FX(i-1,j,k)- &
321!^ & tl_FX(i ,j,k)+ &
322!^ & tl_FC(i,k )- &
323!^ & tl_FC(i,k-1)))*on_u(i,j)
324!^
325 adfac=on_u(i,j)*tl_ru(i,j,k,nrhs)
326 adfac1=adfac*cff
327 adfac2=adfac1*(z_w(i-1,j,n(ng))- &
328 & z_w(i ,j,n(ng)))
329 adfac3=adfac1*(hz(i-1,j,k)+ &
330 & hz(i ,j,k))
331 adfac4=adfac*cff1
332 ad_hz(i-1,j,k)=ad_hz(i-1,j,k)+adfac2
333 ad_hz(i ,j,k)=ad_hz(i ,j,k)+adfac2
334 ad_z_w(i-1,j,n(ng))=ad_z_w(i-1,j,n(ng))+adfac3
335 ad_z_w(i ,j,n(ng))=ad_z_w(i ,j,n(ng))-adfac3
336 ad_fx(i-1,j,k)=ad_fx(i-1,j,k)+adfac4
337 ad_fx(i ,j,k)=ad_fx(i ,j,k)-adfac4
338 ad_fc(i,k )=ad_fc(i,k )+adfac4
339 ad_fc(i,k-1)=ad_fc(i,k-1)-adfac4
340 ad_ru(i,j,k,nrhs)=0.0_r8
341!^ tl_FC(i,k-1)=0.5_r8* &
342!^ & (tl_dh*(P(i,j,k-1)+P(i-1,j,k-1))+ &
343!^ & dh*(tl_P(i,j,k-1)+tl_P(i-1,j,k-1)))
344!^
345 adfac=0.5_r8*tl_fc(i,k-1)
346 adfac1=adfac*dh
347 ad_dh=ad_dh+(p(i,j,k-1)+p(i-1,j,k-1))*adfac
348 ad_p(i-1,j,k-1)=ad_p(i-1,j,k-1)+adfac1
349 ad_p(i ,j,k-1)=ad_p(i ,j,k-1)+adfac1
350 ad_fc(i,k-1)=0.0_r8
351!^ tl_dh=tl_z_w(i,j,k-1)-tl_z_w(i-1,j,k-1)
352!^
353 ad_z_w(i-1,j,k-1)=ad_z_w(i-1,j,k-1)-ad_dh
354 ad_z_w(i ,j,k-1)=ad_z_w(i ,j,k-1)+ad_dh
355 ad_dh=0.0_r8
356 END DO
357 END DO
358 DO i=istru,iend
359!^ tl_FC(i,N(ng))=0.0_r8
360 ad_fc(i,n(ng))=0.0_r8
361 END DO
362 END IF
363!
364! Compute adjoint pressure and its vertical integral.
365!
366 DO k=1,n(ng)
367 DO i=istru-1,iend
368!^ tl_FX(i,j,k)=0.5_r8* &
369!^ & (tl_Hz(i,j,k)*(P(i,j,k)+P(i,j,k-1))+ &
370!^ & Hz(i,j,k)*(tl_P(i,j,k)+tl_P(i,j,k-1)))
371!^
372 adfac=0.5_r8*ad_fx(i,j,k)
373 adfac1=adfac*hz(i,j,k)
374 ad_hz(i,j,k)=ad_hz(i,j,k)+(p(i,j,k)+p(i,j,k-1))*adfac
375 ad_p(i,j,k-1)=ad_p(i,j,k-1)+adfac1
376 ad_p(i,j,k )=ad_p(i,j,k )+adfac1
377 ad_fx(i,j,k)=0.0_r8
378!^ tl_P(i,j,k-1)=tl_P(i,j,k)+ &
379!^ & tl_Hz(i,j,k)*rho(i,j,k)+ &
380!^ & Hz(i,j,k)*tl_rho(i,j,k)
381!^
382 ad_p(i,j,k)=ad_p(i,j,k)+ad_p(i,j,k-1)
383 ad_hz(i,j,k)=ad_hz(i,j,k)+rho(i,j,k)*ad_p(i,j,k-1)
384 ad_rho(i,j,k)=ad_rho(i,j,k)+hz(i,j,k)*ad_p(i,j,k-1)
385 END DO
386 END DO
387 DO i=istru-1,iend
388#ifdef TIDE_GENERATING_FORCES
389!^ tl_P(i,j,N(ng))=tl_P(i,j,N(ng))-g*tl_eq_tide(i,j)
390!^
391 ad_eq_tide(i,j)=ad_eq_tide(i,j)-g*ad_p(i,j,n(ng))
392#endif
393!^ tl_P(i,j,N(ng))=0.0_r8
394!^
395 ad_p(i,j,n(ng))=0.0_r8
396 END DO
397 END DO j_loop
398!
399 RETURN
400 END SUBROUTINE ad_prsgrd40_tile
401
402 END MODULE ad_prsgrd_mod
subroutine, public ad_prsgrd(ng, tile)
Definition ad_prsgrd31.h:37
subroutine ad_prsgrd40_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, om_v, on_u, hz, ad_hz, 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