/[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.15 - (hide annotations) (download)
Wed Oct 20 16:41:21 2004 UTC (19 years, 7 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint58l_post, checkpoint64z, checkpoint57t_post, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64i, checkpoint64h, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint64a, checkpoint57o_post, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint58e_post, checkpoint57v_post, checkpoint58u_post, checkpoint58w_post, checkpoint57m_post, checkpoint57s_post, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint57k_post, checkpoint57d_post, checkpoint57g_post, checkpoint64, checkpoint65, checkpoint60, checkpoint61, checkpoint62, checkpoint63, checkpoint57b_post, checkpoint57c_pre, checkpoint58r_post, checkpoint55j_post, checkpoint56b_post, checkpoint57i_post, checkpoint57y_post, checkpoint57e_post, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint58n_post, checkpoint58x_post, checkpoint57g_pre, checkpoint58t_post, checkpoint58h_post, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, checkpoint56c_post, checkpoint57y_pre, checkpoint57f_pre, checkpoint57a_post, checkpoint58q_post, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint58j_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint57r_post, checkpoint59, checkpoint58, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, eckpoint57e_pre, checkpoint57h_done, checkpoint58f_post, checkpoint57x_post, checkpoint57n_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint57f_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58g_post, checkpoint58o_post, checkpoint57z_post, checkpoint62c, checkpoint62b, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62w, checkpoint62v, checkpoint62u, checkpoint62t, checkpoint57c_post, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint58y_post, checkpoint58k_post, checkpoint58v_post, checkpoint58s_post, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint58p_post, checkpoint61a, checkpoint61n, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61q, checkpoint57j_post, checkpoint61z, checkpoint61x, checkpoint61y, checkpoint58b_post, checkpoint57h_pre, checkpoint58m_post, checkpoint57l_post, checkpoint57h_post, checkpoint56a_post, HEAD
Changes since 1.14: +3 -3 lines
Change file names for input data fields

1 molod 1.15 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_init_chem.F,v 1.14 2004/09/01 18:35:33 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 molod 1.14 integer time0
72 molod 1.1 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.15 open(ku,file='data.sage',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 molod 1.14 do time0=1,ntime
94 molod 1.1 read (ku,1001)
95     do lat=1,nlat
96 molod 1.14 read(ku,1000) (qz(lat,lev,time0),lev=1,nlev)
97 molod 1.1 enddo
98     enddo
99    
100     c Convert from Volume Mixing Ratio to Mass Mixing Ratio
101     c -----------------------------------------------------
102 molod 1.14 do time0 = 1,ntime
103 molod 1.1 do lev = 1,nlev
104     do lat = 1,nlat
105 molod 1.14 qz(lat,lev,time0) = qz(lat,lev,time0)*voltomas
106 molod 1.1 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 molod 1.14 integer time0
140 molod 1.1 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.15 open(ku,file='data.gcmo3',form='unformatted',access='direct',
166 molod 1.7 . recl=nlat*4)
167 molod 1.1
168 molod 1.14 do time0=1,ntime
169 molod 1.1 do lev=1,nlev
170 molod 1.5 C Note: 2 quantities in Ozone Dataset
171 molod 1.14 nrec = lev+(time0-1)*nlev*2
172 molod 1.1 read(ku,rec=nrec) o3
173 molod 1.13 #if defined( _BYTESWAPIO )
174 molod 1.12 call mds_byteswapr4(nlat,o3)
175 molod 1.13 #endif
176 molod 1.1 do lat=1,nlat
177 molod 1.14 oz(lat,nlev-lev+1,time0) = o3(lat)
178 molod 1.1 enddo
179     enddo
180     enddo
181    
182     close (ku)
183     return
184     end
185    
186     subroutine get_methane_n2o (pres,Nrphys,n2o,methane)
187     C***********************************************************************
188     C PURPOSE
189     C Compute methane and n2o
190     C
191     C ARGUMENTS DESCRIPTION
192     C
193     C***********************************************************************
194     C* Climatological Annual and Global Mean Height Data *
195     C***********************************************************************
196    
197     implicit none
198    
199 molod 1.6 integer Nrphys
200 molod 1.1 _RL n2o(Nrphys),methane(Nrphys)
201     _RL pres(Nrphys)
202 molod 1.3 _RL hght(Nrphys), slope,pr1,pr2,hpr1,hpr2
203 molod 1.1 integer L,L1,L2,lup,ldn
204    
205 molod 1.3 _RL plevc (46), plevz(46)
206     _RL hghtc (46), hghtz(46)
207 molod 1.1
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