/[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.5 - (show annotations) (download)
Thu Jun 10 21:50:33 2004 UTC (20 years ago) by molod
Branch: MAIN
CVS Tags: checkpoint53d_post, checkpoint54a_pre, checkpoint54a_post, checkpoint54b_post, checkpoint54, checkpoint53g_post, checkpoint53f_post
Changes since 1.4: +4 -3 lines
Developing

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

  ViewVC Help
Powered by ViewVC 1.1.22