/[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.8 - (show annotations) (download)
Fri Jul 16 19:37:04 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint54c_post
Changes since 1.7: +46 -55 lines
Debugging
Add Land Surface Model (Koster-Suarez) code to fizhi

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

  ViewVC Help
Powered by ViewVC 1.1.22