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

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

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


Revision 1.2 - (show annotations) (download)
Tue Jun 8 15:09:01 2004 UTC (20 years ago) by molod
Branch: MAIN
Changes since 1.1: +0 -0 lines
Developing input for fizhi

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

  ViewVC Help
Powered by ViewVC 1.1.22