/[MITgcm]/MITgcm/pkg/fizhi/update_chemistry_exports.F
ViewVC logotype

Contents of /MITgcm/pkg/fizhi/update_chemistry_exports.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.10 - (show annotations) (download)
Wed Jul 28 01:25:07 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint55c_post, checkpoint54e_post, checkpoint55d_pre, checkpoint55h_post, checkpoint55b_post, checkpoint55, checkpoint54f_post, checkpoint55g_post, checkpoint55f_post, checkpoint55e_post, checkpoint55a_post, checkpoint55d_post
Changes since 1.9: +15 -15 lines
debugging

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

  ViewVC Help
Powered by ViewVC 1.1.22