ROMS
Loading...
Searching...
No Matches
rp_prsgrd32.h
Go to the documentation of this file.
1 MODULE rp_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 routine evaluates the representers tangent linear baroclinic, !
11! hydrostatic pressure gradient term using a nonconservative Density !
12! Jacobian scheme, based on cubic polynomial fits for "rho" and !
13! "z_r" as functions of nondimensional coordinates (XI,ETA,s), that !
14! is, its respective array indices. The cubic polynomials are !
15! monotonized by using harmonic mean instead of linear averages to !
16! interpolate 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 :: rp_prsgrd
36!
37 CONTAINS
38!
39!***********************************************************************
40 SUBROUTINE rp_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, irpm, 23, __line__, myfile)
67#endif
68 CALL rp_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) % tl_Hz, &
80 & grid(ng) % z_r, &
81 & grid(ng) % tl_z_r, &
82 & grid(ng) % z_w, &
83 & grid(ng) % tl_z_w, &
84 & ocean(ng) % rho, &
85 & ocean(ng) % tl_rho, &
86#ifdef TIDE_GENERATING_FORCES
87 & ocean(ng) % eq_tide, &
88 & ocean(ng) % tl_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) % tl_ru, &
98 & ocean(ng) % tl_rv)
99#ifdef PROFILE
100 CALL wclock_off (ng, irpm, 23, __line__, myfile)
101#endif
102!
103 RETURN
104 END SUBROUTINE rp_prsgrd
105!
106!***********************************************************************
107 SUBROUTINE rp_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, tl_Hz, &
116 & z_r, tl_z_r, &
117 & z_w, tl_z_w, &
118 & rho, tl_rho, &
119#ifdef TIDE_GENERATING_FORCES
120 & eq_tide, tl_eq_tide, &
121#endif
122#ifdef ATM_PRESS
123 & Pair, &
124#endif
125#ifdef DIAGNOSTICS_UV
126!! & DiaRU, DiaRV, &
127#endif
128 & tl_ru, tl_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
153 real(r8), intent(in) :: tl_Hz(LBi:,LBj:,:)
154 real(r8), intent(in) :: tl_z_r(LBi:,LBj:,:)
155 real(r8), intent(in) :: tl_z_w(LBi:,LBj:,0:)
156 real(r8), intent(in) :: tl_rho(LBi:,LBj:,:)
157# ifdef TIDE_GENERATING_FORCES
158 real(r8), intent(in) :: eq_tide(LBi:,LBj:)
159 real(r8), intent(in) :: tl_eq_tide(LBi:,LBj:)
160# endif
161# ifdef ATM_PRESS
162 real(r8), intent(in) :: Pair(LBi:,LBj:)
163# endif
164# ifdef DIAGNOSTICS_UV
165!! real(r8), intent(inout) :: DiaRU(LBi:,LBj:,:,:,:)
166!! real(r8), intent(inout) :: DiaRV(LBi:,LBj:,:,:,:)
167# endif
168 real(r8), intent(inout) :: tl_ru(LBi:,LBj:,0:,:)
169 real(r8), intent(inout) :: tl_rv(LBi:,LBj:,0:,:)
170#else
171# ifdef MASKING
172 real(r8), intent(in) :: umask(LBi:UBi,LBj:UBj)
173 real(r8), intent(in) :: vmask(LBi:UBi,LBj:UBj)
174# endif
175 real(r8), intent(in) :: om_v(LBi:UBi,LBj:UBj)
176 real(r8), intent(in) :: on_u(LBi:UBi,LBj:UBj)
177 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
178 real(r8), intent(in) :: z_r(LBi:UBi,LBj:UBj,N(ng))
179 real(r8), intent(in) :: z_w(LBi:UBi,LBj:UBj,0:N(ng))
180 real(r8), intent(in) :: rho(LBi:UBi,LBj:UBj,N(ng))
181
182 real(r8), intent(in) :: tl_Hz(LBi:UBi,LBj:UBj,N(ng))
183 real(r8), intent(in) :: tl_z_r(LBi:UBi,LBj:UBj,N(ng))
184 real(r8), intent(in) :: tl_z_w(LBi:UBi,LBj:UBj,0:N(ng))
185 real(r8), intent(in) :: tl_rho(LBi:UBi,LBj:UBj,N(ng))
186# ifdef TIDE_GENERATING_FORCES
187 real(r8), intent(in) :: eq_tide(LBi:UBi,LBj:UBj)
188 real(r8), intent(in) :: tl_eq_tide(LBi:UBi,LBj:UBj)
189# endif
190# ifdef ATM_PRESS
191 real(r8), intent(in) :: Pair(LBi:UBi,LBj:UBj)
192# endif
193# ifdef DIAGNOSTICS_UV
194!! real(r8), intent(inout) :: DiaRU(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
195!! real(r8), intent(inout) :: DiaRV(LBi:UBi,LBj:UBj,N(ng),2,NDrhs)
196# endif
197 real(r8), intent(inout) :: tl_ru(LBi:UBi,LBj:UBj,0:N(ng),2)
198 real(r8), intent(inout) :: tl_rv(LBi:UBi,LBj:UBj,0:N(ng),2)
199#endif
200!
201! Local variable declarations.
202!
203 integer :: i, j, k
204
205 real(r8), parameter :: OneFifth = 0.2_r8
206 real(r8), parameter :: OneTwelfth = 1.0_r8/12.0_r8
207 real(r8), parameter :: eps = 1.0e-10_r8
208
209 real(r8) :: GRho, GRho0, HalfGRho
210 real(r8) :: cff, cff1, cff2
211 real(r8) :: tl_cff, tl_cff1, tl_cff2
212#ifdef ATM_PRESS
213 real(r8) :: OneAtm, fac
214#endif
215 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: P
216
217 real(r8), dimension(IminS:ImaxS,JminS:JmaxS,N(ng)) :: tl_P
218
219 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dR
220 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dR1
221 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dZ
222 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: dZ1
223
224 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_dR
225 real(r8), dimension(IminS:ImaxS,0:N(ng)) :: tl_dZ
226
227 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: FC
228 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: aux
229 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: dRx
230 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: dZx
231
232 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_FC
233 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_aux
234 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_dRx
235 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tl_dZx
236
237#include "set_bounds.h"
238!
239!-----------------------------------------------------------------------
240! Preliminary step (same for XI- and ETA-components:
241!-----------------------------------------------------------------------
242!
243 grho=g/rho0
244 grho0=1000.0_r8*grho
245 halfgrho=0.5_r8*grho
246#ifdef ATM_PRESS
247 oneatm=1013.25_r8 ! 1 atm = 1013.25 mb
248 fac=100.0_r8/rho0
249#endif
250!
251 DO j=jstrv-1,jend
252 DO k=1,n(ng)-1
253 DO i=istru-1,iend
254 dr(i,k)=rho(i,j,k+1)-rho(i,j,k)
255 tl_dr(i,k)=tl_rho(i,j,k+1)-tl_rho(i,j,k)
256 dz(i,k)=z_r(i,j,k+1)-z_r(i,j,k)
257 tl_dz(i,k)=tl_z_r(i,j,k+1)-tl_z_r(i,j,k)
258 END DO
259 END DO
260 DO i=istru-1,iend
261 dr(i,n(ng))=dr(i,n(ng)-1)
262 tl_dr(i,n(ng))=tl_dr(i,n(ng)-1)
263 dz(i,n(ng))=dz(i,n(ng)-1)
264 tl_dz(i,n(ng))=tl_dz(i,n(ng)-1)
265 dr(i,0)=dr(i,1)
266 tl_dr(i,0)=tl_dr(i,1)
267 dz(i,0)=dz(i,1)
268 tl_dz(i,0)=tl_dz(i,1)
269 END DO
270 DO k=n(ng),1,-1
271 DO i=istru-1,iend
272 cff=2.0_r8*dr(i,k)*dr(i,k-1)
273 tl_cff=2.0_r8*(tl_dr(i,k)*dr(i,k-1)+ &
274 & dr(i,k)*tl_dr(i,k-1))- &
275#ifdef TL_IOMS
276 & cff
277#endif
278 dr1(i,k)=dr(i,k)
279 IF (cff.gt.eps) THEN
280 dr(i,k)=cff/(dr(i,k)+dr(i,k-1))
281 tl_dr(i,k)=(tl_cff-dr(i,k)*(tl_dr(i,k)+tl_dr(i,k-1)))/ &
282 & (dr1(i,k)+dr(i,k-1))+ &
283#ifdef TL_IOMS
284 & dr(i,k)
285#endif
286 ELSE
287 dr(i,k)=0.0_r8
288 tl_dr(i,k)=0.0_r8
289 END IF
290 dz1(i,k)=dz(i,k)
291 dz(i,k)=2.0_r8*dz(i,k)*dz(i,k-1)/(dz(i,k)+dz(i,k-1))
292 tl_dz(i,k)=(2.0_r8*(tl_dz(i,k)*dz(i,k-1)+ &
293 & dz1(i,k)*tl_dz(i,k-1))- &
294 & dz(i,k)*(tl_dz(i,k)+tl_dz(i,k-1)))/ &
295 & (dz1(i,k)+dz(i,k-1))
296 END DO
297 END DO
298 DO i=istru-1,iend
299 cff1=1.0_r8/(z_r(i,j,n(ng))-z_r(i,j,n(ng)-1))
300 tl_cff1=-cff1*cff1*(tl_z_r(i,j,n(ng))-tl_z_r(i,j,n(ng)-1))+ &
301#ifdef TL_IOMS
302 & 2.0_r8*cff1
303#endif
304 cff2=0.5_r8*(rho(i,j,n(ng))-rho(i,j,n(ng)-1))* &
305 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))*cff1
306 tl_cff2=0.5_r8*((tl_rho(i,j,n(ng))-tl_rho(i,j,n(ng)-1))* &
307 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))*cff1+ &
308 & (rho(i,j,n(ng))-rho(i,j,n(ng)-1))* &
309 & ((tl_z_w(i,j,n(ng))-tl_z_r(i,j,n(ng)))*cff1+ &
310 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))*tl_cff1))- &
311#ifdef TL_IOMS
312 & 2.0_r8*cff2
313#endif
314 p(i,j,n(ng))=g*z_w(i,j,n(ng))+ &
315#ifdef ATM_PRESS
316 & fac*(pair(i,j)-oneatm)+ &
317#endif
318 & grho*(rho(i,j,n(ng))+cff2)* &
319 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))
320 tl_p(i,j,n(ng))=g*tl_z_w(i,j,n(ng))+ &
321 & grho*((tl_rho(i,j,n(ng))+tl_cff2)* &
322 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))+ &
323 & (rho(i,j,n(ng))+cff2)* &
324 & (tl_z_w(i,j,n(ng))-tl_z_r(i,j,n(ng))))- &
325#ifdef TL_IOMS
326 & grho*(rho(i,j,n(ng))+cff2)* &
327 & (z_w(i,j,n(ng))-z_r(i,j,n(ng)))
328#endif
329#ifdef TIDE_GENERATING_FORCES
330 p(i,j,n(ng))=p(i,j,n(ng))-g*eq_tide(i,j)
331 tl_p(i,j,n(ng))=tl_p(i,j,n(ng))-g*tl_eq_tide(i,j)
332#endif
333 END DO
334 DO k=n(ng)-1,1,-1
335 DO i=istru-1,iend
336 cff=halfgrho*((rho(i,j,k+1)+rho(i,j,k))* &
337 & (z_r(i,j,k+1)-z_r(i,j,k))- &
338 & onefifth* &
339 & ((dr(i,k+1)-dr(i,k))* &
340 & (z_r(i,j,k+1)-z_r(i,j,k)- &
341 & onetwelfth* &
342 & (dz(i,k+1)+dz(i,k)))- &
343 & (dz(i,k+1)-dz(i,k))* &
344 & (rho(i,j,k+1)-rho(i,j,k)- &
345 & onetwelfth* &
346 & (dr(i,k+1)+dr(i,k)))))
347 tl_cff=halfgrho*((tl_rho(i,j,k+1)+tl_rho(i,j,k))* &
348 & (z_r(i,j,k+1)-z_r(i,j,k))+ &
349 & (rho(i,j,k+1)+rho(i,j,k))* &
350 & (tl_z_r(i,j,k+1)-tl_z_r(i,j,k))- &
351 & onefifth* &
352 & ((tl_dr(i,k+1)-tl_dr(i,k))* &
353 & (z_r(i,j,k+1)-z_r(i,j,k)- &
354 & onetwelfth* &
355 & (dz(i,k+1)+dz(i,k)))+ &
356 & (dr(i,k+1)-dr(i,k))* &
357 & (tl_z_r(i,j,k+1)-tl_z_r(i,j,k)- &
358 & onetwelfth* &
359 & (tl_dz(i,k+1)+tl_dz(i,k)))- &
360 & (tl_dz(i,k+1)-tl_dz(i,k))* &
361 & (rho(i,j,k+1)-rho(i,j,k)- &
362 & onetwelfth* &
363 & (dr(i,k+1)+dr(i,k)))- &
364 & (dz(i,k+1)-dz(i,k))* &
365 & (tl_rho(i,j,k+1)-tl_rho(i,j,k)- &
366 & onetwelfth* &
367 & (tl_dr(i,k+1)+tl_dr(i,k)))))- &
368#ifdef TL_IOMS
369 & cff
370#endif
371 p(i,j,k)=p(i,j,k+1)+cff
372 tl_p(i,j,k)=tl_p(i,j,k+1)+tl_cff
373 END DO
374 END DO
375 END DO
376!
377!-----------------------------------------------------------------------
378! Compute XI-component pressure gradient term.
379!-----------------------------------------------------------------------
380!
381 DO k=n(ng),1,-1
382 DO j=jstr,jend
383 DO i=istru-1,iend+1
384 aux(i,j)=z_r(i,j,k)-z_r(i-1,j,k)
385 tl_aux(i,j)=tl_z_r(i,j,k)-tl_z_r(i-1,j,k)
386#ifdef MASKING
387 aux(i,j)=aux(i,j)*umask(i,j)
388 tl_aux(i,j)=tl_aux(i,j)*umask(i,j)
389#endif
390 fc(i,j)=rho(i,j,k)-rho(i-1,j,k)
391 tl_fc(i,j)=tl_rho(i,j,k)-tl_rho(i-1,j,k)
392#ifdef MASKING
393 fc(i,j)=fc(i,j)*umask(i,j)
394 tl_fc(i,j)=tl_fc(i,j)*umask(i,j)
395#endif
396 END DO
397 END DO
398!
399 DO j=jstr,jend
400 DO i=istru-1,iend
401 cff=2.0_r8*aux(i,j)*aux(i+1,j)
402 tl_cff=2.0_r8*(tl_aux(i,j)*aux(i+1,j)+ &
403 & aux(i,j)*tl_aux(i+1,j))- &
404#ifdef TL_IOMS
405 & cff
406#endif
407 IF (cff.gt.eps) THEN
408 cff1=1.0_r8/(aux(i,j)+aux(i+1,j))
409 tl_cff1=-cff1*cff1*(tl_aux(i,j)+tl_aux(i+1,j))+ &
410#ifdef TL_IOMS
411 & 2.0_r8*cff1
412#endif
413 dzx(i,j)=cff*cff1
414 tl_dzx(i,j)=tl_cff*cff1+cff*tl_cff1- &
415#ifdef TL_IOMS
416 & dzx(i,j)
417#endif
418 ELSE
419 dzx(i,j)=0.0_r8
420 tl_dzx(i,j)=0.0_r8
421 END IF
422 cff1=2.0_r8*fc(i,j)*fc(i+1,j)
423 tl_cff1=2.0_r8*(tl_fc(i,j)*fc(i+1,j)+ &
424 & fc(i,j)*tl_fc(i+1,j))- &
425#ifdef TL_IOMS
426 & cff1
427#endif
428 IF (cff1.gt.eps) THEN
429 cff2=1.0_r8/(fc(i,j)+fc(i+1,j))
430 tl_cff2=-cff2*cff2*(tl_fc(i,j)+tl_fc(i+1,j))+ &
431#ifdef TL_IOMS
432 & 2.0_r8*cff2
433#endif
434 drx(i,j)=cff1*cff2
435 tl_drx(i,j)=tl_cff1*cff2+cff1*tl_cff2- &
436#ifdef TL_IOMS
437 & drx(i,j)
438#endif
439 ELSE
440 drx(i,j)=0.0_r8
441 tl_drx(i,j)=0.0_r8
442 END IF
443 END DO
444 END DO
445!
446 DO j=jstr,jend
447 DO i=istru,iend
448!^ ru(i,j,k,nrhs)=on_u(i,j)*0.5_r8* &
449!^ & (Hz(i,j,k)+Hz(i-1,j,k))* &
450!^ & (P(i-1,j,k)-P(i,j,k)- &
451!^ & HalfGRho* &
452!^ & ((rho(i,j,k)+rho(i-1,j,k))* &
453!^ & (z_r(i,j,k)-z_r(i-1,j,k))- &
454!^ & OneFifth* &
455!^ & ((dRx(i,j)-dRx(i-1,j))* &
456!^ & (z_r(i,j,k)-z_r(i-1,j,k)- &
457!^ & OneTwelfth* &
458!^ & (dZx(i,j)+dZx(i-1,j)))- &
459!^ & (dZx(i,j)-dZx(i-1,j))* &
460!^ & (rho(i,j,k)-rho(i-1,j,k)- &
461!^ & OneTwelfth* &
462!^ & (dRx(i,j)+dRx(i-1,j))))))
463!^
464 tl_ru(i,j,k,nrhs)=on_u(i,j)*0.5_r8* &
465 & ((tl_hz(i,j,k)+tl_hz(i-1,j,k))* &
466 & (p(i-1,j,k)-p(i,j,k)- &
467 & halfgrho* &
468 & ((rho(i,j,k)+rho(i-1,j,k))* &
469 & (z_r(i,j,k)-z_r(i-1,j,k))- &
470 & onefifth* &
471 & ((drx(i,j)-drx(i-1,j))* &
472 & (z_r(i,j,k)-z_r(i-1,j,k)- &
473 & onetwelfth* &
474 & (dzx(i,j)+dzx(i-1,j)))- &
475 & (dzx(i,j)-dzx(i-1,j))* &
476 & (rho(i,j,k)-rho(i-1,j,k)- &
477 & onetwelfth* &
478 & (drx(i,j)+drx(i-1,j))))))+ &
479 & (hz(i,j,k)+hz(i-1,j,k))* &
480 & (tl_p(i-1,j,k)-tl_p(i,j,k)- &
481 & halfgrho* &
482 & ((tl_rho(i,j,k)+tl_rho(i-1,j,k))* &
483 & (z_r(i,j,k)-z_r(i-1,j,k))+ &
484 & (rho(i,j,k)+rho(i-1,j,k))* &
485 & (tl_z_r(i,j,k)-tl_z_r(i-1,j,k))- &
486 & onefifth* &
487 & ((tl_drx(i,j)-tl_drx(i-1,j))* &
488 & (z_r(i,j,k)-z_r(i-1,j,k)- &
489 & onetwelfth* &
490 & (dzx(i,j)+dzx(i-1,j)))+ &
491 & (drx(i,j)-drx(i-1,j))* &
492 & (tl_z_r(i,j,k)-tl_z_r(i-1,j,k)- &
493 & onetwelfth* &
494 & (tl_dzx(i,j)+tl_dzx(i-1,j)))- &
495 & (tl_dzx(i,j)-tl_dzx(i-1,j))* &
496 & (rho(i,j,k)-rho(i-1,j,k)- &
497 & onetwelfth* &
498 & (drx(i,j)+drx(i-1,j)))- &
499 & (dzx(i,j)-dzx(i-1,j))* &
500 & (tl_rho(i,j,k)-tl_rho(i-1,j,k)- &
501 & onetwelfth* &
502 & (tl_drx(i,j)+tl_drx(i-1,j)))))))- &
503#ifdef TL_IOMS
504 & on_u(i,j)*0.5_r8* &
505 & (hz(i,j,k)+hz(i-1,j,k))* &
506 & (p(i-1,j,k)-p(i,j,k)- &
507 & 2.0_r8*halfgrho* &
508 & ((rho(i,j,k)+rho(i-1,j,k))* &
509 & (z_r(i,j,k)-z_r(i-1,j,k))- &
510 & onefifth* &
511 & ((drx(i,j)-drx(i-1,j))* &
512 & (z_r(i,j,k)-z_r(i-1,j,k)- &
513 & onetwelfth* &
514 & (dzx(i,j)+dzx(i-1,j)))- &
515 & (dzx(i,j)-dzx(i-1,j))* &
516 & (rho(i,j,k)-rho(i-1,j,k)- &
517 & onetwelfth* &
518 & (drx(i,j)+drx(i-1,j))))))
519#endif
520#ifdef DIAGNOSTICS_UV
521!! DiaRU(i,j,k,nrhs,M3pgrd)=ru(i,j,k,nrhs)
522#endif
523 END DO
524 END DO
525 END DO
526!
527!-----------------------------------------------------------------------
528! ETA-component pressure gradient term.
529!-----------------------------------------------------------------------
530!
531 DO k=n(ng),1,-1
532 DO j=jstrv-1,jend+1
533 DO i=istr,iend
534 aux(i,j)=z_r(i,j,k)-z_r(i,j-1,k)
535 tl_aux(i,j)=tl_z_r(i,j,k)-tl_z_r(i,j-1,k)
536#ifdef MASKING
537 aux(i,j)=aux(i,j)*vmask(i,j)
538 tl_aux(i,j)=tl_aux(i,j)*vmask(i,j)
539#endif
540 fc(i,j)=rho(i,j,k)-rho(i,j-1,k)
541 tl_fc(i,j)=tl_rho(i,j,k)-tl_rho(i,j-1,k)
542#ifdef MASKING
543 fc(i,j)=fc(i,j)*vmask(i,j)
544 tl_fc(i,j)=tl_fc(i,j)*vmask(i,j)
545#endif
546 END DO
547 END DO
548!
549 DO j=jstrv-1,jend
550 DO i=istr,iend
551 cff=2.0_r8*aux(i,j)*aux(i,j+1)
552 tl_cff=2.0_r8*(tl_aux(i,j)*aux(i,j+1)+ &
553 & aux(i,j)*tl_aux(i,j+1))- &
554#ifdef TL_IOMS
555 & cff
556#endif
557 IF (cff.gt.eps) THEN
558 cff1=1.0_r8/(aux(i,j)+aux(i,j+1))
559 tl_cff1=-cff1*cff1*(tl_aux(i,j)+tl_aux(i,j+1))+ &
560#ifdef TL_IOMS
561 & 2.0_r8*cff1
562#endif
563 dzx(i,j)=cff*cff1
564 tl_dzx(i,j)=tl_cff*cff1+cff*tl_cff1- &
565#ifdef TL_IOMS
566 & dzx(i,j)
567#endif
568 ELSE
569 dzx(i,j)=0.0_r8
570 tl_dzx(i,j)=0.0_r8
571 END IF
572 cff1=2.0_r8*fc(i,j)*fc(i,j+1)
573 tl_cff1=2.0_r8*(tl_fc(i,j)*fc(i,j+1)+ &
574 & fc(i,j)*tl_fc(i,j+1))- &
575#ifdef TL_IOMS
576 & cff1
577#endif
578 IF (cff1.gt.eps) THEN
579 cff2=1.0_r8/(fc(i,j)+fc(i,j+1))
580 tl_cff2=-cff2*cff2*(tl_fc(i,j)+tl_fc(i,j+1))+ &
581#ifdef TL_IOMS
582 & 2.0_r8*cff2
583#endif
584 drx(i,j)=cff1*cff2
585 tl_drx(i,j)=tl_cff1*cff2+cff1*tl_cff2- &
586#ifdef TL_IOMS
587 & drx(i,j)
588#endif
589 ELSE
590 drx(i,j)=0.0_r8
591 tl_drx(i,j)=0.0_r8
592 END IF
593 END DO
594 END DO
595!
596 DO j=jstrv,jend
597 DO i=istr,iend
598!^ rv(i,j,k,nrhs)=om_v(i,j)*0.5_r8* &
599!^ & (Hz(i,j,k)+Hz(i,j-1,k))* &
600!^ & (P(i,j-1,k)-P(i,j,k)- &
601!^ & HalfGRho* &
602!^ & ((rho(i,j,k)+rho(i,j-1,k))* &
603!^ & (z_r(i,j,k)-z_r(i,j-1,k))- &
604!^ & OneFifth* &
605!^ & ((dRx(i,j)-dRx(i,j-1))* &
606!^ & (z_r(i,j,k)-z_r(i,j-1,k)- &
607!^ & OneTwelfth* &
608!^ & (dZx(i,j)+dZx(i,j-1)))- &
609!^ & (dZx(i,j)-dZx(i,j-1))* &
610!^ & (rho(i,j,k)-rho(i,j-1,k)- &
611!^ & OneTwelfth* &
612!^ & (dRx(i,j)+dRx(i,j-1))))))
613!^
614 tl_rv(i,j,k,nrhs)=om_v(i,j)*0.5_r8* &
615 & ((tl_hz(i,j,k)+tl_hz(i,j-1,k))* &
616 & (p(i,j-1,k)-p(i,j,k)- &
617 & halfgrho* &
618 & ((rho(i,j,k)+rho(i,j-1,k))* &
619 & (z_r(i,j,k)-z_r(i,j-1,k))- &
620 & onefifth* &
621 & ((drx(i,j)-drx(i,j-1))* &
622 & (z_r(i,j,k)-z_r(i,j-1,k)- &
623 & onetwelfth* &
624 & (dzx(i,j)+dzx(i,j-1)))- &
625 & (dzx(i,j)-dzx(i,j-1))* &
626 & (rho(i,j,k)-rho(i,j-1,k)- &
627 & onetwelfth* &
628 & (drx(i,j)+drx(i,j-1))))))+ &
629 & (hz(i,j,k)+hz(i,j-1,k))* &
630 & (tl_p(i,j-1,k)-tl_p(i,j,k)- &
631 & halfgrho* &
632 & ((tl_rho(i,j,k)+tl_rho(i,j-1,k))* &
633 & (z_r(i,j,k)-z_r(i,j-1,k))+ &
634 & (rho(i,j,k)+rho(i,j-1,k))* &
635 & (tl_z_r(i,j,k)-tl_z_r(i,j-1,k))- &
636 & onefifth* &
637 & ((tl_drx(i,j)-tl_drx(i,j-1))* &
638 & (z_r(i,j,k)-z_r(i,j-1,k)- &
639 & onetwelfth* &
640 & (dzx(i,j)+dzx(i,j-1)))+ &
641 & (drx(i,j)-drx(i,j-1))* &
642 & (tl_z_r(i,j,k)-tl_z_r(i,j-1,k)- &
643 & onetwelfth* &
644 & (tl_dzx(i,j)+tl_dzx(i,j-1)))- &
645 & (tl_dzx(i,j)-tl_dzx(i,j-1))* &
646 & (rho(i,j,k)-rho(i,j-1,k)- &
647 & onetwelfth* &
648 & (drx(i,j)+drx(i,j-1)))- &
649 & (dzx(i,j)-dzx(i,j-1))* &
650 & (tl_rho(i,j,k)-tl_rho(i,j-1,k)- &
651 & onetwelfth* &
652 & (tl_drx(i,j)+tl_drx(i,j-1)))))))- &
653#ifdef TL_IOMS
654 & om_v(i,j)*0.5_r8* &
655 & (hz(i,j,k)+hz(i,j-1,k))* &
656 & (p(i,j-1,k)-p(i,j,k)- &
657 & 2.0_r8*halfgrho* &
658 & ((rho(i,j,k)+rho(i,j-1,k))* &
659 & (z_r(i,j,k)-z_r(i,j-1,k))- &
660 & onefifth* &
661 & ((drx(i,j)-drx(i,j-1))* &
662 & (z_r(i,j,k)-z_r(i,j-1,k)- &
663 & onetwelfth* &
664 & (dzx(i,j)+dzx(i,j-1)))- &
665 & (dzx(i,j)-dzx(i,j-1))* &
666 & (rho(i,j,k)-rho(i,j-1,k)- &
667 & onetwelfth* &
668 & (drx(i,j)+drx(i,j-1))))))
669#endif
670#ifdef DIAGNOSTICS_UV
671!! DiaRV(i,j,k,nrhs,M3pgrd)=rv(i,j,k,nrhs)
672#endif
673 END DO
674 END DO
675 END DO
676!
677 RETURN
678 END SUBROUTINE rp_prsgrd32_tile
679
680 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_prsgrd32_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nrhs, umask, vmask, 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