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

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

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


Revision 1.9 - (hide annotations) (download)
Mon Jul 26 18:45:17 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint54d_post
Changes since 1.8: +25 -28 lines
Went to use of FIZHI_OPTIONS and _RL in all routines

1 molod 1.9 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/update_chemistry_exports.F,v 1.8 2004/07/16 19:37:04 molod Exp $
2 molod 1.1 C $Name: $
3    
4 molod 1.9 #include "FIZHI_OPTIONS.h"
5 molod 1.1 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 molod 1.8 #include "fizhi_land_SIZE.h"
18 molod 1.1 #include "GRID.h"
19     #include "DYNVARS.h"
20 molod 1.2 #include "fizhi_chemistry_coms.h"
21 molod 1.4 #include "fizhi_coms.h"
22 molod 1.1 #include "gridalt_mapping.h"
23     #include "EEPARAMS.h"
24 molod 1.4 #include "chronos.h"
25 molod 1.1
26 molod 1.4 integer myTime, myIter, myThid
27 molod 1.1
28     c pe on physics grid refers to bottom edge
29 molod 1.9 _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 molod 1.4
34 molod 1.5 integer i, j, L, LL, bi, bj
35 molod 1.4 integer im1, im2, jm1, jm2, idim1, idim2, jdim1, jdim2
36     integer nhms1,nymd1,nhms2,nymd2,imns,ipls
37 molod 1.9 _RL facm, facp
38 molod 1.6 logical alarm
39     external alarm
40 molod 1.4
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 molod 1.1
52     do bj = myByLo(myThid), myByHi(myThid)
53     do bi = myBxLo(myThid), myBxHi(myThid)
54    
55 molod 1.5 c Construct the physics grid pressures - count pephy levels top down
56     c (even though dpphy counted bottom up)
57 molod 1.1 do j = 1,sNy
58     do i = 1,sNx
59 molod 1.7 pephy(i,j,Nrphys+1,bi,bj)=(Ro_surf(i,j,bi,bj)+etaH(i,j,bi,bj))
60 molod 1.1 do L = 2,Nrphys+1
61 molod 1.5 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 molod 1.1 enddo
64     enddo
65     enddo
66 molod 1.4 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 molod 1.1
74 molod 1.4 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 molod 1.1
102     enddo
103     enddo
104    
105 molod 1.4 endif
106    
107     return
108     end
109 molod 1.3
110     subroutine interp_chemistry (stratq,nwatlevs,nwatlats,watlevs,
111 molod 1.8 . watlats,ozone,nozlevs,nozlats,ozlevs,ozlats,
112     . qz,plz,xlat,im,jm,lm,nSx,nSy,bi,bj,ozrad,qzrad)
113 molod 1.3
114     implicit none
115    
116     c Input Variables
117     c ---------------
118 molod 1.8 integer nwatlevs,nwatlats,nozlevs,nozlats,nSx,nSy,bi,bj
119 molod 1.9 _RL stratq(nwatlats,nwatlevs),ozone(nozlats,nozlevs)
120 molod 1.8 _RL watlevs(nwatlevs),watlats(nwatlats)
121     _RL ozlevs(nozlevs),ozlats(nozlats)
122     integer im,jm,lm
123 molod 1.9 _RL qz(im,jm,lm),plz(im,jm,lm)
124     _RL xlat(im,jm)
125 molod 1.8 _RL ozrad(im,jm,lm,nSx,nSy)
126     _RL qzrad(im,jm,lm,nSx,nSy)
127 molod 1.3
128     C **********************************************************************
129     C **** Get Ozone and Stratospheric Moisture Data ****
130     C **********************************************************************
131    
132     call interp_qz (stratq,nwatlevs,nwatlats,watlevs,watlats,im*jm,
133 molod 1.4 . xlat,lm,plz,qz,qzrad(1,1,1,bi,bj))
134 molod 1.8 call interp_oz (ozone ,nozlevs,nozlats,ozlevs,ozlats,im*jm,
135 molod 1.4 . xlat,lm,plz,ozrad(1,1,1,bi,bj))
136 molod 1.3 return
137     end
138    
139     subroutine interp_qz(stratq,nwatlevs,nwatlats,watlevs,watlats,
140     . irun,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 molod 1.8 C qz_out .... Combination of Chemistry Moisture and Model Moisture
153     C (kg/kg mass mixing ratio)
154 molod 1.3 C
155     C***********************************************************************
156    
157     implicit none
158     integer nwatlevs,nwatlats
159 molod 1.9 _RL stratq ( nwatlats,nwatlevs )
160 molod 1.8 _RL watlats (nwatlats)
161     _RL watlevs (nwatlevs)
162 molod 1.3
163     integer irun,nlevs
164 molod 1.9 _RL xlat (irun)
165     _RL pres (irun,nlevs)
166     _RL qz_in (irun,nlevs)
167 molod 1.8 _RL qz_out(irun,nlevs)
168 molod 1.3
169     c Local Variables
170     c ---------------
171     integer pqu,pql,dpq
172     parameter ( pqu = 100. )
173     parameter ( pql = 300. )
174     parameter ( dpq = pql-pqu )
175    
176     integer i,k,L1,L2,LM,LP
177 molod 1.9 _RL h2o_time_lat (irun,nwatlevs)
178     _RL qz_clim(irun,nlevs)
179 molod 1.3
180 molod 1.9 _RL qpr1(irun), qpr2(irun), slope(irun)
181     _RL pr1(irun), pr2(irun)
182 molod 1.3
183     integer jlat,jlatm,jlatp
184    
185     C **********************************************************************
186     C **** Interpolate Moisture data to model latitudes ***
187     C **********************************************************************
188    
189     DO 32 k = 1, nwatlevs
190     DO 34 i = 1,irun
191    
192     DO 36 jlat = 1, nwatlats
193     IF( watlats(jlat).gt.xlat(i) ) THEN
194     IF( jlat.EQ.1 ) THEN
195     jlatm = 1
196     jlatp = 1
197     slope(i) = 0
198     ELSE
199     jlatm = jlat -1
200     jlatp = jlat
201     slope(i) = ( xlat(i) -watlats(jlat-1) )
202     . / ( watlats(jlat)-watlats(jlat-1) )
203     ENDIF
204     GOTO 37
205     ENDIF
206     36 CONTINUE
207     jlatm = nwatlats
208     jlatp = nwatlats
209     slope(i) = 1
210     37 CONTINUE
211     QPR1(i) = stratq(jlatm,k)
212     QPR2(i) = stratq(jlatp,k)
213     34 CONTINUE
214    
215     do i = 1,irun
216     h2o_time_lat(i,k) = qpr1(i) + slope(i)*(qpr2(i)-qpr1(i))
217     enddo
218    
219     32 CONTINUE
220    
221     C **********************************************************************
222     C **** Interpolate Latitude Moisture data to model pressures ***
223     C **********************************************************************
224    
225     DO 40 L2 = 1,nlevs
226    
227     DO 44 i= 1, irun
228     DO 46 L1 = 1,nwatlevs
229     IF( watlevs(L1).GT.pres(i,L2) ) THEN
230     IF( L1.EQ.1 ) THEN
231     LM = 1
232     LP = 2
233     ELSE
234     LM = L1-1
235     LP = L1
236     ENDIF
237     GOTO 47
238     ENDIF
239     46 CONTINUE
240     LM = nwatlevs-1
241     LP = nwatlevs
242     47 CONTINUE
243     PR1(i) = watlevs (LM)
244     PR2(i) = watlevs (LP)
245     QPR1(i) = h2o_time_lat(i,LM)
246     QPR2(i) = h2o_time_lat(i,LP)
247     44 CONTINUE
248    
249     do i= 1, irun
250     slope(i) =(QPR1(i)-QPR2(i)) / (PR1(i)-PR2(i))
251     qz_clim(i,L2) = QPR2(i) + (pres(i,L2)-PR2(i))*SLOPE(i)
252     enddo
253    
254     40 CONTINUE
255    
256     c
257     c ... Above 100 mb, using climatological water data set ...................
258     c ... Below 300 mb, using model predicted water data set ...................
259     c ... In between, using linear interpolation ...............................
260     c
261     do k= 1, nlevs
262     do i= 1, irun
263     if( pres(i,k).ge.pqu .and. pres(i,k).le. pql) then
264     qz_out(i,k) = qz_clim(i,k)+(qz_in(i,k)-
265     1 qz_clim(i,k))*(pres(i,k)-pqu)/dpq
266     else if( pres(i,k) .gt. pql ) then
267     qz_out(i,k) = qz_in (i,k)
268     else
269     qz_out(i,k) = qz_clim(i,k)
270     endif
271     enddo
272     enddo
273    
274     return
275     end
276    
277 molod 1.8 subroutine interp_oz (ozone,nozlevs,nozlats,ozlevs,ozlats,
278 molod 1.3 . irun,xlat,nlevs,plevs,ozrad)
279     C***********************************************************************
280     C Purpose
281     C To Interpolate Chemistry Ozone from Chemistry Grid to Physics Grid
282     C
283     C INPUT Argument Description
284     C ozone ..... Climatological Ozone
285     C chemistry .. Chemistry State Data Structure
286     C irun ....... Number of Columns to be filled
287     C xlat ....... Latitude in Degrees
288     C nlevs ...... Vertical Dimension
289     C pres ....... Three-dimensional array of pressures
290     C ozrad ...... Ozone on Physics Grid (kg/kg mass mixing radtio)
291     C
292     C***********************************************************************
293     implicit none
294 molod 1.8 integer nozlevs,nozlats,irun,nlevs
295 molod 1.9 _RL ozone(nozlats,nozlevs)
296     _RL xlat(irun)
297     _RL plevs(irun,nlevs)
298 molod 1.8 _RL ozrad(irun,nlevs)
299     _RL ozlevs(nozlevs),ozlats(nozlats)
300 molod 1.3
301     c Local Variables
302     c ---------------
303 molod 1.9 _RL zero,one,o3min,voltomas
304 molod 1.3 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 molod 1.9 _RL O3INT1(IRUN,nozlevs)
312     _RL QPR1(IRUN), QPR2(IRUN), SLOPE(IRUN)
313     _RL PR1(IRUN), PR2(IRUN)
314 molod 1.3
315     C **********************************************************************
316     C **** INTERPOLATE ozone data to model latitudes ***
317     C **********************************************************************
318    
319 molod 1.8 DO 32 K=1,nozlevs
320 molod 1.3 DO 34 I=1,IRUN
321    
322 molod 1.8 DO 36 jlat = 1,nozlats
323     IF( ozlats(jlat).gt.xlat(i) ) THEN
324 molod 1.3 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 molod 1.8 slope(i) = ( XLAT(I) -ozlats(jlat-1) )
332     . / ( ozlats(jlat)-ozlats(jlat-1) )
333 molod 1.3 ENDIF
334     GOTO 37
335     ENDIF
336     36 CONTINUE
337 molod 1.8 jlatm = nozlats
338     jlatp = nozlats
339 molod 1.3 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 molod 1.8 DO 46 L1 = 1,nozlevs
359     IF( ozlevs(L1).GT.PLEVS(I,L2) ) THEN
360 molod 1.3 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 molod 1.8 LM = nozlevs-1
371     LP = nozlevs
372 molod 1.3 47 CONTINUE
373 molod 1.8 PR1(I) = ozlevs (LM)
374     PR2(I) = ozlevs (LP)
375 molod 1.3 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