1 |
C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_gwdrag.F,v 1.1 2005/05/20 23:50:52 molod Exp $ |
2 |
C $Name: $ |
3 |
#include "FIZHI_OPTIONS.h" |
4 |
subroutine gwdrag (myid,pz,pl,ple,dpres,pkz,uz,vz,tz,qz,phis_var, |
5 |
. dudt,dvdt,dtdt,im,jm,lm,bi,bj,istrip,npcs,imglobal) |
6 |
C*********************************************************************** |
7 |
C |
8 |
C PURPOSE: |
9 |
C ======== |
10 |
C Driver Routine for Gravity Wave Drag |
11 |
C |
12 |
C INPUT: |
13 |
C ====== |
14 |
C myid ....... Process ID |
15 |
C pz ....... Surface Pressure [im,jm] |
16 |
C pl ....... 3D pressure field [im,jm,lm] |
17 |
C ple ....... 3d pressure at model level edges [im,jm,lm+1] |
18 |
C dpres ....... pressure difference across level [im,jm,lm] |
19 |
C pkz ....... pressure**kappa [im,jm,lm] |
20 |
C uz ....... zonal velocity [im,jm,lm] |
21 |
C vz ....... meridional velocity [im,jm,lm] |
22 |
C tz ....... temperature [im,jm,lm] |
23 |
C qz ....... specific humidity [im,jm,lm] |
24 |
C phis_var .... topography variance |
25 |
C im ....... number of grid points in x direction |
26 |
C jm ....... number of grid points in y direction |
27 |
C lm ....... number of grid points in vertical |
28 |
C istrip ...... 'strip' length for cache size control |
29 |
C npcs ....... number of strips |
30 |
C imglobal .... (avg) number of longitude points around the globe |
31 |
C |
32 |
C INPUT/OUTPUT: |
33 |
C ============ |
34 |
C dudt ....... Updated U-Wind Tendency including Gravity Wave Drag |
35 |
C dvdt ....... Updated V-Wind Tendency including Gravity Wave Drag |
36 |
C dtdt ....... Updated Pi*Theta Tendency including Gravity Wave Drag |
37 |
C |
38 |
C*********************************************************************** |
39 |
implicit none |
40 |
|
41 |
c Input Variables |
42 |
c --------------- |
43 |
integer myid,im,jm,lm,bi,bj,istrip,npcs,imglobal |
44 |
_RL pz(im,jm) |
45 |
_RL pl(im,jm,lm) |
46 |
_RL ple(im,jm,lm+1) |
47 |
_RL dpres(im,jm,lm) |
48 |
_RL pkz(im,jm,lm) |
49 |
_RL uz(im,jm,lm) |
50 |
_RL vz(im,jm,lm) |
51 |
_RL tz(im,jm,lm) |
52 |
_RL qz(im,jm,lm) |
53 |
_RL phis_var(im,jm) |
54 |
|
55 |
_RL dudt(im,jm,lm) |
56 |
_RL dvdt(im,jm,lm) |
57 |
_RL dtdt(im,jm,lm) |
58 |
|
59 |
c Local Variables |
60 |
c --------------- |
61 |
_RL tv(im,jm,lm) |
62 |
_RL dragu(im,jm,lm), dragv(im,jm,lm) |
63 |
_RL dragt(im,jm,lm) |
64 |
_RL dragx(im,jm), dragy(im,jm) |
65 |
_RL sumu(im,jm) |
66 |
integer nthin(im,jm),nbase(im,jm) |
67 |
integer nthini, nbasei |
68 |
|
69 |
_RL phis_std(im,jm) |
70 |
|
71 |
_RL std(istrip), ps(istrip) |
72 |
_RL us(istrip,lm), vs(istrip,lm), ts(istrip,lm) |
73 |
_RL dragus(istrip,lm), dragvs(istrip,lm) |
74 |
_RL dragxs(istrip), dragys(istrip) |
75 |
_RL plstr(istrip,lm),plestr(istrip,lm),dpresstr(istrip,lm) |
76 |
integer nthinstr(istrip),nbasestr(istrip) |
77 |
|
78 |
integer n,i,j,L |
79 |
_RL getcon, pi |
80 |
_RL grav, rgas, cp, cpinv, lstar |
81 |
#ifdef ALLOW_DIAGNOSTICS |
82 |
logical diagnostics_is_on |
83 |
external diagnostics_is_on |
84 |
_RL tmpdiag(im,jm) |
85 |
#endif |
86 |
|
87 |
c Initialization |
88 |
c -------------- |
89 |
pi = 4.0*atan(1.0) |
90 |
grav = getcon('GRAVITY') |
91 |
rgas = getcon('RGAS') |
92 |
cp = getcon('CP') |
93 |
cpinv = 1.0/cp |
94 |
lstar = 2*getcon('EARTH RADIUS')*cos(pi/3.0)/imglobal |
95 |
|
96 |
c Compute NTHIN and NBASE |
97 |
c ----------------------- |
98 |
do j=1,jm |
99 |
do i=1,im |
100 |
|
101 |
do nthini = 1,lm+1 |
102 |
if( 1000.0-ple(i,j,lm+2-nthini).gt.25. ) then |
103 |
nthin(i,j) = nthini |
104 |
goto 10 |
105 |
endif |
106 |
enddo |
107 |
10 continue |
108 |
do nbasei = 1,lm+1 |
109 |
if( ple(i,j,lm+2-nbasei).lt.666.7 ) then |
110 |
nbase(i,j) = nbasei |
111 |
goto 20 |
112 |
endif |
113 |
enddo |
114 |
20 continue |
115 |
if( 666.7-ple(i,j,lm+2-nbase(i,j)) .gt. |
116 |
. ple(i,j,lm+3-nbase(i,j))-666.7 ) then |
117 |
nbase(i,j) = nbase(i,j)-1 |
118 |
endif |
119 |
|
120 |
enddo |
121 |
enddo |
122 |
c Compute Topography Sub-Grid Standard Deviation |
123 |
c ---------------------------------------------- |
124 |
do j=1,jm |
125 |
do i=1,im |
126 |
phis_std(i,j) = min( 400.0, sqrt( max(0.0,phis_var(i,j)) )/grav ) |
127 |
enddo |
128 |
enddo |
129 |
|
130 |
c Compute Virtual Temperatures |
131 |
c ---------------------------- |
132 |
do L = 1,lm |
133 |
do j = 1,jm |
134 |
do i = 1,im |
135 |
tv(i,j,L) = tz(i,j,L)*pkz(i,j,L)*(1.+.609*qz(i,j,L)) |
136 |
enddo |
137 |
enddo |
138 |
enddo |
139 |
|
140 |
c Call Gravity Wave Drag Paramterization on A-Grid |
141 |
c ------------------------------------------------ |
142 |
|
143 |
do n=1,npcs |
144 |
|
145 |
call strip ( phis_std,std,im*jm,istrip,1,n ) |
146 |
|
147 |
call strip ( pz,ps,im*jm,istrip,1 ,n ) |
148 |
call strip ( uz,us,im*jm,istrip,lm,n ) |
149 |
call strip ( vz,vs,im*jm,istrip,lm,n ) |
150 |
call strip ( tv,ts,im*jm,istrip,lm,n ) |
151 |
call strip ( pl,plstr,im*jm,istrip,lm,n ) |
152 |
call strip ( ple,plestr,im*jm,istrip,lm,n ) |
153 |
call strip ( dpres,dpresstr,im*jm,istrip,lm,n ) |
154 |
call stripint ( nthin,nthinstr,im*jm,istrip,lm,n ) |
155 |
call stripint ( nbase,nbasestr,im*jm,istrip,lm,n ) |
156 |
|
157 |
call GWDD ( ps,us,vs,ts, |
158 |
. dragus,dragvs,dragxs,dragys,std, |
159 |
. plstr,plestr,dpresstr,grav,rgas,cp, |
160 |
. istrip,lm,nthinstr,nbasestr,lstar ) |
161 |
|
162 |
call paste ( dragus,dragu,istrip,im*jm,lm,n ) |
163 |
call paste ( dragvs,dragv,istrip,im*jm,lm,n ) |
164 |
call paste ( dragxs,dragx,istrip,im*jm,1 ,n ) |
165 |
call paste ( dragys,dragy,istrip,im*jm,1 ,n ) |
166 |
|
167 |
enddo |
168 |
|
169 |
c Add Gravity-Wave Drag to Wind and Theta Tendencies |
170 |
c -------------------------------------------------- |
171 |
do L = 1,lm |
172 |
do j = 1,jm |
173 |
do i = 1,im |
174 |
dragu(i,j,L) = sign( min(0.006,abs(dragu(i,j,L))),dragu(i,j,L) ) |
175 |
dragv(i,j,L) = sign( min(0.006,abs(dragv(i,j,L))),dragv(i,j,L) ) |
176 |
dragt(i,j,L) = -( uz(i,j,L)*dragu(i,j,L)+vz(i,j,L)*dragv(i,j,L) ) |
177 |
. *cpinv |
178 |
dudt(i,j,L) = dudt(i,j,L) + dragu(i,j,L) |
179 |
dvdt(i,j,L) = dvdt(i,j,L) + dragv(i,j,L) |
180 |
dtdt(i,j,L) = dtdt(i,j,L) + dragt(i,j,L)*pz(i,j)/pkz(i,j,L) |
181 |
enddo |
182 |
enddo |
183 |
enddo |
184 |
|
185 |
c Compute Diagnostics |
186 |
c ------------------- |
187 |
#ifdef ALLOW_DIAGNOSTICS |
188 |
do L = 1,lm |
189 |
|
190 |
if(diagnostics_is_on('GWDU ',myid) ) then |
191 |
do j=1,jm |
192 |
do i=1,im |
193 |
tmpdiag(i,j) = dragu(i,j,L)*86400 |
194 |
enddo |
195 |
enddo |
196 |
call diagnostics_fill(tmpdiag,'GWDU ',L,1,3,bi,bj,myid) |
197 |
endif |
198 |
|
199 |
if(diagnostics_is_on('GWDV ',myid) ) then |
200 |
do j=1,jm |
201 |
do i=1,im |
202 |
tmpdiag(i,j) = dragv(i,j,L)*86400 |
203 |
enddo |
204 |
enddo |
205 |
call diagnostics_fill(tmpdiag,'GWDV ',L,1,3,bi,bj,myid) |
206 |
endif |
207 |
|
208 |
if(diagnostics_is_on('GWDT ',myid) ) then |
209 |
do j=1,jm |
210 |
do i=1,im |
211 |
tmpdiag(i,j) = dragt(i,j,L)*86400 |
212 |
enddo |
213 |
enddo |
214 |
call diagnostics_fill(tmpdiag,'GWDT ',L,1,3,bi,bj,myid) |
215 |
endif |
216 |
|
217 |
enddo |
218 |
|
219 |
c Gravity Wave Drag at Surface (U-Wind) |
220 |
c ------------------------------------- |
221 |
if(diagnostics_is_on('GWDUS ',myid) ) then |
222 |
call diagnostics_fill(dragx,'GWDUS ',0,1,3,bi,bj,myid) |
223 |
endif |
224 |
|
225 |
c Gravity Wave Drag at Surface (V-Wind) |
226 |
c ------------------------------------- |
227 |
if(diagnostics_is_on('GWDVS ',myid) ) then |
228 |
call diagnostics_fill(dragy,'GWDVS ',0,1,3,bi,bj,myid) |
229 |
endif |
230 |
|
231 |
c Gravity Wave Drag at Model Top (U-Wind) |
232 |
c --------------------------------------- |
233 |
if(diagnostics_is_on('GWDUT ',myid) ) then |
234 |
do j = 1,jm |
235 |
do i = 1,im |
236 |
sumu(i,j) = 0.0 |
237 |
enddo |
238 |
enddo |
239 |
do L = 1,lm |
240 |
do j = 1,jm |
241 |
do i = 1,im |
242 |
sumu(i,j) = sumu(i,j) + dragu(i,j,L)*dpres(i,j,L)/pz(i,j) |
243 |
enddo |
244 |
enddo |
245 |
enddo |
246 |
do j=1,jm |
247 |
do i=1,im |
248 |
tmpdiag(i,j) = dragx(i,j) + sumu(i,j)*pz(i,j)/grav*100 |
249 |
enddo |
250 |
enddo |
251 |
call diagnostics_fill(tmpdiag,'GWDUT ',0,1,3,bi,bj,myid) |
252 |
endif |
253 |
|
254 |
c Gravity Wave Drag at Model Top (V-Wind) |
255 |
c --------------------------------------- |
256 |
if(diagnostics_is_on('GWDVT ',myid) ) then |
257 |
do j = 1,jm |
258 |
do i = 1,im |
259 |
sumu(i,j) = 0.0 |
260 |
enddo |
261 |
enddo |
262 |
do L = 1,lm |
263 |
do j = 1,jm |
264 |
do i = 1,im |
265 |
sumu(i,j) = sumu(i,j) + dragv(i,j,L)*dpres(i,j,L)/pz(i,j) |
266 |
enddo |
267 |
enddo |
268 |
enddo |
269 |
do j=1,jm |
270 |
do i=1,im |
271 |
tmpdiag(i,j) = dragy(i,j) + sumu(i,j)*pz(i,j)/grav*100 |
272 |
enddo |
273 |
enddo |
274 |
call diagnostics_fill(tmpdiag,'GWDVT ',0,1,3,bi,bj,myid) |
275 |
endif |
276 |
#endif |
277 |
|
278 |
return |
279 |
end |
280 |
SUBROUTINE GWDD ( ps,u,v,t,dudt,dvdt,xdrag,ydrag, |
281 |
. std,pl,ple,dpres, |
282 |
. grav,rgas,cp,irun,lm,nthin,nbase,lstar ) |
283 |
C*********************************************************************** |
284 |
C |
285 |
C Description: |
286 |
C ============ |
287 |
C Parameterization to introduce a Gravity Wave Drag |
288 |
C due to sub-grid scale orographic forcing |
289 |
C |
290 |
C Input: |
291 |
C ====== |
292 |
C ps ......... Surface Pressure |
293 |
C u .......... Zonal Wind (m/sec) |
294 |
C v .......... Meridional Wind (m/sec) |
295 |
C t .......... Virtual Temperature (deg K) |
296 |
C std ........ Standard Deviation of sub-grid Orography (m) |
297 |
C ple ....... Model pressure Edge Values |
298 |
C pl ........ Model pressure Values |
299 |
C dpres....... Model Delta pressure Values |
300 |
C grav ....... Gravitational constant (m/sec**2) |
301 |
C rgas ....... Gas constant |
302 |
C cp ......... Specific Heat at constant pressure |
303 |
C irun ....... Number of grid-points in horizontal dimension |
304 |
C lm ......... Number of grid-points in vertical dimension |
305 |
C lstar ...... Monochromatic Wavelength/(2*pi) |
306 |
C |
307 |
C Output: |
308 |
C ======= |
309 |
C dudt ....... Zonal Acceleration due to GW Drag (m/sec**2) |
310 |
C dvdt ....... Meridional Acceleration due to GW Drag (m/sec**2) |
311 |
C xdrag ...... Zonal Surface and Base Layer Stress (Pa) |
312 |
C ydrag ...... Meridional Surface and Base Layer Stress (Pa) |
313 |
C |
314 |
C*********************************************************************** |
315 |
|
316 |
implicit none |
317 |
|
318 |
c Input Variables |
319 |
c --------------- |
320 |
integer irun,lm |
321 |
_RL ps(irun) |
322 |
_RL u(irun,lm), v(irun,lm), t(irun,lm) |
323 |
_RL dudt(irun,lm), dvdt(irun,lm) |
324 |
_RL xdrag(irun), ydrag(irun) |
325 |
_RL std(irun) |
326 |
_RL ple(irun,lm+1), pl(irun,lm), dpres(irun,lm) |
327 |
_RL grav, rgas, cp |
328 |
integer nthin(irun),nbase(irun) |
329 |
_RL lstar |
330 |
|
331 |
c Dynamic Allocation Variables |
332 |
c ---------------------------- |
333 |
_RL ubar(irun), vbar(irun), robar(irun) |
334 |
_RL speed(irun), ang(irun) |
335 |
_RL bv(irun,lm) |
336 |
_RL nbar(irun) |
337 |
|
338 |
_RL tstd(irun) |
339 |
_RL XTENS(irun,lm+1), YTENS(irun,lm+1) |
340 |
_RL TENSIO(irun,lm+1) |
341 |
_RL DRAGSF(irun) |
342 |
_RL RO(irun,lm), DZ(irun,lm) |
343 |
|
344 |
integer icrilv(irun) |
345 |
|
346 |
c Local Variables |
347 |
c --------------- |
348 |
integer i,l |
349 |
_RL a,g,stdmax,agrav,akwnmb |
350 |
_RL gocp,roave,roiave,frsf,gstar,vai1,vai2 |
351 |
_RL vaisd,velco,deluu,delvv,delve2,delz,vsqua |
352 |
_RL richsn,crifro,crif2,fro2,coef |
353 |
|
354 |
c Initialization |
355 |
c -------------- |
356 |
a = 1.0 |
357 |
g = 1.0 |
358 |
agrav = 1.0/GRAV |
359 |
akwnmb = 1.0/lstar |
360 |
gocp = GRAV/CP |
361 |
|
362 |
c Constrain the Maximum Value of the Standard Deviation |
363 |
c ----------------------------------------------------- |
364 |
stdmax = 400. |
365 |
do i = 1,irun |
366 |
tstd(i) = std(i) |
367 |
if( std(i).gt.stdmax ) tstd(i) = stdmax |
368 |
enddo |
369 |
|
370 |
c Compute Atmospheric Density |
371 |
c --------------------------- |
372 |
do l = 1,lm |
373 |
do i = 1,irun |
374 |
ro(i,l) = pl(i,l)/(rgas*t(i,lm+1-l)) |
375 |
enddo |
376 |
enddo |
377 |
|
378 |
c Compute Layer Thicknesses |
379 |
c ------------------------- |
380 |
do l = 2,lm |
381 |
do i = 1,irun |
382 |
roiave = ( 1./ro(i,l-1) + 1./ro(i,l) )*0.5 |
383 |
dz(i,l) = agrav*roiave*( pl(i,l-1)-pl(i,l) ) |
384 |
enddo |
385 |
enddo |
386 |
|
387 |
|
388 |
c****************************************************** |
389 |
c Surface and Base Layer Stress * |
390 |
c****************************************************** |
391 |
|
392 |
c Definition of Surface Wind Vector |
393 |
c --------------------------------- |
394 |
do i = 1,irun |
395 |
robar(i) = 0.0 |
396 |
ubar(i) = 0.0 |
397 |
vbar(i) = 0.0 |
398 |
enddo |
399 |
|
400 |
do i = 1,irun |
401 |
do L = 1,nbase(i)-1 |
402 |
robar(i) = robar(i) + ro(i,L) *(ple(i,L)-ple(i,L+1)) |
403 |
ubar(i) = ubar(i) + u(i,lm+1-L)*(ple(i,L)-ple(i,L+1)) |
404 |
vbar(i) = vbar(i) + v(i,lm+1-L)*(ple(i,L)-ple(i,L+1)) |
405 |
enddo |
406 |
enddo |
407 |
|
408 |
do i = 1,irun |
409 |
robar(i) = robar(i)/(ple(i,1)-ple(i,nbase(i))) * 100.0 |
410 |
ubar(i) = ubar(i)/(ple(i,1)-ple(i,nbase(i))) |
411 |
vbar(i) = vbar(i)/(ple(i,1)-ple(i,nbase(i))) |
412 |
|
413 |
speed(i) = SQRT( ubar(i)*ubar(i) + vbar(i)*vbar(i) ) |
414 |
ang(i) = ATAN2(vbar(i),ubar(i)) |
415 |
|
416 |
enddo |
417 |
|
418 |
c Brunt Vaisala Frequency |
419 |
c ----------------------- |
420 |
do i = 1,irun |
421 |
do l = 2,nbase(i) |
422 |
VAI1 = (T(i,lm+1-l)-T(i,lm+2-l))/DZ(i,l)+GOCP |
423 |
if( VAI1.LT.0.0 ) then |
424 |
VAI1 = 0.0 |
425 |
endif |
426 |
VAI2 = 2.0*GRAV/( T(i,lm+1-l)+T(i,lm+2-l) ) |
427 |
VSQUA = VAI1*VAI2 |
428 |
BV(i,l) = SQRT(VSQUA) |
429 |
enddo |
430 |
enddo |
431 |
|
432 |
c Stress at the Surface Level |
433 |
c --------------------------- |
434 |
do i = 1,irun |
435 |
nbar(i) = 0.0 |
436 |
enddo |
437 |
do i = 1,irun |
438 |
do l = 2,nbase(i) |
439 |
NBAR(i) = NBAR(i) + BV(i,l)*(pl(i,l-1)-pl(i,l)) |
440 |
enddo |
441 |
enddo |
442 |
|
443 |
do i = 1,irun |
444 |
NBAR(i) = NBAR(i)/(pl(i,1)-pl(i,nbase(i))) |
445 |
FRSF = NBAR(i)*tstd(i)/speed(i) |
446 |
|
447 |
if( speed(i).eq.0.0 .or. nbar(i).eq.0.0 ) then |
448 |
TENSIO(i,1) = 0.0 |
449 |
else |
450 |
GSTAR = G*FRSF*FRSF/(FRSF*FRSF+A*A) |
451 |
TENSIO(i,1) = GSTAR*(ROBAR(i)*speed(i)*speed(i)*speed(i)) |
452 |
. / (NBAR(i)*LSTAR) |
453 |
endif |
454 |
|
455 |
XTENS(i,1) = TENSIO(i,1) * cos(ang(i)) |
456 |
YTENS(i,1) = TENSIO(i,1) * sin(ang(i)) |
457 |
DRAGSF(i) = TENSIO(i,1) |
458 |
XDRAG(i) = XTENS(i,1) |
459 |
YDRAG(i) = YTENS(i,1) |
460 |
enddo |
461 |
|
462 |
c Check for Very thin lowest layer |
463 |
c -------------------------------- |
464 |
do i = 1,irun |
465 |
if( NTHIN(i).gt.1 ) then |
466 |
do l = 1,nthin(i) |
467 |
TENSIO(i,l) = TENSIO(i,1) |
468 |
XTENS(i,l) = XTENS(i,1) |
469 |
YTENS(i,l) = YTENS(i,1) |
470 |
enddo |
471 |
endif |
472 |
enddo |
473 |
|
474 |
c****************************************************** |
475 |
c Compute Gravity Wave Stress from NTHIN+1 to NBASE * |
476 |
c****************************************************** |
477 |
|
478 |
do i = 1,irun |
479 |
do l = nthin(i)+1,nbase(i) |
480 |
|
481 |
velco = 0.5*( (u(i,lm+1-l)*ubar(i) + v(i,lm+1-l)*vbar(i)) |
482 |
. + (u(i,lm+2-l)*ubar(i) + v(i,lm+2-l)*vbar(i)) ) |
483 |
. / speed(i) |
484 |
|
485 |
C Convert to Newton/m**2 |
486 |
roave = 0.5*(ro(i,l-1)+ro(i,l)) * 100.0 |
487 |
|
488 |
if( VELCO.le.0.0 ) then |
489 |
TENSIO(i,l) = TENSIO(i,l-1) |
490 |
goto 1500 |
491 |
endif |
492 |
|
493 |
c Froude number squared |
494 |
c --------------------- |
495 |
FRO2 = bv(i,l)/(AKWNMB*ROAVE*VELCO*VELCO*VELCO)*TENSIO(i,l-1) |
496 |
DELUU = u(i,lm+1-l)-u(i,lm+2-l) |
497 |
DELVV = v(i,lm+1-l)-v(i,lm+2-l) |
498 |
DELVE2 = ( DELUU*DELUU + DELVV*DELVV ) |
499 |
|
500 |
c Compute Richarson Number |
501 |
c ------------------------ |
502 |
if( DELVE2.ne.0.0 ) then |
503 |
DELZ = DZ(i,l) |
504 |
VSQUA = BV(i,l)*BV(i,l) |
505 |
RICHSN = DELZ*DELZ*VSQUA/DELVE2 |
506 |
else |
507 |
RICHSN = 99999.0 |
508 |
endif |
509 |
|
510 |
if( RICHSN.le.0.25 ) then |
511 |
TENSIO(i,l) = TENSIO(i,l-1) |
512 |
goto 1500 |
513 |
endif |
514 |
|
515 |
c Stress in the Base Layer changes if the local Froude number |
516 |
c exceeds the Critical Froude number |
517 |
c ---------------------------------- |
518 |
CRIFRO = 1.0 - 0.25/RICHSN |
519 |
CRIF2 = CRIFRO*CRIFRO |
520 |
if( l.eq.2 ) CRIF2 = MIN(0.7,CRIF2) |
521 |
|
522 |
if( FRO2.gt.CRIF2 ) then |
523 |
TENSIO(i,l) = CRIF2/FRO2*TENSIO(i,l-1) |
524 |
else |
525 |
TENSIO(i,l) = TENSIO(i,l-1) |
526 |
endif |
527 |
|
528 |
1500 CONTINUE |
529 |
XTENS(i,l) = TENSIO(i,l)*COS(ang(i)) |
530 |
YTENS(i,l) = TENSIO(i,l)*SIN(ang(i)) |
531 |
|
532 |
enddo |
533 |
enddo |
534 |
|
535 |
c****************************************************** |
536 |
c Compute Gravity Wave Stress from Base+1 to Top * |
537 |
c****************************************************** |
538 |
|
539 |
do i = 1,irun |
540 |
icrilv(i) = 0 |
541 |
enddo |
542 |
|
543 |
do i = 1,irun |
544 |
do l = nbase(i)+1,lm+1 |
545 |
|
546 |
TENSIO(i,l) = 0.0 |
547 |
|
548 |
c Check for Critical Level Absorption |
549 |
c ----------------------------------- |
550 |
if( icrilv(i).eq.1 ) goto 130 |
551 |
|
552 |
c Let Remaining Stress escape out the top edge of model |
553 |
c ----------------------------------------------------- |
554 |
if( l.eq.lm+1 ) then |
555 |
TENSIO(i,l) = TENSIO(i,l-1) |
556 |
goto 130 |
557 |
endif |
558 |
|
559 |
ROAVE = 0.5*(ro(i,l-1)+ro(i,l)) * 100.0 |
560 |
VAI1 = (T(i,lm+1-l)-T(i,lm+2-l))/DZ(i,l)+GOCP |
561 |
|
562 |
if( VAI1.lt.0.0 ) then |
563 |
icrilv(i) = 1 |
564 |
TENSIO(i,l) = 0.0 |
565 |
goto 130 |
566 |
endif |
567 |
|
568 |
VAI2 = 2.0*GRAV/(T(i,lm+1-l)+T(i,lm+2-l)) |
569 |
VSQUA = VAI1*VAI2 |
570 |
VAISD = SQRT(VSQUA) |
571 |
|
572 |
velco = 0.5*( (u(i,lm+1-l)*ubar(i) + v(i,lm+1-l)*vbar(i)) |
573 |
. + (u(i,lm+2-l)*ubar(i) + v(i,lm+2-l)*vbar(i)) ) |
574 |
. / speed(i) |
575 |
|
576 |
if( velco.lt.0.0 ) then |
577 |
icrilv(i) = 1 |
578 |
TENSIO(i,l) = 0.0 |
579 |
goto 130 |
580 |
endif |
581 |
|
582 |
c Froude number squared |
583 |
c --------------------- |
584 |
FRO2 = vaisd/(AKWNMB*ROAVE*VELCO*VELCO*VELCO)*TENSIO(i,l-1) |
585 |
DELUU = u(i,lm+1-l)-u(i,lm+2-l) |
586 |
DELVV = v(i,lm+1-l)-v(i,lm+2-l) |
587 |
DELVE2 = ( DELUU*DELUU + DELVV*DELVV ) |
588 |
|
589 |
c Compute Richarson Number |
590 |
c ------------------------ |
591 |
if( DELVE2.ne.0.0 ) then |
592 |
DELZ = DZ(i,l) |
593 |
RICHSN = DELZ*DELZ*VSQUA/DELVE2 |
594 |
else |
595 |
RICHSN = 99999.0 |
596 |
endif |
597 |
|
598 |
if( RICHSN.le.0.25 ) then |
599 |
TENSIO(i,l) = 0.0 |
600 |
icrilv(i) = 1 |
601 |
goto 130 |
602 |
endif |
603 |
|
604 |
c Stress in Layer changes if the local Froude number |
605 |
c exceeds the Critical Froude number |
606 |
c ---------------------------------- |
607 |
CRIFRO = 1.0 - 0.25/RICHSN |
608 |
CRIF2 = CRIFRO*CRIFRO |
609 |
|
610 |
if( FRO2.ge.CRIF2 ) then |
611 |
TENSIO(i,l) = CRIF2/FRO2*TENSIO(i,l-1) |
612 |
else |
613 |
TENSIO(i,l) = TENSIO(i,l-1) |
614 |
endif |
615 |
|
616 |
130 continue |
617 |
XTENS(i,l) = TENSIO(i,l)*COS(ang(i)) |
618 |
YTENS(i,l) = TENSIO(i,l)*SIN(ang(i)) |
619 |
enddo |
620 |
enddo |
621 |
|
622 |
C ****************************************************** |
623 |
C MOMENTUM CHANGE FOR FREE ATMOSPHERE * |
624 |
C ****************************************************** |
625 |
|
626 |
do i = 1,irun |
627 |
do l = nthin(i)+1,lm |
628 |
coef = -grav*ple(i,lm+1)/dpres(i,lm+1-l) |
629 |
dudt(i,lm+1-l) = coef*(XTENS(i,l+1)-XTENS(i,l)) |
630 |
dvdt(i,lm+1-l) = coef*(YTENS(i,l+1)-YTENS(i,l)) |
631 |
enddo |
632 |
enddo |
633 |
|
634 |
c Momentum change near the surface |
635 |
c -------------------------------- |
636 |
do i = 1,irun |
637 |
coef = grav*ple(i,lm+1)/(ple(i,lm+1-nthin(i))-ple(i,lm+1)) |
638 |
dudt(i,lm) = coef*(XTENS(i,nthin(i)+1)-XTENS(i,1)) |
639 |
dvdt(i,lm) = coef*(YTENS(i,nthin(i)+1)-YTENS(i,1)) |
640 |
enddo |
641 |
|
642 |
c If Lowest layer is very thin, it is strapped to next layer |
643 |
c ---------------------------------------------------------- |
644 |
do i = 1,irun |
645 |
if( nthin(i).gt.1 ) then |
646 |
do l = 2,nthin(i) |
647 |
dudt(i,lm+1-l) = dudt(i,lm) |
648 |
dvdt(i,lm+1-l) = dvdt(i,lm) |
649 |
enddo |
650 |
endif |
651 |
enddo |
652 |
|
653 |
c Convert Units to (m/sec**2) |
654 |
c --------------------------- |
655 |
do l = 1,lm |
656 |
do i = 1,irun |
657 |
dudt(i,l) = - dudt(i,l)/ps(i)*0.01 |
658 |
dvdt(i,l) = - dvdt(i,l)/ps(i)*0.01 |
659 |
enddo |
660 |
enddo |
661 |
|
662 |
return |
663 |
end |