ROMS
Loading...
Searching...
No Matches
rp_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 representers tangent linear baroclinic !
11! hydrostatic pressure gradient term using the standard and !
12! weighted Jacobians scheme of Song and Wright (1997). Note that !
13! horizontal gradients are computed before of the vertical !
14! integration. !
15! !
16! The pressure gradient terms (m4/s2) are loaded into right-hand- !
17! side arrays "tl_ru" and "tl_rv". !
18! !
19! Reference: !
20! !
21! Song, Y.T., 1998: A general pressure gradient formulation for !
22! numerical ocean models. Part I: Scheme design and diagnostic !
23! analysis, Monthly Weather Rev., 126, 3213-3230. !
24! !
25!=======================================================================
26!
27 implicit none
28!
29 PRIVATE
30 PUBLIC :: rp_prsgrd
31!
32 CONTAINS
33!
34!***********************************************************************
35 SUBROUTINE rp_prsgrd (ng, tile)
36!***********************************************************************
37!
38 USE mod_param
39#ifdef DIAGNOSTICS
40!! USE mod_diags
41#endif
42#ifdef ATM_PRESS
43 USE mod_forces
44#endif
45 USE mod_grid
46 USE mod_ocean
47 USE mod_stepping
48!
49! Imported variable declarations.
50!
51 integer, intent(in) :: ng, tile
52!
53! Local variable declarations.
54!
55 character (len=*), parameter :: myfile = &
56 & __FILE__
57!
58#include "tile.h"
59!
60#ifdef PROFILE
61 CALL wclock_on (ng, irpm, 23, __line__, myfile)
62#endif
63 CALL rp_prsgrd31_tile (ng, tile, &
64 & lbi, ubi, lbj, ubj, &
65 & imins, imaxs, jmins, jmaxs, &
66 & nrhs(ng), &
67 & grid(ng) % om_v, &
68 & grid(ng) % on_u, &
69 & grid(ng) % Hz, &
70 & grid(ng) % tl_Hz, &
71 & grid(ng) % z_r, &
72 & grid(ng) % tl_z_r, &
73 & grid(ng) % z_w, &
74 & grid(ng) % tl_z_w, &
75 & ocean(ng) % rho, &
76 & ocean(ng) % tl_rho, &
77#ifdef TIDE_GENERATING_FORCES
78 & ocean(ng) % eq_tide, &
79 & ocean(ng) % tl_eq_tide, &
80#endif
81#ifdef ATM_PRESS
82 & forces(ng) % Pair, &
83#endif
84#ifdef DIAGNOSTICS_UV
85!! & DIAGS(ng) % DiaRU, &
86!! & DIAGS(ng) % DiaRV, &
87#endif
88 & ocean(ng) % tl_ru, &
89 & ocean(ng) % tl_rv)
90#ifdef PROFILE
91 CALL wclock_off (ng, irpm, 23, __line__, myfile)
92#endif
93!
94 RETURN
95 END SUBROUTINE rp_prsgrd
96!
97!***********************************************************************
98 SUBROUTINE rp_prsgrd31_tile (ng, tile, &
99 & LBi, UBi, LBj, UBj, &
100 & IminS, ImaxS, JminS, JmaxS, &
101 & nrhs, &
102 & om_v, on_u, &
103 & Hz, tl_Hz, &
104 & z_r, tl_z_r, &
105 & z_w, tl_z_w, &
106 & rho, tl_rho, &
107#ifdef TIDE_GENERATING_FORCES
108 & eq_tide, tl_eq_tide, &
109#endif
110#ifdef ATM_PRESS
111 & Pair, &
112#endif
113#ifdef DIAGNOSTICS_UV
114!! & DiaRU, DiaRV, &
115#endif
116 & tl_ru, tl_rv)
117!***********************************************************************
118!
119 USE mod_param
120 USE mod_scalars
121!
122! Imported variable declarations.
123!
124 integer, intent(in) :: ng, tile
125 integer, intent(in) :: LBi, UBi, LBj, UBj
126 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
127 integer, intent(in) :: nrhs
128
129#ifdef ASSUMED_SHAPE
130 real(r8), intent(in) :: om_v(LBi:,LBj:)
131 real(r8), intent(in) :: on_u(LBi:,LBj:)
132 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
133 real(r8), intent(in) :: z_r(LBi:,LBj:,:)
134 real(r8), intent(in) :: z_w(LBi:,LBj:,0:)
135 real(r8), intent(in) :: rho(LBi:,LBj:,:)
136
137 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
138 real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:)
139 real(r8), intent(in) :: tl_z_w(LBi:,LBj:,0:)
140 real(r8), intent(in) :: tl_rho(LBi:,LBj:,:)
141
142# ifdef TIDE_GENERATING_FORCES
143 real(r8), intent(in) :: eq_tide(LBi:,LBj:)
144 real(r8), intent(in) :: tl_eq_tide(LBi:,LBj:)
145# endif
146# ifdef ATM_PRESS
147 real(r8), intent(in) :: Pair(LBi:,LBj:)
148# endif
149# ifdef DIAGNOSTICS_UV
150!! real(r8), intent(inout) :: DiaRU(LBi:,LBj:,:,:,:)
151!! real(r8), intent(inout) :: DiaRV(LBi:,LBj:,:,:,:)
152# endif
153 real(r8), intent(inout) :: tl_ru(LBi:,LBj:,0:,:)
154 real(r8), intent(inout) :: tl_rv(LBi:,LBj:,0:,:)
155#else
156 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
157 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
158 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
159 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
160 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
161 real(r8), intent(in) :: rho(LBi:UBi,LBj:UBj,N(ng))
162
163 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
164 real(r8), intent(in) :: tl_z_r(LBi:UBi,LBj:UBj,N(ng))
165 real(r8), intent(in) :: tl_z_w(LBi:UBi,LBj:UBj,0:N(ng))
166 real(r8), intent(in) :: tl_rho(LBi:UBi,LBj:UBj,N(ng))
167
168# ifdef TIDE_GENERATING_FORCES
169 real(r8), intent(in) :: eq_tide(LBi:,LBj:)
170 real(r8), intent(in) :: tl_eq_tide(LBi:,LBj:)
171# endif
172# ifdef ATM_PRESS
173 real(r8), intent(in) :: Pair(LBi:UBi,LBj:UBj)
174# endif
175# ifdef DIAGNOSTICS_UV
176!! real(r8), intent(inout) :: DiaRU(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
177!! real(r8), intent(inout) :: DiaRV(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
178# endif
179 real(r8), intent(inout) :: tl_ru(LBi:UBi,LBj:UBj,0:N(ng),2)
180 real(r8), intent(inout) :: tl_rv(LBi:UBi,LBj:UBj,0:N(ng),2)
181#endif
182!
183! Local variable declarations.
184!
185 integer :: i, j, k
186
187 real(r8) :: fac, fac1, fac2, fac3
188 real(r8) :: cff1, cff2, cff3, cff4
189 real(r8) :: tl_cff1, tl_cff2, tl_cff3, tl_cff4
190#ifdef WJ_GRADP
191 real(r8) :: gamma, tl_gamma
192#endif
193
194 real(r8), dimension(IminS:ImaxS) :: phie
195 real(r8), dimension(IminS:ImaxS) :: phix
196
197 real(r8), dimension(IminS:ImaxS) :: tl_phie
198 real(r8), dimension(IminS:ImaxS) :: tl_phix
199
200#include "set_bounds.h"
201!
202!-----------------------------------------------------------------------
203! Calculate pressure gradient in the XI-direction (m4/s2).
204!-----------------------------------------------------------------------
205!
206! Compute surface baroclinic pressure gradient.
207!
208#ifdef ATM_PRESS
209 fac=100.0_r8/rho0
210#endif
211 fac1=0.5_r8*g/rho0
212 fac2=1000.0_r8*g/rho0
213 fac3=0.25_r8*g/rho0
214
215 DO j=jstr,jend
216 DO i=istru,iend
217#ifdef TIDE_GENERATING_FORCES
218 cff1=z_w(i ,j,n(ng))-eq_tide(i ,j)-z_r(i ,j,n(ng))+ &
219 & z_w(i-1,j,n(ng))-eq_tide(i-1,j)-z_r(i-1,j,n(ng))
220 tl_cff1=tl_z_w(i ,j,n(ng))-tl_eq_tide(i ,j)- &
221 & tl_z_r(i ,j,n(ng))+ &
222 & tl_z_w(i-1,j,n(ng))-tl_eq_tide(i-1,j)- &
223 & tl_z_r(i-1,j,n(ng))
224#else
225 cff1=z_w(i ,j,n(ng))-z_r(i ,j,n(ng))+ &
226 & z_w(i-1,j,n(ng))-z_r(i-1,j,n(ng))
227 tl_cff1=tl_z_w(i ,j,n(ng))-tl_z_r(i ,j,n(ng))+ &
228 & tl_z_w(i-1,j,n(ng))-tl_z_r(i-1,j,n(ng))
229#endif
230 phix(i)=fac1*(rho(i,j,n(ng))-rho(i-1,j,n(ng)))*cff1
231 tl_phix(i)=fac1* &
232 & ((tl_rho(i,j,n(ng))-tl_rho(i-1,j,n(ng)))*cff1+ &
233 & (rho(i,j,n(ng))-rho(i-1,j,n(ng)))*tl_cff1)- &
234#ifdef TL_IOMS
235 & phix(i)
236#endif
237#ifdef ATM_PRESS
238 phix(i)=phix(i)+fac*(pair(i,j)-pair(i-1,j))
239#endif
240#ifdef RHO_SURF
241 phix(i)=phix(i)+ &
242 & (fac2+fac1*(rho(i,j,n(ng))+rho(i-1,j,n(ng))))* &
243 & (z_w(i,j,n(ng))-z_w(i-1,j,n(ng)))
244 tl_phix(i)=tl_phix(i)+ &
245 & (fac1*(tl_rho(i,j,n(ng))+tl_rho(i-1,j,n(ng))))* &
246 & (z_w(i,j,n(ng))-z_w(i-1,j,n(ng)))+ &
247 & (fac2+fac1*(rho(i,j,n(ng))+rho(i-1,j,n(ng))))* &
248 & (tl_z_w(i,j,n(ng))-tl_z_w(i-1,j,n(ng)))- &
249# ifdef TL_IOMS
250 & fac1*(rho(i,j,n(ng))+rho(i-1,j,n(ng)))* &
251 & (z_w(i,j,n(ng))-z_w(i-1,j,n(ng)))
252# endif
253#endif
254!^ ru(i,j,N(ng),nrhs)=-0.5_r8*(Hz(i,j,N(ng))+Hz(i-1,j,N(ng)))* &
255!^ & phix(i)*on_u(i,j)
256!^
257 tl_ru(i,j,n(ng),nrhs)=-0.5_r8*on_u(i,j)* &
258 & ((tl_hz(i ,j,n(ng))+ &
259 & tl_hz(i-1,j,n(ng)))*phix(i)+ &
260 & (hz(i ,j,n(ng))+ &
261 & hz(i-1,j,n(ng)))*tl_phix(i))+ &
262#ifdef TL_IOMS
263 & 0.5_r8*on_u(i,j)* &
264 & (hz(i ,j,n(ng))+ &
265 & hz(i-1,j,n(ng)))*phix(i)
266#endif
267#ifdef DIAGNOSTICS_UV
268!! DiaRU(i,j,N(ng),nrhs,M3pgrd)=ru(i,j,N(ng),nrhs)
269#endif
270 END DO
271!
272! Compute interior baroclinic pressure gradient. Differentiate and
273! then vertically integrate.
274!
275 DO k=n(ng)-1,1,-1
276 DO i=istru,iend
277#ifdef WJ_GRADP
278 cff1=1.0_r8/((z_r(i ,j,k+1)-z_r(i ,j,k))* &
279 & (z_r(i-1,j,k+1)-z_r(i-1,j,k)))
280 tl_cff1=-cff1*cff1*((tl_z_r(i ,j,k+1)-tl_z_r(i ,j,k))* &
281 & (z_r(i-1,j,k+1)-z_r(i-1,j,k))+ &
282 & (z_r(i ,j,k+1)-z_r(i ,j,k))* &
283 & (tl_z_r(i-1,j,k+1)-tl_z_r(i-1,j,k)))+ &
284# ifdef TL_IOMS
285 & 3.0_r8*cff1
286# endif
287 cff2=z_r(i ,j,k )-z_r(i-1,j,k )+ &
288 & z_r(i ,j,k+1)-z_r(i-1,j,k+1)
289 tl_cff2=tl_z_r(i ,j,k )-tl_z_r(i-1,j,k )+ &
290 & tl_z_r(i ,j,k+1)-tl_z_r(i-1,j,k+1)
291 cff3=z_r(i ,j,k+1)-z_r(i ,j,k )- &
292 & z_r(i-1,j,k+1)+z_r(i-1,j,k )
293 tl_cff3=tl_z_r(i ,j,k+1)-tl_z_r(i ,j,k )- &
294 & tl_z_r(i-1,j,k+1)+tl_z_r(i-1,j,k )
295 gamma=0.125_r8*cff1*cff2*cff3
296 tl_gamma=0.125_r8*(tl_cff1*cff2*cff3+ &
297 & cff1*(tl_cff2*cff3+ &
298 & cff2*tl_cff3))- &
299# ifdef TL_IOMS
300 & 2.0_r8*gamma
301# endif
302
303 cff1=(1.0_r8+gamma)*(rho(i,j,k+1)-rho(i-1,j,k+1))+ &
304 & (1.0_r8-gamma)*(rho(i,j,k )-rho(i-1,j,k ))
305 tl_cff1=tl_gamma*(rho(i,j,k+1)-rho(i-1,j,k+1)- &
306 & rho(i,j,k )+rho(i-1,j,k ))+ &
307 & (1.0_r8+gamma)*(tl_rho(i ,j,k+1)- &
308 & tl_rho(i-1,j,k+1))+ &
309 & (1.0_r8-gamma)*(tl_rho(i ,j,k )- &
310 & tl_rho(i-1,j,k ))- &
311# ifdef TL_IOMS
312 & gamma*((rho(i,j,k+1)-rho(i-1,j,k+1))- &
313 & (rho(i,j,k )-rho(i-1,j,k )))
314# endif
315 cff2=rho(i,j,k+1)+rho(i-1,j,k+1)- &
316 & rho(i,j,k )-rho(i-1,j,k )
317 tl_cff2=tl_rho(i,j,k+1)+tl_rho(i-1,j,k+1)- &
318 & tl_rho(i,j,k )-tl_rho(i-1,j,k )
319 cff3=z_r(i,j,k+1)+z_r(i-1,j,k+1)- &
320 & z_r(i,j,k )-z_r(i-1,j,k )
321 tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i-1,j,k+1)- &
322 & tl_z_r(i,j,k )-tl_z_r(i-1,j,k )
323 cff4=(1.0_r8+gamma)*(z_r(i,j,k+1)-z_r(i-1,j,k+1))+ &
324 & (1.0_r8-gamma)*(z_r(i,j,k )-z_r(i-1,j,k ))
325 tl_cff4=tl_gamma*(z_r(i,j,k+1)-z_r(i-1,j,k+1)- &
326 & z_r(i,j,k )+z_r(i-1,j,k ))+ &
327 & (1.0_r8+gamma)*(tl_z_r(i ,j,k+1)- &
328 & tl_z_r(i-1,j,k+1))+ &
329 & (1.0_r8-gamma)*(tl_z_r(i ,j,k )- &
330 & tl_z_r(i-1,j,k ))- &
331# ifdef TL_IOMS
332 & gamma*((z_r(i,j,k+1)-z_r(i-1,j,k+1))- &
333 & (z_r(i,j,k )-z_r(i-1,j,k )))
334# endif
335 phix(i)=phix(i)+ &
336 & fac3*(cff1*cff3-cff2*cff4)
337 tl_phix(i)=tl_phix(i)+ &
338 & fac3*(tl_cff1*cff3+ &
339 & cff1*tl_cff3- &
340 & tl_cff2*cff4- &
341 & cff2*tl_cff4)- &
342# ifdef TL_IOMS
343 & fac3*(cff1*cff3- &
344 & cff2*cff4)
345# endif
346#else
347 cff1=rho(i,j,k+1)-rho(i-1,j,k+1)+ &
348 & rho(i,j,k )-rho(i-1,j,k )
349 cff2=rho(i,j,k+1)+rho(i-1,j,k+1)- &
350 & rho(i,j,k )-rho(i-1,j,k )
351 tl_cff1=tl_rho(i,j,k+1)-tl_rho(i-1,j,k+1)+ &
352 & tl_rho(i,j,k )-tl_rho(i-1,j,k )
353 tl_cff2=tl_rho(i,j,k+1)+tl_rho(i-1,j,k+1)- &
354 & tl_rho(i,j,k )-tl_rho(i-1,j,k )
355 cff3=z_r(i,j,k+1)+z_r(i-1,j,k+1)- &
356 & z_r(i,j,k )-z_r(i-1,j,k )
357 cff4=z_r(i,j,k+1)-z_r(i-1,j,k+1)+ &
358 & z_r(i,j,k )-z_r(i-1,j,k )
359 tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i-1,j,k+1)- &
360 & tl_z_r(i,j,k )-tl_z_r(i-1,j,k )
361 tl_cff4=tl_z_r(i,j,k+1)-tl_z_r(i-1,j,k+1)+ &
362 & tl_z_r(i,j,k )-tl_z_r(i-1,j,k )
363 phix(i)=phix(i)+ &
364 & fac3*(cff1*cff3-cff2*cff4)
365 tl_phix(i)=tl_phix(i)+ &
366 & fac3*(tl_cff1*cff3+ &
367 & cff1*tl_cff3- &
368 & tl_cff2*cff4- &
369 & cff2*tl_cff4)- &
370# ifdef TL_IOMS
371 & fac3*(cff1*cff3- &
372 & cff2*cff4)
373# endif
374#endif
375!^ ru(i,j,k,nrhs)=-0.5_r8*(Hz(i,j,k)+Hz(i-1,j,k))* &
376!^ & phix(i)*on_u(i,j)
377!^
378 tl_ru(i,j,k,nrhs)=-0.5_r8*on_u(i,j)* &
379 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
380 & phix(i)+ &
381 & (hz(i,j,k)+hz(i-1,j,k))* &
382 & tl_phix(i))+ &
383#ifdef TL_IOMS
384 & 0.5_r8*on_u(i,j)* &
385 & (hz(i,j,k)+hz(i-1,j,k))*phix(i)
386#endif
387#ifdef DIAGNOSTICS_UV
388!! DiaRU(i,j,k,nrhs,M3pgrd)=ru(i,j,k,nrhs)
389#endif
390 END DO
391 END DO
392!
393!-----------------------------------------------------------------------
394! Calculate pressure gradient in the ETA-direction (m4/s2).
395!-----------------------------------------------------------------------
396!
397! Compute surface baroclinic pressure gradient.
398!
399 IF (j.ge.jstrv) THEN
400 DO i=istr,iend
401 cff1=z_w(i,j ,n(ng))-z_r(i,j ,n(ng))+ &
402 & z_w(i,j-1,n(ng))-z_r(i,j-1,n(ng))
403 tl_cff1=tl_z_w(i,j ,n(ng))-tl_z_r(i,j ,n(ng))+ &
404 & tl_z_w(i,j-1,n(ng))-tl_z_r(i,j-1,n(ng))
405 phie(i)=fac1*(rho(i,j,n(ng))-rho(i,j-1,n(ng)))*cff1
406 tl_phie(i)=fac1* &
407 & ((tl_rho(i,j,n(ng))-tl_rho(i,j-1,n(ng)))*cff1+ &
408 & (rho(i,j,n(ng))-rho(i,j-1,n(ng)))*tl_cff1)- &
409#ifdef TL_IOMS
410 & phie(i)
411#endif
412#ifdef ATM_PRESS
413 phie(i)=phie(i)+fac*(pair(i,j)-pair(i,j-1))
414#endif
415#ifdef RHO_SURF
416 phie(i)=phie(i)+ &
417 & (fac2+fac1*(rho(i,j,n(ng))+rho(i,j-1,n(ng))))* &
418 & (z_w(i,j,n(ng))-z_w(i,j-1,n(ng)))
419 tl_phie(i)=tl_phie(i)+ &
420 & (fac1*(tl_rho(i,j,n(ng))+tl_rho(i,j-1,n(ng))))* &
421 & (z_w(i,j,n(ng))-z_w(i,j-1,n(ng)))+ &
422 & (fac2+fac1*(rho(i,j,n(ng))+rho(i,j-1,n(ng))))* &
423 & (tl_z_w(i,j,n(ng))-tl_z_w(i,j-1,n(ng)))- &
424# ifdef TL_IOMS
425 & fac1*(rho(i,j,n(ng))+rho(i,j-1,n(ng)))* &
426 & (z_w(i,j,n(ng))-z_w(i,j-1,n(ng)))
427# endif
428#endif
429!^ rv(i,j,N(ng),nrhs)=-0.5_r8*(Hz(i,j,N(ng))+Hz(i,j-1,N(ng)))* &
430!^ & phie(i)*om_v(i,j)
431!^
432 tl_rv(i,j,n(ng),nrhs)=-0.5_r8*om_v(i,j)* &
433 & ((tl_hz(i,j ,n(ng))+ &
434 & tl_hz(i,j-1,n(ng)))*phie(i)+ &
435 & (hz(i,j ,n(ng))+ &
436 & hz(i,j-1,n(ng)))*tl_phie(i))+ &
437#ifdef TL_IOMS
438 & 0.5_r8*om_v(i,j)* &
439 & (hz(i,j ,n(ng))+ &
440 & hz(i,j-1,n(ng)))*phie(i)
441#endif
442#ifdef DIAGNOSTICS_UV
443!! DiaRV(i,j,N(ng),nrhs,M3pgrd)=rv(i,j,N(ng),nrhs)
444#endif
445 END DO
446!
447! Compute interior baroclinic pressure gradient. Differentiate and
448! then vertically integrate.
449!
450 DO k=n(ng)-1,1,-1
451 DO i=istr,iend
452#ifdef WJ_GRADP
453 cff1=1.0_r8/((z_r(i,j ,k+1)-z_r(i,j ,k))* &
454 & (z_r(i,j-1,k+1)-z_r(i,j-1,k)))
455 tl_cff1=-cff1*cff1*((tl_z_r(i,j ,k+1)-tl_z_r(i,j ,k))* &
456 & (z_r(i,j-1,k+1)-z_r(i,j-1,k))+ &
457 & (z_r(i,j ,k+1)-z_r(i,j ,k))* &
458 & (tl_z_r(i,j-1,k+1)-tl_z_r(i,j-1,k)))+ &
459# ifdef TL_IOMS
460 & 3.0_r8*cff1
461# endif
462 cff2=z_r(i,j ,k )-z_r(i,j-1,k )+ &
463 & z_r(i,j ,k+1)-z_r(i,j-1,k+1)
464 tl_cff2=tl_z_r(i,j ,k )-tl_z_r(i,j-1,k )+ &
465 & tl_z_r(i,j ,k+1)-tl_z_r(i,j-1,k+1)
466 cff3=z_r(i,j ,k+1)-z_r(i,j ,k )- &
467 & z_r(i,j-1,k+1)+z_r(i,j-1,k )
468 tl_cff3=tl_z_r(i,j ,k+1)-tl_z_r(i,j ,k )- &
469 & tl_z_r(i,j-1,k+1)+tl_z_r(i,j-1,k )
470 gamma=0.125_r8*cff1*cff2*cff3
471 tl_gamma=0.125_r8*(tl_cff1*cff2*cff3+ &
472 & cff1*(tl_cff2*cff3+ &
473 & cff2*tl_cff3))- &
474# ifdef TL_IOMS
475 & 2.0_r8*gamma
476# endif
477
478 cff1=(1.0_r8+gamma)*(rho(i,j,k+1)-rho(i,j-1,k+1))+ &
479 & (1.0_r8-gamma)*(rho(i,j,k )-rho(i,j-1,k ))
480 tl_cff1=tl_gamma*(rho(i,j,k+1)-rho(i,j-1,k+1)- &
481 & rho(i,j,k )+rho(i,j-1,k ))+ &
482 & (1.0_r8+gamma)*(tl_rho(i,j ,k+1)- &
483 & tl_rho(i,j-1,k+1))+ &
484 & (1.0_r8-gamma)*(tl_rho(i,j ,k )- &
485 & tl_rho(i,j-1,k ))- &
486# ifdef TL_IOMS
487 & gamma*((rho(i,j,k+1)-rho(i,j-1,k+1))- &
488 & (rho(i,j,k )-rho(i,j-1,k )))
489# endif
490 cff2=rho(i,j,k+1)+rho(i,j-1,k+1)- &
491 & rho(i,j,k )-rho(i,j-1,k )
492 tl_cff2=tl_rho(i,j,k+1)+tl_rho(i,j-1,k+1)- &
493 & tl_rho(i,j,k )-tl_rho(i,j-1,k )
494 cff3=z_r(i,j,k+1)+z_r(i,j-1,k+1)- &
495 & z_r(i,j,k )-z_r(i,j-1,k )
496 tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i,j-1,k+1)- &
497 & tl_z_r(i,j,k )-tl_z_r(i,j-1,k )
498 cff4=(1.0_r8+gamma)*(z_r(i,j,k+1)-z_r(i,j-1,k+1))+ &
499 & (1.0_r8-gamma)*(z_r(i,j,k )-z_r(i,j-1,k ))
500 tl_cff4=tl_gamma*(z_r(i,j,k+1)-z_r(i,j-1,k+1)- &
501 & z_r(i,j,k )+z_r(i,j-1,k ))+ &
502 & (1.0_r8+gamma)*(tl_z_r(i,j ,k+1)- &
503 & tl_z_r(i,j-1,k+1))+ &
504 & (1.0_r8-gamma)*(tl_z_r(i,j ,k )- &
505 & tl_z_r(i,j-1,k ))- &
506# ifdef TL_IOMS
507 & gamma*((z_r(i,j,k+1)-z_r(i,j-1,k+1))- &
508 & (z_r(i,j,k )-z_r(i,j-1,k )))
509# endif
510 phie(i)=phie(i)+ &
511 & fac3*(cff1*cff3-cff2*cff4)
512 tl_phie(i)=tl_phie(i)+ &
513 & fac3*(tl_cff1*cff3+ &
514 & cff1*tl_cff3- &
515 & tl_cff2*cff4- &
516 & cff2*tl_cff4)- &
517# ifdef TL_IOMS
518 & fac3*(cff1*cff3- &
519 & cff2*cff4)
520# endif
521#else
522 cff1=rho(i,j,k+1)-rho(i,j-1,k+1)+ &
523 & rho(i,j,k )-rho(i,j-1,k )
524 cff2=rho(i,j,k+1)+rho(i,j-1,k+1)- &
525 & rho(i,j,k )-rho(i,j-1,k )
526 tl_cff1=tl_rho(i,j,k+1)-tl_rho(i,j-1,k+1)+ &
527 & tl_rho(i,j,k )-tl_rho(i,j-1,k )
528 tl_cff2=tl_rho(i,j,k+1)+tl_rho(i,j-1,k+1)- &
529 & tl_rho(i,j,k )-tl_rho(i,j-1,k )
530 cff3=z_r(i,j,k+1)+z_r(i,j-1,k+1)- &
531 & z_r(i,j,k )-z_r(i,j-1,k )
532 cff4=z_r(i,j,k+1)-z_r(i,j-1,k+1)+ &
533 & z_r(i,j,k )-z_r(i,j-1,k )
534 tl_cff3=tl_z_r(i,j,k+1)+tl_z_r(i,j-1,k+1)- &
535 & tl_z_r(i,j,k )-tl_z_r(i,j-1,k )
536 tl_cff4=tl_z_r(i,j,k+1)-tl_z_r(i,j-1,k+1)+ &
537 & tl_z_r(i,j,k )-tl_z_r(i,j-1,k )
538 phie(i)=phie(i)+ &
539 & fac3*(cff1*cff3-cff2*cff4)
540 tl_phie(i)=tl_phie(i)+ &
541 & fac3*(tl_cff1*cff3+ &
542 & cff1*tl_cff3- &
543 & tl_cff2*cff4- &
544 & cff2*tl_cff4)- &
545# ifdef TL_IOMS
546 & fac3*(cff1*cff3- &
547 & cff2*cff4)
548# endif
549#endif
550!^ rv(i,j,k,nrhs)=-0.5_r8*(Hz(i,j,k)+Hz(i,j-1,k))* &
551!^ & phie(i)*om_v(i,j)
552!^
553 tl_rv(i,j,k,nrhs)=-0.5_r8*om_v(i,j)* &
554 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
555 & phie(i)+ &
556 & (hz(i,j,k)+hz(i,j-1,k))* &
557 & tl_phie(i))+ &
558#ifdef TL_IOMS
559 & 0.5_r8*om_v(i,j)* &
560 & (hz(i,j,k)+hz(i,j-1,k))*phie(i)
561#endif
562#ifdef DIAGNOSTICS_UV
563!! DiaRV(i,j,k,nrhs,M3pgrd)=rv(i,j,k,nrhs)
564#endif
565 END DO
566 END DO
567 END IF
568 END DO
569!
570 RETURN
571 END SUBROUTINE rp_prsgrd31_tile
572
573 END MODULE rp_prsgrd_mod
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 irpm
Definition mod_param.F:664
real(dp) g
real(dp) rho0
integer, dimension(:), allocatable nrhs
subroutine, public rp_prsgrd(ng, tile)
Definition rp_prsgrd31.h:36
subroutine rp_prsgrd31_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, om_v, on_u, hz, tl_hz, z_r, tl_z_r, z_w, tl_z_w, rho, tl_rho, eq_tide, tl_eq_tide, pair, tl_ru, tl_rv)
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