1 | SUBROUTINE concentrations(ngrid,nlayer,nq, |
---|
2 | & pplay,pt,pdt,pq,pdq,ptimestep) |
---|
3 | |
---|
4 | use tracer_mod, only: igcm_co2, igcm_co, igcm_o, igcm_o1d, |
---|
5 | & igcm_o2, igcm_o3, igcm_h, igcm_h2, |
---|
6 | & igcm_oh, igcm_ho2, igcm_n2, igcm_ar, |
---|
7 | & igcm_h2o_vap, igcm_n, igcm_no, igcm_no2, |
---|
8 | & igcm_n2d, igcm_co2plus, igcm_oplus, |
---|
9 | & igcm_o2plus, igcm_coplus, igcm_cplus, |
---|
10 | & igcm_nplus, igcm_noplus, igcm_n2plus, |
---|
11 | & igcm_hplus, igcm_hco2plus, mmol |
---|
12 | use conc_mod, only: mmean, Akknew, rnew, cpnew |
---|
13 | USE comcstfi_h |
---|
14 | implicit none |
---|
15 | |
---|
16 | !======================================================================= |
---|
17 | ! CALCULATION OF MEAN MOLECULAR MASS, Cp, Akk and R |
---|
18 | ! |
---|
19 | ! mmean(ngrid,nlayer) amu |
---|
20 | ! cpnew(ngrid,nlayer) J/kg/K |
---|
21 | ! rnew(ngrid,nlayer) J/kg/K |
---|
22 | ! akknew(ngrid,nlayer) coefficient of thermal concduction |
---|
23 | ! |
---|
24 | ! version: April 2012 - Franck Lefevre |
---|
25 | !======================================================================= |
---|
26 | |
---|
27 | ! declarations |
---|
28 | |
---|
29 | #include "callkeys.h" |
---|
30 | #include "chimiedata.h" |
---|
31 | |
---|
32 | ! input/output |
---|
33 | |
---|
34 | integer,intent(in) :: ngrid ! number of atmospheric columns |
---|
35 | integer,intent(in) :: nlayer ! number of atmospheric layers |
---|
36 | integer,intent(in) :: nq ! number of tracers |
---|
37 | real,intent(in) :: pplay(ngrid,nlayer) |
---|
38 | real,intent(in) :: pt(ngrid,nlayer) |
---|
39 | real,intent(in) :: pdt(ngrid,nlayer) |
---|
40 | real,intent(in) :: pq(ngrid,nlayer,nq) |
---|
41 | real,intent(in) :: pdq(ngrid,nlayer,nq) |
---|
42 | real,intent(in) :: ptimestep |
---|
43 | |
---|
44 | ! local variables |
---|
45 | |
---|
46 | integer :: i, l, ig, iq |
---|
47 | integer, save :: nbq |
---|
48 | integer,allocatable,save :: niq(:) |
---|
49 | real :: ni(nq), ntot |
---|
50 | real :: zq(ngrid, nlayer, nq) |
---|
51 | real :: zt(ngrid, nlayer) |
---|
52 | real,allocatable,save :: aki(:) |
---|
53 | real,allocatable,save :: cpi(:) |
---|
54 | |
---|
55 | logical, save :: firstcall = .true. |
---|
56 | |
---|
57 | if (firstcall) then |
---|
58 | |
---|
59 | ! allocate local saved arrays: |
---|
60 | allocate(aki(nq)) |
---|
61 | allocate(cpi(nq)) |
---|
62 | allocate(niq(nq)) |
---|
63 | ! find index of chemical tracers to use |
---|
64 | ! initialize thermal conductivity and specific heat coefficients |
---|
65 | ! !? values are estimated |
---|
66 | |
---|
67 | nbq = 0 ! to count number of tracers used in this subroutine |
---|
68 | |
---|
69 | if (igcm_co2 /= 0) then |
---|
70 | nbq = nbq + 1 |
---|
71 | niq(nbq) = igcm_co2 |
---|
72 | aki(nbq) = 3.072e-4 |
---|
73 | cpi(nbq) = 0.834e3 |
---|
74 | end if |
---|
75 | if (igcm_co /= 0) then |
---|
76 | nbq = nbq + 1 |
---|
77 | niq(nbq) = igcm_co |
---|
78 | aki(nbq) = 4.87e-4 |
---|
79 | cpi(nbq) = 1.034e3 |
---|
80 | end if |
---|
81 | if (igcm_o /= 0) then |
---|
82 | nbq = nbq + 1 |
---|
83 | niq(nbq) = igcm_o |
---|
84 | aki(nbq) = 7.59e-4 |
---|
85 | cpi(nbq) = 1.3e3 |
---|
86 | end if |
---|
87 | if (igcm_o1d /= 0) then |
---|
88 | nbq = nbq + 1 |
---|
89 | niq(nbq) = igcm_o1d |
---|
90 | aki(nbq) = 7.59e-4 !? |
---|
91 | cpi(nbq) = 1.3e3 !? |
---|
92 | end if |
---|
93 | if (igcm_o2 /= 0) then |
---|
94 | nbq = nbq + 1 |
---|
95 | niq(nbq) = igcm_o2 |
---|
96 | aki(nbq) = 5.68e-4 |
---|
97 | cpi(nbq) = 0.9194e3 |
---|
98 | end if |
---|
99 | if (igcm_o3 /= 0) then |
---|
100 | nbq = nbq + 1 |
---|
101 | niq(nbq) = igcm_o3 |
---|
102 | aki(nbq) = 3.00e-4 !? |
---|
103 | cpi(nbq) = 0.800e3 !? |
---|
104 | end if |
---|
105 | if (igcm_h /= 0) then |
---|
106 | nbq = nbq + 1 |
---|
107 | niq(nbq) = igcm_h |
---|
108 | aki(nbq) = 0.0 |
---|
109 | cpi(nbq) = 20.780e3 |
---|
110 | end if |
---|
111 | if (igcm_h2 /= 0) then |
---|
112 | nbq = nbq + 1 |
---|
113 | niq(nbq) = igcm_h2 |
---|
114 | aki(nbq) = 36.314e-4 |
---|
115 | cpi(nbq) = 14.266e3 |
---|
116 | end if |
---|
117 | if (igcm_oh /= 0) then |
---|
118 | nbq = nbq + 1 |
---|
119 | niq(nbq) = igcm_oh |
---|
120 | aki(nbq) = 7.00e-4 !? |
---|
121 | cpi(nbq) = 1.045e3 |
---|
122 | end if |
---|
123 | if (igcm_ho2 /= 0) then |
---|
124 | nbq = nbq + 1 |
---|
125 | niq(nbq) = igcm_ho2 |
---|
126 | aki(nbq) = 0.0 |
---|
127 | cpi(nbq) = 1.065e3 !? |
---|
128 | end if |
---|
129 | if (igcm_n2 /= 0) then |
---|
130 | nbq = nbq + 1 |
---|
131 | niq(nbq) = igcm_n2 |
---|
132 | aki(nbq) = 5.6e-4 |
---|
133 | cpi(nbq) = 1.034e3 |
---|
134 | end if |
---|
135 | if (igcm_ar /= 0) then |
---|
136 | nbq = nbq + 1 |
---|
137 | niq(nbq) = igcm_ar |
---|
138 | aki(nbq) = 0.0 !? |
---|
139 | cpi(nbq) = 1.000e3 !? |
---|
140 | end if |
---|
141 | if (igcm_h2o_vap /= 0) then |
---|
142 | nbq = nbq + 1 |
---|
143 | niq(nbq) = igcm_h2o_vap |
---|
144 | aki(nbq) = 0.0 |
---|
145 | cpi(nbq) = 1.870e3 |
---|
146 | end if |
---|
147 | if (igcm_n /= 0) then |
---|
148 | nbq = nbq + 1 |
---|
149 | niq(nbq) = igcm_n |
---|
150 | aki(nbq) = 0.0 |
---|
151 | cpi(nbq) = 0.0 |
---|
152 | endif |
---|
153 | if(igcm_no /= 0) then |
---|
154 | nbq = nbq + 1 |
---|
155 | niq(nbq) = igcm_no |
---|
156 | aki(nbq) = 0.0 |
---|
157 | cpi(nbq) = 0.0 |
---|
158 | endif |
---|
159 | if(igcm_no2 /= 0) then |
---|
160 | nbq = nbq + 1 |
---|
161 | niq(nbq) = igcm_no2 |
---|
162 | aki(nbq) = 0.0 |
---|
163 | cpi(nbq) = 0.0 |
---|
164 | endif |
---|
165 | if(igcm_n2d /= 0) then |
---|
166 | nbq = nbq + 1 |
---|
167 | niq(nbq) = igcm_n2d |
---|
168 | aki(nbq) = 0.0 |
---|
169 | cpi(nbq) = 0.0 |
---|
170 | endif |
---|
171 | if(igcm_co2plus /= 0) then |
---|
172 | nbq = nbq + 1 |
---|
173 | niq(nbq) = igcm_co2plus |
---|
174 | aki(nbq) = 0.0 |
---|
175 | cpi(nbq) = 0.0 |
---|
176 | endif |
---|
177 | if(igcm_oplus /= 0) then |
---|
178 | nbq = nbq + 1 |
---|
179 | niq(nbq) = igcm_oplus |
---|
180 | aki(nbq) = 0.0 |
---|
181 | cpi(nbq) = 0.0 |
---|
182 | endif |
---|
183 | if(igcm_o2plus /= 0) then |
---|
184 | nbq = nbq + 1 |
---|
185 | niq(nbq) = igcm_o2plus |
---|
186 | aki(nbq) = 0.0 |
---|
187 | cpi(nbq) = 0.0 |
---|
188 | endif |
---|
189 | if(igcm_coplus /= 0) then |
---|
190 | nbq = nbq + 1 |
---|
191 | niq(nbq) = igcm_coplus |
---|
192 | aki(nbq) = 0.0 |
---|
193 | cpi(nbq) = 0.0 |
---|
194 | endif |
---|
195 | if(igcm_cplus /= 0) then |
---|
196 | nbq = nbq + 1 |
---|
197 | niq(nbq) = igcm_cplus |
---|
198 | aki(nbq) = 0.0 |
---|
199 | cpi(nbq) = 0.0 |
---|
200 | endif |
---|
201 | if(igcm_nplus /= 0) then |
---|
202 | nbq = nbq + 1 |
---|
203 | niq(nbq) = igcm_nplus |
---|
204 | aki(nbq) = 0.0 |
---|
205 | cpi(nbq) = 0.0 |
---|
206 | endif |
---|
207 | if(igcm_noplus /= 0) then |
---|
208 | nbq = nbq + 1 |
---|
209 | niq(nbq) = igcm_noplus |
---|
210 | aki(nbq) = 0.0 |
---|
211 | cpi(nbq) = 0.0 |
---|
212 | endif |
---|
213 | if(igcm_n2plus /= 0) then |
---|
214 | nbq = nbq + 1 |
---|
215 | niq(nbq) = igcm_n2plus |
---|
216 | aki(nbq) = 0.0 |
---|
217 | cpi(nbq) = 0.0 |
---|
218 | endif |
---|
219 | if(igcm_hplus /= 0) then |
---|
220 | nbq = nbq + 1 |
---|
221 | niq(nbq) = igcm_hplus |
---|
222 | aki(nbq) = 0.0 |
---|
223 | cpi(nbq) = 0.0 |
---|
224 | endif |
---|
225 | if(igcm_hco2plus /= 0) then |
---|
226 | nbq = nbq + 1 |
---|
227 | niq(nbq) = igcm_hco2plus |
---|
228 | aki(nbq) = 0.0 |
---|
229 | cpi(nbq) = 0.0 |
---|
230 | endif |
---|
231 | |
---|
232 | ! tell the world about it: |
---|
233 | write(*,*) "concentrations: firstcall, nbq=",nbq |
---|
234 | write(*,*) " niq(1:nbq)=",niq(1:nbq) |
---|
235 | write(*,*) " aki(1:nbq)=",aki(1:nbq) |
---|
236 | write(*,*) " cpi(1:nbq)=",cpi(1:nbq) |
---|
237 | |
---|
238 | firstcall = .false. |
---|
239 | |
---|
240 | end if ! if (firstcall) |
---|
241 | |
---|
242 | ! update temperature |
---|
243 | |
---|
244 | do l = 1,nlayer |
---|
245 | do ig = 1,ngrid |
---|
246 | zt(ig,l) = pt(ig,l) + pdt(ig,l)*ptimestep |
---|
247 | end do |
---|
248 | end do |
---|
249 | |
---|
250 | ! update tracers |
---|
251 | |
---|
252 | do l = 1,nlayer |
---|
253 | do ig = 1,ngrid |
---|
254 | do i = 1,nbq |
---|
255 | iq = niq(i) |
---|
256 | zq(ig,l,iq) = max(1.e-30, pq(ig,l,iq) |
---|
257 | $ + pdq(ig,l,iq)*ptimestep) |
---|
258 | end do |
---|
259 | end do |
---|
260 | end do |
---|
261 | |
---|
262 | ! mmean : mean molecular mass |
---|
263 | ! rnew : specific gas constant |
---|
264 | |
---|
265 | mmean(:,:) = 0. |
---|
266 | |
---|
267 | do l = 1,nlayer |
---|
268 | do ig = 1,ngrid |
---|
269 | do i = 1,nbq |
---|
270 | iq = niq(i) |
---|
271 | mmean(ig,l) = mmean(ig,l) + zq(ig,l,iq)/mmol(iq) |
---|
272 | end do |
---|
273 | mmean(ig,l) = 1./mmean(ig,l) |
---|
274 | rnew(ig,l) = 8.314/mmean(ig,l)*1.e3 ! J/kg/K |
---|
275 | end do |
---|
276 | end do |
---|
277 | |
---|
278 | ! cpnew : specicic heat |
---|
279 | ! akknew : thermal conductivity cofficient |
---|
280 | |
---|
281 | cpnew(:,:) = 0. |
---|
282 | akknew(:,:) = 0. |
---|
283 | |
---|
284 | do l = 1,nlayer |
---|
285 | do ig = 1,ngrid |
---|
286 | ntot = pplay(ig,l)/(1.381e-23*zt(ig,l))*1.e-6 ! in #/cm3 |
---|
287 | do i = 1,nbq |
---|
288 | iq = niq(i) |
---|
289 | ni(iq) = ntot*zq(ig,l,iq)*mmean(ig,l)/mmol(iq) |
---|
290 | cpnew(ig,l) = cpnew(ig,l) + ni(iq)*cpi(i) |
---|
291 | akknew(ig,l) = akknew(ig,l) + ni(iq)*aki(i) |
---|
292 | end do |
---|
293 | cpnew(ig,l) = cpnew(ig,l)/ntot |
---|
294 | akknew(ig,l) = akknew(ig,l)/ntot |
---|
295 | end do |
---|
296 | ! print*, l, mmean(1,l), cpnew(1,l), rnew(1,l) |
---|
297 | end do |
---|
298 | |
---|
299 | return |
---|
300 | end |
---|