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

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

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


Revision 1.12 - (hide annotations) (download)
Thu Aug 5 20:09:48 2004 UTC (19 years, 10 months ago) by molod
Branch: MAIN
Changes since 1.11: +3 -3 lines
Do byteswapping for ozone data read

1 molod 1.12 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_init_chem.F,v 1.11 2004/07/28 01:25:07 molod Exp $
2 molod 1.1 C $Name: $
3    
4 molod 1.9 #include "FIZHI_OPTIONS.h"
5 molod 1.7 subroutine fizhi_init_chem(mythid,nozlats,nozlevs,ntimesoz,
6     . ozlats,ozlevs,o3,
7     . nwatlats,nwatlevs,ntimesq,watlats,watlevs,water,
8 molod 1.6 . Nrphys,pressure,n2o,methane,co2,cfc11,cfc12,cfc22)
9 molod 1.1 C***********************************************************************
10     C Subroutine fizhi_init_chem - routine to read in the ozone and upper
11     C atmosphere water vapor, and set the methane, n2o and cfc values
12     C
13     C INPUT:
14     C
15     C mythid - thread number (processor number)
16     C
17     C OUTPUT:
18     C
19     C***********************************************************************
20     implicit none
21     integer mythid,nozlevs,nozlats,nwatlevs,nwatlats,Nrphys
22 molod 1.7 integer ntimesoz,ntimesq
23     _RL o3(nozlats,nozlevs,ntimesoz)
24     _RL water(nwatlats,nwatlevs,ntimesq)
25 molod 1.1 _RL ozlats(nozlats), ozlevs(nozlevs)
26     _RL watlats(nwatlats), watlevs(nwatlevs)
27     _RL pressure(Nrphys),methane(Nrphys),n2o(Nrphys)
28     _RL co2,cfc11,cfc12,cfc22
29    
30 molod 1.9 _RL getcon
31 molod 1.6 integer kqz,koz
32     call mdsfindunit( kqz, myThid )
33     call mdsfindunit( koz, myThid )
34 molod 1.1
35 molod 1.7 call read_qz (kqz,water,watlats,watlevs,nwatlats,nwatlevs,ntimesq)
36     call read_oz (koz,o3,ozlats,ozlevs,nozlats,nozlevs,ntimesoz)
37 molod 1.6
38 molod 1.1 call get_methane_n2o (pressure,Nrphys,n2o,methane)
39 molod 1.6
40 molod 1.1 co2 = getcon('CO2' )*1.e-6
41     cfc11 = getcon('CFC11')*1.e-9
42     cfc12 = getcon('CFC12')*1.e-9
43     cfc22 = getcon('CFC22')*1.e-9
44    
45     RETURN
46     END
47    
48     subroutine read_qz (ku,qz,lats,levs,nlat,nlev,ntime)
49     C***********************************************************************
50     C PURPOSE
51     C To Read Stratospheric Moisture Data
52     C
53     C ARGUMENTS DESCRIPTION
54     C ku ...... Unit to Read Moisture Data
55     C qz ...... Stratospheric Moisture Data
56     C lats .... Stratospheric Moisture Data Latitudes (degrees)
57     C levs .... Stratospheric Moisture Data Levels (mb)
58     C nlat .... Number of ozone latitudes
59     C nlev .... Number of ozone levels
60     C ntime ... Number of ozone time values
61     C
62     C***********************************************************************
63    
64     implicit none
65     integer ku,nlat,nlev,ntime
66    
67 molod 1.6 _RL qz(nlat,nlev,ntime)
68 molod 1.3 _RL lats(nlat)
69     _RL levs(nlev)
70 molod 1.1
71     integer time
72     integer lat
73     integer lev
74    
75 molod 1.9 _RL voltomas
76 molod 1.1 parameter ( voltomas = 0.622e-6 )
77    
78 molod 1.7 open(ku,file='sage.data',form='formatted')
79 molod 1.1 rewind ku
80    
81     c Set Moisture Data Latitudes
82     c ---------------------------
83     do lat = 1,nlat
84     lats(lat) = -85. + (lat-1)*10.
85     enddo
86    
87     c Read Moisture Pressure Levels
88     c -----------------------------
89     read(ku,1000) (levs(lev),lev=1,nlev)
90    
91     c Read Moisture Amounts by Month and Level
92     c ----------------------------------------
93     do time=1,ntime
94     read (ku,1001)
95     do lat=1,nlat
96     read(ku,1000) (qz(lat,lev,time),lev=1,nlev)
97     enddo
98     enddo
99    
100     c Convert from Volume Mixing Ratio to Mass Mixing Ratio
101     c -----------------------------------------------------
102     do time = 1,ntime
103     do lev = 1,nlev
104     do lat = 1,nlat
105     qz(lat,lev,time) = qz(lat,lev,time)*voltomas
106     enddo
107     enddo
108     enddo
109    
110     1000 format (3(5x,7(2x,f6.1)/))
111     1001 format (1x)
112     return
113     end
114    
115     subroutine read_oz (ku,oz,lats,levs,nlat,nlev,ntime)
116     C***********************************************************************
117     C PURPOSE
118     C To Read Ozone Value
119     C
120     C ARGUMENTS DESCRIPTION
121     C ku ...... Unit to Read Ozone Data
122     C oz ...... Ozone Data
123     C lats .... Ozone Data Latitudes (degrees)
124     C levs .... Ozone Data Levels (mb)
125     C nlat .... Number of ozone latitudes
126     C nlev .... Number of ozone levels
127     C ntime ... Number of ozone time values
128     C
129     C***********************************************************************
130    
131     implicit none
132     integer ku,nlat,nlev,ntime
133    
134     _RL oz(nlat,nlev,ntime)
135 molod 1.10 real*4 o3(nlat)
136 molod 1.3 _RL lats(nlat)
137     _RL levs(nlev)
138 molod 1.1
139     integer time
140     integer lat
141     integer lev
142     integer nrec
143    
144 molod 1.3 _RL plevs(34)
145 molod 1.1 data plevs/ 0.003, 0.005, 0.007, 0.01, 0.015, 0.02, 0.03, 0.05,
146     . 0.07, 0.1, 0.15, 0.2, 0.3, 0.5, 0.7, 1.0, 1.5, 2.0,
147     . 3.0, 5.0, 7.0, 10.0, 15.0, 20.0, 30.0, 50.0, 70.0,
148     . 100.0, 150.0, 200.0, 300.0, 500.0, 700.0, 1000.0 /
149    
150     c Set Ozone Data Latitudes
151     c ------------------------
152     do lat = 1,nlat
153     lats(lat) = -90. + (lat-1)*5.
154     enddo
155    
156     c Set Ozone Data Levels
157     c ------------------------
158     do lev = 1,nlev
159     levs(lev) = plevs(lev)
160     enddo
161    
162     c Read Ozone Amounts by Month and Level
163     c -------------------------------------
164     close (ku)
165 molod 1.7 open(ku,file='gcmo3.data',form='unformatted',access='direct',
166     . recl=nlat*4)
167 molod 1.1
168     do time=1,ntime
169     do lev=1,nlev
170 molod 1.5 C Note: 2 quantities in Ozone Dataset
171     nrec = lev+(time-1)*nlev*2
172 molod 1.1 read(ku,rec=nrec) o3
173 molod 1.12 call mds_byteswapr4(nlat,o3)
174 molod 1.1 do lat=1,nlat
175 molod 1.12 oz(lat,nlev-lev+1,time) = o3(lat)
176 molod 1.1 enddo
177     enddo
178     enddo
179    
180     close (ku)
181     return
182     end
183    
184     subroutine get_methane_n2o (pres,Nrphys,n2o,methane)
185     C***********************************************************************
186     C PURPOSE
187     C Compute methane and n2o
188     C
189     C ARGUMENTS DESCRIPTION
190     C
191     C***********************************************************************
192     C* Climatological Annual and Global Mean Height Data *
193     C***********************************************************************
194    
195     implicit none
196    
197 molod 1.6 integer Nrphys
198 molod 1.1 _RL n2o(Nrphys),methane(Nrphys)
199     _RL pres(Nrphys)
200 molod 1.3 _RL hght(Nrphys), slope,pr1,pr2,hpr1,hpr2
201 molod 1.1 integer L,L1,L2,lup,ldn
202    
203 molod 1.3 _RL plevc (46), plevz(46)
204     _RL hghtc (46), hghtz(46)
205 molod 1.1
206     data plevc /1000.00, 975.00, 950.00, 925.00, 900.00,
207     . 875.00, 850.00, 825.00, 800.00, 750.00,
208     . 700.00, 650.00, 600.00, 550.00, 500.00,
209     . 450.00, 400.00, 350.00, 300.00, 250.00,
210     . 200.00, 150.00, 100.00, 70.00, 50.00,
211     . 40.00, 30.00, 20.00, 10.00, 7.00,
212     . 5.00, 4.00, 3.00, 2.00, 1.00,
213     . 0.70, 0.50, 0.40, 0.30, 0.20,
214     . 0.10, 0.07, 0.05, 0.04, 0.03,
215     . 0.02 /
216    
217     data hghtc/ 0.128733 , 0.316985 , 0.528275 , 0.749515 , 0.976471 ,
218     . 1.208910 , 1.446800 , 1.690980 , 1.941630 , 2.463530 ,
219     . 3.016200 , 3.603490 , 4.229410 , 4.899870 , 5.622320 ,
220     . 6.405940 , 7.263450 , 8.211920 , 9.275540 , 10.49150 ,
221     . 11.92420 , 13.70200 , 16.12980 , 18.24120 , 20.26480 ,
222     . 21.63100 , 23.41250 , 25.96570 , 30.45890 , 32.85240 ,
223     . 35.17360 , 36.75040 , 38.82900 , 41.84600 , 47.15580 ,
224     . 49.90100 , 52.46230 , 54.13890 , 56.26340 , 59.17640 ,
225     . 63.89980 , 66.20240 , 68.29210 , 69.63550 , 71.32330 ,
226     . 73.62110 /
227    
228     do L=1,46
229     plevz(L) = plevc(47-L)
230     hghtz(L) = hghtc(47-L)
231     enddo
232    
233     C **********************************************************************
234     C Interpolate Heights to Model Pressures ****
235     C **********************************************************************
236    
237     do L2 = 1,Nrphys
238    
239     do L1 = 1,46
240     if( plevz(L1).gt.pres(L2) ) then
241     if( L1.eq.1 ) then
242     lup = 1
243     ldn = 2
244     else
245     lup = L1-1
246     ldn = L1
247     endif
248     goto 10
249     endif
250     enddo
251     lup = 45
252     ldn = 46
253    
254     10 continue
255     pr1 = plevz(lup)
256     pr2 = plevz(ldn)
257     hpr1 = hghtz(lup)
258     hpr2 = hghtz(ldn)
259    
260     slope = ( hpr1-hpr2 )/( pr1-pr2 )
261     hght(L2) = hpr2 + ( pres(L2)-pr2 )*slope
262    
263     enddo
264    
265     C **********************************************************************
266     C Set the profiles of N2O and CH4 based on Bresser and Pawson 1996 ****
267     C **********************************************************************
268    
269     do L = 1,Nrphys
270     if( hght(L).gt.26. ) then
271     n2o(L) = 120.* exp( (26.- hght(L)) / 5.69 ) * 1.e-9
272     else if( hght(L).gt.16. ) then
273     n2o(L) = 307.* exp( (16.- hght(L)) /10.47 ) * 1.e-9
274     else
275     n2o(L) = 307.e-9
276     endif
277     enddo
278    
279     do L = 1,Nrphys
280     if( hght(L).gt.55. ) then
281     methane(L) = 0.2e-6
282     else if( hght(L).gt.14. ) then
283     methane(L) = 1.7* exp( (14.- hght(L)) /19.16 ) * 1.e-6
284     else
285     methane(L) = 1.7e-6
286     endif
287     enddo
288    
289     return
290     end

  ViewVC Help
Powered by ViewVC 1.1.22