1 |
C $Header: /u/gcmpack/MITgcm/pkg/fizhi/update_chemistry_exports.F,v 1.11 2004/10/22 14:52:14 molod Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "FIZHI_OPTIONS.h" |
5 |
subroutine update_chemistry_exports (myTime, myIter, myThid) |
6 |
c---------------------------------------------------------------------- |
7 |
c Subroutine update_chemistry_exports - 'Wrapper' routine to update |
8 |
c the fields related to the earth's chemistry that are needed |
9 |
c by fizhi. |
10 |
c Also: Set up "bi, bj loop" and some timers and clocks here. |
11 |
c |
12 |
c Call: interp_chemistry |
13 |
c----------------------------------------------------------------------- |
14 |
implicit none |
15 |
#include "SIZE.h" |
16 |
#include "fizhi_SIZE.h" |
17 |
#include "fizhi_land_SIZE.h" |
18 |
#include "GRID.h" |
19 |
#include "DYNVARS.h" |
20 |
#include "fizhi_chemistry_coms.h" |
21 |
#include "fizhi_coms.h" |
22 |
#include "gridalt_mapping.h" |
23 |
#include "EEPARAMS.h" |
24 |
#include "chronos.h" |
25 |
|
26 |
integer myIter, myThid |
27 |
_RL myTime |
28 |
|
29 |
c pe on physics grid refers to bottom edge |
30 |
_RL pephy(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nrphys+1,nSx,nSy) |
31 |
_RL pphy(sNx,sNy,Nrphys,nSx,nSy) |
32 |
_RL oz1(nlatsoz,nlevsoz), strq1(nlatsq,nlevsq) |
33 |
_RL waterin(sNx,sNy,Nrphys), xlat(sNx,sNy) |
34 |
|
35 |
integer i, j, L, LL, bi, bj |
36 |
integer im1, im2, jm1, jm2 |
37 |
integer nhms1,nymd1,nhms2,nymd2,imns,ipls |
38 |
_RL facm, facp |
39 |
logical alarm |
40 |
external alarm |
41 |
|
42 |
im1 = 1 |
43 |
im2 = sNx |
44 |
jm1 = 1 |
45 |
jm2 = sNy |
46 |
|
47 |
if( alarm('radsw').or.alarm('radlw') ) then |
48 |
|
49 |
do bj = myByLo(myThid), myByHi(myThid) |
50 |
do bi = myBxLo(myThid), myBxHi(myThid) |
51 |
|
52 |
c Construct the physics grid pressures - count pephy levels top down |
53 |
c (even though dpphy counted bottom up) |
54 |
do j = 1,sNy |
55 |
do i = 1,sNx |
56 |
pephy(i,j,Nrphys+1,bi,bj)=(Ro_surf(i,j,bi,bj)+etaH(i,j,bi,bj)) |
57 |
do L = 2,Nrphys+1 |
58 |
LL = Nrphys+2-L |
59 |
pephy(i,j,LL,bi,bj)=pephy(i,j,LL+1,bi,bj)-dpphys(i,j,L-1,bi,bj) |
60 |
enddo |
61 |
enddo |
62 |
enddo |
63 |
do j = 1,sNy |
64 |
do i = 1,sNx |
65 |
do L = 1,Nrphys |
66 |
pphy(i,j,L,bi,bj)=(pephy(i,j,L+1,bi,bj)+pephy(i,j,L,bi,bj)) |
67 |
. /200. |
68 |
enddo |
69 |
enddo |
70 |
enddo |
71 |
|
72 |
do j = 1,sNy |
73 |
do i = 1,sNx |
74 |
xlat(i,j) = yC(i,j,bi,bj) |
75 |
do L = 1,Nrphys |
76 |
waterin(i,j,L) = sphy(i,j,L,bi,bj) |
77 |
enddo |
78 |
enddo |
79 |
enddo |
80 |
|
81 |
call time_bound(nymd,nhms,nymd1,nhms1,nymd2,nhms2,imns,ipls) |
82 |
call interp_time(nymd,nhms,nymd1,nhms1,nymd2,nhms2,facm,facp) |
83 |
|
84 |
do L = 1,nlevsoz |
85 |
do j = 1,nlatsoz |
86 |
oz1(j,L) = ozone(j,L,imns)*facm + ozone(j,L,ipls)*facp |
87 |
enddo |
88 |
enddo |
89 |
|
90 |
do L = 1,nlevsq |
91 |
do j = 1,nlatsq |
92 |
strq1(j,L) = stratq(j,L,imns)*facm + stratq(j,L,ipls)*facp |
93 |
enddo |
94 |
enddo |
95 |
|
96 |
call interp_chemistry(strq1,nlevsq,nlatsq,levsq,latsq, |
97 |
. oz1,nlevsoz,nlatsoz,levsoz,latsoz, |
98 |
. waterin,pphy(1,1,1,bi,bj),xlat, |
99 |
. im2,jm2,Nrphys,nSx,nSy,bi,bj,o3,qstr) |
100 |
|
101 |
enddo |
102 |
enddo |
103 |
|
104 |
endif |
105 |
|
106 |
return |
107 |
end |
108 |
|
109 |
subroutine interp_chemistry (stratq,nwatlevs,nwatlats,watlevs, |
110 |
. watlats,ozone,nozlevs,nozlats,ozlevs,ozlats, |
111 |
. qz,plz,xlat,im,jm,lm,nSx,nSy,bi,bj,ozrad,qzrad) |
112 |
|
113 |
implicit none |
114 |
|
115 |
c Input Variables |
116 |
c --------------- |
117 |
integer nwatlevs,nwatlats,nozlevs,nozlats,nSx,nSy,bi,bj |
118 |
_RL stratq(nwatlats,nwatlevs),ozone(nozlats,nozlevs) |
119 |
_RL watlevs(nwatlevs),watlats(nwatlats) |
120 |
_RL ozlevs(nozlevs),ozlats(nozlats) |
121 |
integer im,jm,lm |
122 |
_RL qz(im,jm,lm),plz(im,jm,lm) |
123 |
_RL xlat(im,jm) |
124 |
_RL ozrad(im,jm,lm,nSx,nSy) |
125 |
_RL qzrad(im,jm,lm,nSx,nSy) |
126 |
|
127 |
C ********************************************************************** |
128 |
C **** Get Ozone and Stratospheric Moisture Data **** |
129 |
C ********************************************************************** |
130 |
|
131 |
call interp_qz (stratq,nwatlevs,nwatlats,watlevs,watlats,im*jm, |
132 |
. bi,bj, xlat,lm,plz,qz,qzrad(1,1,1,bi,bj)) |
133 |
call interp_oz (ozone ,nozlevs,nozlats,ozlevs,ozlats,im*jm, |
134 |
. bi,bj, xlat,lm,plz,ozrad(1,1,1,bi,bj)) |
135 |
|
136 |
return |
137 |
end |
138 |
|
139 |
subroutine interp_qz(stratq,nwatlevs,nwatlats,watlevs,watlats, |
140 |
. irun,bi,bj,xlat,nlevs,pres,qz_in,qz_out ) |
141 |
C*********************************************************************** |
142 |
C Purpose |
143 |
C To Interpolate Chemistry Moisture from Chemistry Grid to Physics Grid |
144 |
C |
145 |
C INPUT Argument Description |
146 |
C stratq .... Climatological SAGE Stratospheric Moisture |
147 |
C irun ...... Number of Columns to be filled |
148 |
C xlat ...... Latitude in Degrees |
149 |
C nlevs ..... Vertical Dimension |
150 |
C pres ...... PRES (IM,JM,nlevs) Three-dimensional array of pressures |
151 |
C qz_in ..... Model Moisture (kg/kg mass mixing radtio) |
152 |
C qz_out .... Combination of Chemistry Moisture and Model Moisture |
153 |
C (kg/kg mass mixing ratio) |
154 |
C |
155 |
C*********************************************************************** |
156 |
|
157 |
implicit none |
158 |
integer nwatlevs,nwatlats |
159 |
integer bi,bj |
160 |
_RL stratq ( nwatlats,nwatlevs ) |
161 |
_RL watlats (nwatlats) |
162 |
_RL watlevs (nwatlevs) |
163 |
|
164 |
integer irun,nlevs |
165 |
_RL xlat (irun) |
166 |
_RL pres (irun,nlevs) |
167 |
_RL qz_in (irun,nlevs) |
168 |
_RL qz_out(irun,nlevs) |
169 |
|
170 |
c Local Variables |
171 |
c --------------- |
172 |
integer pqu,pql,dpq |
173 |
parameter ( pqu = 100. ) |
174 |
parameter ( pql = 300. ) |
175 |
parameter ( dpq = pql-pqu ) |
176 |
|
177 |
integer i,k,L1,L2,LM,LP |
178 |
_RL h2o_time_lat (irun,nwatlevs) |
179 |
_RL qz_clim(irun,nlevs) |
180 |
|
181 |
_RL qpr1(irun), qpr2(irun), slope(irun) |
182 |
_RL pr1(irun), pr2(irun) |
183 |
|
184 |
integer jlat,jlatm,jlatp |
185 |
|
186 |
C ********************************************************************** |
187 |
C **** Interpolate Moisture data to model latitudes *** |
188 |
C ********************************************************************** |
189 |
|
190 |
DO 32 k = 1, nwatlevs |
191 |
DO 34 i = 1,irun |
192 |
|
193 |
DO 36 jlat = 1, nwatlats |
194 |
IF( watlats(jlat).gt.xlat(i) ) THEN |
195 |
IF( jlat.EQ.1 ) THEN |
196 |
jlatm = 1 |
197 |
jlatp = 1 |
198 |
slope(i) = 0 |
199 |
ELSE |
200 |
jlatm = jlat -1 |
201 |
jlatp = jlat |
202 |
slope(i) = ( xlat(i) -watlats(jlat-1) ) |
203 |
. / ( watlats(jlat)-watlats(jlat-1) ) |
204 |
ENDIF |
205 |
GOTO 37 |
206 |
ENDIF |
207 |
36 CONTINUE |
208 |
jlatm = nwatlats |
209 |
jlatp = nwatlats |
210 |
slope(i) = 1 |
211 |
37 CONTINUE |
212 |
QPR1(i) = stratq(jlatm,k) |
213 |
QPR2(i) = stratq(jlatp,k) |
214 |
34 CONTINUE |
215 |
|
216 |
do i = 1,irun |
217 |
h2o_time_lat(i,k) = qpr1(i) + slope(i)*(qpr2(i)-qpr1(i)) |
218 |
enddo |
219 |
|
220 |
32 CONTINUE |
221 |
|
222 |
C ********************************************************************** |
223 |
C **** Interpolate Latitude Moisture data to model pressures *** |
224 |
C ********************************************************************** |
225 |
|
226 |
DO 40 L2 = 1,nlevs |
227 |
|
228 |
DO 44 i= 1, irun |
229 |
DO 46 L1 = 1,nwatlevs |
230 |
IF( watlevs(L1).GT.pres(i,L2) ) THEN |
231 |
IF( L1.EQ.1 ) THEN |
232 |
LM = 1 |
233 |
LP = 2 |
234 |
ELSE |
235 |
LM = L1-1 |
236 |
LP = L1 |
237 |
ENDIF |
238 |
GOTO 47 |
239 |
ENDIF |
240 |
46 CONTINUE |
241 |
LM = nwatlevs-1 |
242 |
LP = nwatlevs |
243 |
47 CONTINUE |
244 |
PR1(i) = watlevs (LM) |
245 |
PR2(i) = watlevs (LP) |
246 |
QPR1(i) = h2o_time_lat(i,LM) |
247 |
QPR2(i) = h2o_time_lat(i,LP) |
248 |
44 CONTINUE |
249 |
|
250 |
do i= 1, irun |
251 |
slope(i) =(QPR1(i)-QPR2(i)) / (PR1(i)-PR2(i)) |
252 |
qz_clim(i,L2) = QPR2(i) + (pres(i,L2)-PR2(i))*SLOPE(i) |
253 |
enddo |
254 |
|
255 |
40 CONTINUE |
256 |
|
257 |
c |
258 |
c ... Above 100 mb, using climatological water data set ................... |
259 |
c ... Below 300 mb, using model predicted water data set ................... |
260 |
c ... In between, using linear interpolation ............................... |
261 |
c |
262 |
do k= 1, nlevs |
263 |
do i= 1, irun |
264 |
if( pres(i,k).ge.pqu .and. pres(i,k).le. pql) then |
265 |
qz_out(i,k) = qz_clim(i,k)+(qz_in(i,k)- |
266 |
1 qz_clim(i,k))*(pres(i,k)-pqu)/dpq |
267 |
else if( pres(i,k) .gt. pql ) then |
268 |
qz_out(i,k) = qz_in (i,k) |
269 |
else |
270 |
qz_out(i,k) = qz_clim(i,k) |
271 |
endif |
272 |
enddo |
273 |
enddo |
274 |
|
275 |
return |
276 |
end |
277 |
|
278 |
subroutine interp_oz (ozone,nozlevs,nozlats,ozlevs,ozlats, |
279 |
. irun,bi,bj,xlat,nlevs,plevs,ozrad) |
280 |
C*********************************************************************** |
281 |
C Purpose |
282 |
C To Interpolate Chemistry Ozone from Chemistry Grid to Physics Grid |
283 |
C |
284 |
C INPUT Argument Description |
285 |
C ozone ..... Climatological Ozone |
286 |
C chemistry .. Chemistry State Data Structure |
287 |
C irun ....... Number of Columns to be filled |
288 |
C xlat ....... Latitude in Degrees |
289 |
C nlevs ...... Vertical Dimension |
290 |
C pres ....... Three-dimensional array of pressures |
291 |
C ozrad ...... Ozone on Physics Grid (kg/kg mass mixing radtio) |
292 |
C |
293 |
C*********************************************************************** |
294 |
implicit none |
295 |
integer nozlevs,nozlats,irun,nlevs |
296 |
integer bi,bj |
297 |
_RL ozone(nozlats,nozlevs) |
298 |
_RL xlat(irun) |
299 |
_RL plevs(irun,nlevs) |
300 |
_RL ozrad(irun,nlevs) |
301 |
_RL ozlevs(nozlevs),ozlats(nozlats) |
302 |
|
303 |
c Local Variables |
304 |
c --------------- |
305 |
_RL zero,one,o3min,voltomas |
306 |
PARAMETER ( ZERO = 0.0 ) |
307 |
PARAMETER ( ONE = 1.0 ) |
308 |
PARAMETER ( O3MIN = 1.0E-10 ) |
309 |
PARAMETER ( VOLTOMAS = 1.655E-6 ) |
310 |
|
311 |
integer i,k,L1,L2,LM,LP |
312 |
integer jlat,jlatm,jlatp |
313 |
_RL O3INT1(IRUN,nozlevs) |
314 |
_RL QPR1(IRUN), QPR2(IRUN), SLOPE(IRUN) |
315 |
_RL PR1(IRUN), PR2(IRUN) |
316 |
|
317 |
C ********************************************************************** |
318 |
C **** INTERPOLATE ozone data to model latitudes *** |
319 |
C ********************************************************************** |
320 |
|
321 |
DO 32 K=1,nozlevs |
322 |
DO 34 I=1,IRUN |
323 |
|
324 |
DO 36 jlat = 1,nozlats |
325 |
IF( ozlats(jlat).gt.xlat(i) ) THEN |
326 |
IF( jlat.EQ.1 ) THEN |
327 |
jlatm = 1 |
328 |
jlatp = 1 |
329 |
slope(i) = zero |
330 |
ELSE |
331 |
jlatm = jlat-1 |
332 |
jlatp = jlat |
333 |
slope(i) = ( XLAT(I) -ozlats(jlat-1) ) |
334 |
. / ( ozlats(jlat)-ozlats(jlat-1) ) |
335 |
ENDIF |
336 |
GOTO 37 |
337 |
ENDIF |
338 |
36 CONTINUE |
339 |
jlatm = nozlats |
340 |
jlatp = nozlats |
341 |
slope(i) = one |
342 |
37 CONTINUE |
343 |
QPR1(I) = ozone(jlatm,k) |
344 |
QPR2(I) = ozone(jlatp,k) |
345 |
34 CONTINUE |
346 |
|
347 |
DO 38 I=1,IRUN |
348 |
o3int1(i,k) = qpr1(i) + slope(i)*( qpr2(i)-qpr1(i) ) |
349 |
38 CONTINUE |
350 |
|
351 |
32 CONTINUE |
352 |
|
353 |
C ********************************************************************** |
354 |
C **** INTERPOLATE latitude ozone data to model pressures *** |
355 |
C ********************************************************************** |
356 |
|
357 |
DO 40 L2 = 1,NLEVS |
358 |
|
359 |
DO 44 I = 1,IRUN |
360 |
DO 46 L1 = 1,nozlevs |
361 |
IF( ozlevs(L1).GT.PLEVS(I,L2) ) THEN |
362 |
IF( L1.EQ.1 ) THEN |
363 |
LM = 1 |
364 |
LP = 2 |
365 |
ELSE |
366 |
LM = L1-1 |
367 |
LP = L1 |
368 |
ENDIF |
369 |
GOTO 47 |
370 |
ENDIF |
371 |
46 CONTINUE |
372 |
LM = nozlevs-1 |
373 |
LP = nozlevs |
374 |
47 CONTINUE |
375 |
PR1(I) = ozlevs (LM) |
376 |
PR2(I) = ozlevs (LP) |
377 |
QPR1(I) = O3INT1(I,LM) |
378 |
QPR2(I) = O3INT1(I,LP) |
379 |
44 CONTINUE |
380 |
|
381 |
DO 48 I=1,IRUN |
382 |
SLOPE(I) = ( QPR1(I)-QPR2(I) ) |
383 |
. / ( PR1(I)- PR2(I) ) |
384 |
ozrad(I,L2) = QPR2(I) + ( PLEVS(I,L2)-PR2(I) )*SLOPE(I) |
385 |
|
386 |
if( ozrad(i,l2).lt.o3min ) then |
387 |
ozrad(i,l2) = o3min |
388 |
endif |
389 |
|
390 |
48 CONTINUE |
391 |
40 CONTINUE |
392 |
|
393 |
C ********************************************************************** |
394 |
C **** CONVERT FROM VOLUME MIXING RATIO TO MASS MIXING RATIO *** |
395 |
C ********************************************************************** |
396 |
|
397 |
DO 60 I=1,IRUN*NLEVS |
398 |
ozrad (I,1) = ozrad(I,1) * VOLTOMAS |
399 |
60 CONTINUE |
400 |
|
401 |
RETURN |
402 |
END |
403 |
|