1 | subroutine inichim_newstart(ngrid, nq, pq, qsurf, ps, & |
---|
2 | flagh2o, flagthermo) |
---|
3 | |
---|
4 | use tracer_mod |
---|
5 | USE vertical_layers_mod, ONLY: aps,bps |
---|
6 | USE mod_grid_phy_lmdz, ONLY: nbp_lon, nbp_lat, nbp_lev |
---|
7 | USE datafile_mod, ONLY: datadir |
---|
8 | implicit none |
---|
9 | |
---|
10 | !======================================================================= |
---|
11 | ! |
---|
12 | ! Purpose: |
---|
13 | ! -------- |
---|
14 | ! |
---|
15 | ! Initialization of the chemistry for use in the building of a new start file |
---|
16 | ! used by program newstart.F |
---|
17 | ! also used by program testphys1d.F |
---|
18 | ! |
---|
19 | ! Authors: |
---|
20 | ! ------- |
---|
21 | ! Initial version 11/2002 by Sebastien Lebonnois |
---|
22 | ! Revised 07/2003 by Monica Angelats-i-Coll to use input files |
---|
23 | ! Modified 10/2008 Identify tracers by their names Ehouarn Millour |
---|
24 | ! Modified 11/2011 Addition of methane Franck Lefevre |
---|
25 | ! Rewritten 04/2012 Franck Lefevre |
---|
26 | ! |
---|
27 | ! Arguments: |
---|
28 | ! ---------- |
---|
29 | ! |
---|
30 | ! pq(nbp_lon+1,nbp_lat,nbp_lev,nq) Advected fields, ie chemical species here |
---|
31 | ! qsurf(ngrid,nq) Amount of tracer on the surface (kg/m2) |
---|
32 | ! ps(nbp_lon+1,nbp_lat) Surface pressure (Pa) |
---|
33 | ! flagh2o flag for initialisation of h2o (1: yes / 0: no) |
---|
34 | ! flagthermo flag for initialisation of thermosphere only (1: yes / 0: no) |
---|
35 | ! |
---|
36 | !======================================================================= |
---|
37 | |
---|
38 | include "callkeys.h" |
---|
39 | |
---|
40 | ! inputs : |
---|
41 | |
---|
42 | integer,intent(in) :: ngrid ! number of atmospheric columns in the physics |
---|
43 | integer,intent(in) :: nq ! number of tracers |
---|
44 | real,intent(in) :: ps(nbp_lon+1,nbp_lat) ! surface pressure in the gcm (Pa) |
---|
45 | integer,intent(in) :: flagh2o ! flag for h2o initialisation |
---|
46 | integer,intent(in) :: flagthermo ! flag for thermosphere initialisation only |
---|
47 | |
---|
48 | ! outputs : |
---|
49 | |
---|
50 | real,intent(out) :: pq(nbp_lon+1,nbp_lat,nbp_lev,nq) ! advected fields, ie chemical species |
---|
51 | real,intent(out) :: qsurf(ngrid,nq) ! surface values (kg/m2) of tracers |
---|
52 | |
---|
53 | ! local : |
---|
54 | |
---|
55 | integer :: iq, i, j, l, n |
---|
56 | integer :: count, ierr, dummy |
---|
57 | real :: mmean(nbp_lon+1,nbp_lat,nbp_lev) ! mean molecular mass (g) |
---|
58 | real :: pgcm ! pressure at each layer in the gcm (Pa) |
---|
59 | |
---|
60 | integer, parameter :: nalt = 252 ! number of altitudes in the initialization files |
---|
61 | integer :: nspe ! number of species in the initialization files |
---|
62 | integer, allocatable :: niq(:) ! local index of species in initialization files |
---|
63 | real, dimension(nalt) :: tinit, zzfile ! temperature in initialization files |
---|
64 | real, dimension(nalt) :: pinit ! pressure in initialization files |
---|
65 | real, dimension(nalt) :: densinit ! total number density in initialization files |
---|
66 | real, allocatable :: vmrinit(:,:) ! mixing ratios in initialization files |
---|
67 | real, allocatable :: vmrint(:) ! mixing ratio interpolated onto the gcm vertical grid |
---|
68 | real :: vmr |
---|
69 | |
---|
70 | character(len=20) :: txt ! to store some text |
---|
71 | logical :: flagnitro ! checks if N species present |
---|
72 | |
---|
73 | ! 1. identify tracers by their names: (and set corresponding values of mmol) |
---|
74 | |
---|
75 | ! 1.1 initialize tracer indexes to zero: |
---|
76 | nqmx=nq ! initialize value of nqmx |
---|
77 | |
---|
78 | igcm_dustbin(1:nq)=0 |
---|
79 | igcm_co2_ice=0 |
---|
80 | igcm_ccnco2_mass=0 |
---|
81 | igcm_ccnco2_number=0 |
---|
82 | igcm_dust_mass=0 |
---|
83 | igcm_dust_number=0 |
---|
84 | igcm_ccn_mass=0 |
---|
85 | igcm_ccn_number=0 |
---|
86 | igcm_dust_submicron=0 |
---|
87 | igcm_h2o_vap=0 |
---|
88 | igcm_h2o_ice=0 |
---|
89 | igcm_co2=0 |
---|
90 | igcm_co=0 |
---|
91 | igcm_o=0 |
---|
92 | igcm_o1d=0 |
---|
93 | igcm_o2=0 |
---|
94 | igcm_o3=0 |
---|
95 | igcm_h=0 |
---|
96 | igcm_h2=0 |
---|
97 | igcm_oh=0 |
---|
98 | igcm_ho2=0 |
---|
99 | igcm_h2o2=0 |
---|
100 | igcm_ch4=0 |
---|
101 | igcm_n2=0 |
---|
102 | igcm_ar=0 |
---|
103 | igcm_ar_n2=0 |
---|
104 | igcm_n=0 |
---|
105 | igcm_no=0 |
---|
106 | igcm_no2=0 |
---|
107 | igcm_n2d=0 |
---|
108 | igcm_he=0 |
---|
109 | igcm_co2plus=0 |
---|
110 | igcm_oplus=0 |
---|
111 | igcm_o2plus=0 |
---|
112 | igcm_coplus=0 |
---|
113 | igcm_cplus=0 |
---|
114 | igcm_nplus=0 |
---|
115 | igcm_noplus=0 |
---|
116 | igcm_n2plus=0 |
---|
117 | igcm_hplus=0 |
---|
118 | igcm_hco2plus=0 |
---|
119 | igcm_hcoplus=0 |
---|
120 | igcm_h2oplus=0 |
---|
121 | igcm_h3oplus=0 |
---|
122 | igcm_ohplus=0 |
---|
123 | igcm_elec=0 |
---|
124 | |
---|
125 | ! 1.2 find dust tracers |
---|
126 | |
---|
127 | count = 0 |
---|
128 | |
---|
129 | if (dustbin > 0) then |
---|
130 | do iq = 1,nqmx |
---|
131 | txt = " " |
---|
132 | write(txt,'(a4,i2.2)') 'dust', count + 1 |
---|
133 | if (noms(iq) == txt) then |
---|
134 | count = count + 1 |
---|
135 | igcm_dustbin(count) = iq |
---|
136 | mmol(iq) = 100. |
---|
137 | end if |
---|
138 | end do !do iq=1,nqmx |
---|
139 | end if ! of if (dustbin.gt.0) |
---|
140 | |
---|
141 | if (doubleq) then |
---|
142 | do iq = 1,nqmx |
---|
143 | if (noms(iq) == "dust_mass") then |
---|
144 | igcm_dust_mass = iq |
---|
145 | count = count + 1 |
---|
146 | end if |
---|
147 | if (noms(iq) == "dust_number") then |
---|
148 | igcm_dust_number = iq |
---|
149 | count = count + 1 |
---|
150 | end if |
---|
151 | end do |
---|
152 | end if ! of if (doubleq) |
---|
153 | |
---|
154 | if (scavenging) then |
---|
155 | do iq = 1,nqmx |
---|
156 | if (noms(iq) == "ccn_mass") then |
---|
157 | igcm_ccn_mass = iq |
---|
158 | count = count + 1 |
---|
159 | end if |
---|
160 | if (noms(iq) == "ccn_number") then |
---|
161 | igcm_ccn_number = iq |
---|
162 | count = count + 1 |
---|
163 | end if |
---|
164 | end do |
---|
165 | end if ! of if (scavenging) |
---|
166 | |
---|
167 | if (submicron) then |
---|
168 | do iq=1,nqmx |
---|
169 | if (noms(iq) == "dust_submicron") then |
---|
170 | igcm_dust_submicron = iq |
---|
171 | mmol(iq) = 100. |
---|
172 | count = count + 1 |
---|
173 | end if |
---|
174 | end do |
---|
175 | end if ! of if (submicron) |
---|
176 | |
---|
177 | ! 1.3 find chemistry and water tracers |
---|
178 | |
---|
179 | do iq = 1,nqmx |
---|
180 | if (noms(iq) == "co2") then |
---|
181 | igcm_co2 = iq |
---|
182 | mmol(igcm_co2) = 44. |
---|
183 | count = count + 1 |
---|
184 | end if |
---|
185 | if (noms(iq) == "co") then |
---|
186 | igcm_co = iq |
---|
187 | mmol(igcm_co) = 28. |
---|
188 | count = count + 1 |
---|
189 | end if |
---|
190 | if (noms(iq) == "o") then |
---|
191 | igcm_o = iq |
---|
192 | mmol(igcm_o) = 16. |
---|
193 | count = count + 1 |
---|
194 | end if |
---|
195 | if (noms(iq) == "o1d") then |
---|
196 | igcm_o1d = iq |
---|
197 | mmol(igcm_o1d) = 16. |
---|
198 | count = count + 1 |
---|
199 | end if |
---|
200 | if (noms(iq) == "o2") then |
---|
201 | igcm_o2 = iq |
---|
202 | mmol(igcm_o2) = 32. |
---|
203 | count = count + 1 |
---|
204 | end if |
---|
205 | if (noms(iq) == "o3") then |
---|
206 | igcm_o3 = iq |
---|
207 | mmol(igcm_o3) = 48. |
---|
208 | count = count + 1 |
---|
209 | end if |
---|
210 | if (noms(iq) == "h") then |
---|
211 | igcm_h = iq |
---|
212 | mmol(igcm_h) = 1. |
---|
213 | count = count + 1 |
---|
214 | end if |
---|
215 | if (noms(iq) == "h2") then |
---|
216 | igcm_h2 = iq |
---|
217 | mmol(igcm_h2) = 2. |
---|
218 | count = count + 1 |
---|
219 | end if |
---|
220 | if (noms(iq) == "oh") then |
---|
221 | igcm_oh = iq |
---|
222 | mmol(igcm_oh) = 17. |
---|
223 | count = count + 1 |
---|
224 | end if |
---|
225 | if (noms(iq) == "ho2") then |
---|
226 | igcm_ho2 = iq |
---|
227 | mmol(igcm_ho2) = 33. |
---|
228 | count = count + 1 |
---|
229 | end if |
---|
230 | if (noms(iq) == "h2o2") then |
---|
231 | igcm_h2o2 = iq |
---|
232 | mmol(igcm_h2o2) = 34. |
---|
233 | count = count + 1 |
---|
234 | end if |
---|
235 | if (noms(iq) == "ch4") then |
---|
236 | igcm_ch4 = iq |
---|
237 | mmol(igcm_ch4) = 16. |
---|
238 | count = count + 1 |
---|
239 | end if |
---|
240 | if (noms(iq) == "n2") then |
---|
241 | igcm_n2 = iq |
---|
242 | mmol(igcm_n2) = 28. |
---|
243 | count = count + 1 |
---|
244 | end if |
---|
245 | if (noms(iq) == "n") then |
---|
246 | igcm_n = iq |
---|
247 | mmol(igcm_n) = 14. |
---|
248 | count = count + 1 |
---|
249 | end if |
---|
250 | if (noms(iq) == "n2d") then |
---|
251 | igcm_n2d = iq |
---|
252 | mmol(igcm_n2d) = 14. |
---|
253 | count = count + 1 |
---|
254 | end if |
---|
255 | if (noms(iq) == "no") then |
---|
256 | igcm_no = iq |
---|
257 | mmol(igcm_no) = 30. |
---|
258 | count = count + 1 |
---|
259 | end if |
---|
260 | if (noms(iq) == "no2") then |
---|
261 | igcm_no2 = iq |
---|
262 | mmol(igcm_no2) = 46. |
---|
263 | count = count + 1 |
---|
264 | end if |
---|
265 | if (noms(iq) == "ar") then |
---|
266 | igcm_ar = iq |
---|
267 | mmol(igcm_ar) = 40. |
---|
268 | count = count + 1 |
---|
269 | end if |
---|
270 | if (noms(iq) == "h2o_vap") then |
---|
271 | igcm_h2o_vap = iq |
---|
272 | mmol(igcm_h2o_vap) = 18. |
---|
273 | count = count + 1 |
---|
274 | end if |
---|
275 | if (noms(iq) == "h2o_ice") then |
---|
276 | igcm_h2o_ice = iq |
---|
277 | mmol(igcm_h2o_ice) = 18. |
---|
278 | count = count + 1 |
---|
279 | end if |
---|
280 | if (noms(iq).eq."he") then |
---|
281 | igcm_he=iq |
---|
282 | mmol(igcm_he)=4. |
---|
283 | count=count+1 |
---|
284 | endif |
---|
285 | |
---|
286 | ! 1.4 find ions |
---|
287 | |
---|
288 | if (noms(iq) == "co2plus") then |
---|
289 | igcm_co2plus = iq |
---|
290 | mmol(igcm_co2plus) = 44. |
---|
291 | count = count + 1 |
---|
292 | end if |
---|
293 | if (noms(iq) == "oplus") then |
---|
294 | igcm_oplus = iq |
---|
295 | mmol(igcm_oplus) = 16. |
---|
296 | count = count + 1 |
---|
297 | end if |
---|
298 | if (noms(iq) == "o2plus") then |
---|
299 | igcm_o2plus = iq |
---|
300 | mmol(igcm_o2plus) = 32. |
---|
301 | count = count + 1 |
---|
302 | end if |
---|
303 | if (noms(iq) == "coplus") then |
---|
304 | igcm_coplus = iq |
---|
305 | mmol(igcm_coplus) = 28. |
---|
306 | count = count + 1 |
---|
307 | end if |
---|
308 | if (noms(iq) == "cplus") then |
---|
309 | igcm_cplus = iq |
---|
310 | mmol(igcm_cplus) = 12. |
---|
311 | count = count + 1 |
---|
312 | end if |
---|
313 | if (noms(iq) == "nplus") then |
---|
314 | igcm_nplus = iq |
---|
315 | mmol(igcm_nplus) = 14. |
---|
316 | count = count + 1 |
---|
317 | end if |
---|
318 | if (noms(iq) == "noplus") then |
---|
319 | igcm_noplus = iq |
---|
320 | mmol(igcm_noplus) = 30. |
---|
321 | count = count + 1 |
---|
322 | end if |
---|
323 | if (noms(iq) == "n2plus") then |
---|
324 | igcm_n2plus = iq |
---|
325 | mmol(igcm_n2plus) = 28. |
---|
326 | count = count + 1 |
---|
327 | end if |
---|
328 | if (noms(iq) == "hplus") then |
---|
329 | igcm_hplus = iq |
---|
330 | mmol(igcm_hplus) = 1. |
---|
331 | count = count + 1 |
---|
332 | end if |
---|
333 | if (noms(iq) == "hco2plus") then |
---|
334 | igcm_hco2plus = iq |
---|
335 | mmol(igcm_hco2plus) = 45. |
---|
336 | count = count + 1 |
---|
337 | end if |
---|
338 | if (noms(iq) == "hcoplus") then |
---|
339 | igcm_hcoplus = iq |
---|
340 | mmol(igcm_hcoplus) = 29. |
---|
341 | count = count + 1 |
---|
342 | end if |
---|
343 | if (noms(iq) == "h2oplus") then |
---|
344 | igcm_h2oplus = iq |
---|
345 | mmol(igcm_h2oplus) = 18. |
---|
346 | count = count + 1 |
---|
347 | end if |
---|
348 | if (noms(iq) == "h3oplus") then |
---|
349 | igcm_h3oplus = iq |
---|
350 | mmol(igcm_h3oplus) = 19. |
---|
351 | count = count + 1 |
---|
352 | end if |
---|
353 | if (noms(iq) == "ohplus") then |
---|
354 | igcm_ohplus = iq |
---|
355 | mmol(igcm_ohplus) = 17. |
---|
356 | count = count + 1 |
---|
357 | end if |
---|
358 | if (noms(iq) == "elec") then |
---|
359 | igcm_elec = iq |
---|
360 | mmol(igcm_elec) = 1./1822.89 |
---|
361 | count = count + 1 |
---|
362 | end if |
---|
363 | |
---|
364 | ! 1.5 find idealized non-condensible tracer |
---|
365 | |
---|
366 | if (noms(iq) == "Ar_N2") then |
---|
367 | igcm_ar_n2 = iq |
---|
368 | mmol(igcm_ar_n2) = 30. |
---|
369 | count = count + 1 |
---|
370 | end if |
---|
371 | |
---|
372 | end do ! of do iq=1,nqmx |
---|
373 | |
---|
374 | ! 1.6 check that we identified all tracers: |
---|
375 | |
---|
376 | if (count /= nqmx) then |
---|
377 | write(*,*) "inichim_newstart: found only ",count," tracers" |
---|
378 | write(*,*) " expected ",nqmx |
---|
379 | do iq = 1,count |
---|
380 | write(*,*) ' ', iq, ' ', trim(noms(iq)) |
---|
381 | end do |
---|
382 | call abort_physic("inichim_newstart","tracer mismatch",1) |
---|
383 | else |
---|
384 | write(*,*) "inichim_newstart: found all expected tracers" |
---|
385 | do iq = 1,nqmx |
---|
386 | write(*,*) ' ', iq, ' ', trim(noms(iq)) |
---|
387 | end do |
---|
388 | end if |
---|
389 | |
---|
390 | ! 1.7 check if nitrogen species are present: |
---|
391 | |
---|
392 | if(igcm_no == 0) then |
---|
393 | !check that no N species is in traceur.def |
---|
394 | if(igcm_n /= 0 .or. igcm_no2 /= 0 .or. igcm_n2d /= 0) then |
---|
395 | write(*,*)'inichim_newstart error:' |
---|
396 | write(*,*)'N, NO2 and/or N2D are in traceur.def, but not NO' |
---|
397 | write(*,*)'stop' |
---|
398 | call abort_physic("inichim_newstart","missing no tracer",1) |
---|
399 | endif |
---|
400 | flagnitro = .false. |
---|
401 | nspe = 14 |
---|
402 | else |
---|
403 | !check that all N species are in traceur.def |
---|
404 | if(igcm_n == 0 .or. igcm_no2 == 0 .or. igcm_n2d == 0) then |
---|
405 | write(*,*)'inichim_newstart error:' |
---|
406 | write(*,*)'if NO is in traceur.def, N, NO2 and N2D must also be' |
---|
407 | write(*,*)'stop' |
---|
408 | call abort_physic("inichim_newstart","missing n* tracer",1) |
---|
409 | endif |
---|
410 | flagnitro = .true. |
---|
411 | nspe = 18 |
---|
412 | endif |
---|
413 | |
---|
414 | ! 1.8 allocate arrays |
---|
415 | |
---|
416 | allocate(niq(nspe)) |
---|
417 | allocate(vmrinit(nalt,nspe)) |
---|
418 | allocate(vmrint(nspe)) |
---|
419 | |
---|
420 | ! 2. load in chemistry data for initialization: |
---|
421 | |
---|
422 | ! order of major species in initialization file: |
---|
423 | ! |
---|
424 | ! 1: co2 |
---|
425 | ! 2: ar |
---|
426 | ! 3: n2 |
---|
427 | ! 4: o2 |
---|
428 | ! 5: co |
---|
429 | ! 6: o |
---|
430 | ! 7: h2 |
---|
431 | ! |
---|
432 | ! order of minor species in initialization file: |
---|
433 | ! |
---|
434 | ! 1: h |
---|
435 | ! 2: oh |
---|
436 | ! 3: ho2 |
---|
437 | ! 4: h2o |
---|
438 | ! 5: h2o2 |
---|
439 | ! 6: o1d |
---|
440 | ! 7: o3 |
---|
441 | ! |
---|
442 | ! order of nitrogen species in initialization file: |
---|
443 | ! |
---|
444 | ! 1: n |
---|
445 | ! 2: no |
---|
446 | ! 3: no2 |
---|
447 | ! 4: n2d |
---|
448 | |
---|
449 | ! major species: |
---|
450 | |
---|
451 | niq(1) = igcm_co2 |
---|
452 | niq(2) = igcm_ar |
---|
453 | niq(3) = igcm_n2 |
---|
454 | niq(4) = igcm_o2 |
---|
455 | niq(5) = igcm_co |
---|
456 | niq(6) = igcm_o |
---|
457 | niq(7) = igcm_h2 |
---|
458 | |
---|
459 | ! minor species: |
---|
460 | |
---|
461 | niq(8) = igcm_h |
---|
462 | niq(9) = igcm_oh |
---|
463 | niq(10) = igcm_ho2 |
---|
464 | niq(11) = igcm_h2o_vap |
---|
465 | niq(12) = igcm_h2o2 |
---|
466 | niq(13) = igcm_o1d |
---|
467 | niq(14) = igcm_o3 |
---|
468 | |
---|
469 | ! nitrogen species: |
---|
470 | |
---|
471 | if (flagnitro) then |
---|
472 | niq(15) = igcm_n |
---|
473 | niq(16) = igcm_no |
---|
474 | niq(17) = igcm_no2 |
---|
475 | niq(18) = igcm_n2d |
---|
476 | end if |
---|
477 | |
---|
478 | ! 2.1 open initialization files |
---|
479 | |
---|
480 | open(210, iostat=ierr,file=trim(datadir)//'/atmosfera_LMD_may.dat') |
---|
481 | if (ierr /= 0) then |
---|
482 | write(*,*)'Error : cannot open file atmosfera_LMD_may.dat ' |
---|
483 | write(*,*)'(in aeronomars/inichim_newstart.F)' |
---|
484 | write(*,*)'It should be in :', trim(datadir),'/' |
---|
485 | write(*,*)'1) You can change this path in callphys.def with' |
---|
486 | write(*,*)' datadir=/path/to/datafiles/' |
---|
487 | write(*,*)'2) If necessary atmosfera_LMD_may.dat (and others)' |
---|
488 | write(*,*)' can be obtained online on:' |
---|
489 | write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir' |
---|
490 | call abort_physic("inichim_newstart","missing input file",1) |
---|
491 | end if |
---|
492 | open(220, iostat=ierr,file=trim(datadir)//'/atmosfera_LMD_min.dat') |
---|
493 | if (ierr /= 0) then |
---|
494 | write(*,*)'Error : cannot open file atmosfera_LMD_min.dat ' |
---|
495 | write(*,*)'(in aeronomars/inichim_newstart.F)' |
---|
496 | write(*,*)'It should be in :', trim(datadir),'/' |
---|
497 | write(*,*)'1) You can change this path in callphys.def with' |
---|
498 | write(*,*)' datadir=/path/to/datafiles/' |
---|
499 | write(*,*)'2) If necessary atmosfera_LMD_min.dat (and others)' |
---|
500 | write(*,*)' can be obtained online on:' |
---|
501 | write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir' |
---|
502 | call abort_physic("inichim_newstart","missing input file",1) |
---|
503 | end if |
---|
504 | if(flagnitro) then |
---|
505 | open(230, iostat=ierr,file=trim(datadir)//'/atmosfera_LMD_nitr.dat') |
---|
506 | if (ierr.ne.0) then |
---|
507 | write(*,*)'Error : cannot open file atmosfera_LMD_nitr.dat ' |
---|
508 | write(*,*)'(in aeronomars/inichim_newstart.F)' |
---|
509 | write(*,*)'It should be in :', trim(datadir),'/' |
---|
510 | write(*,*)'1) You can change this path in callphys.def with' |
---|
511 | write(*,*)' datadir=/path/to/datafiles/' |
---|
512 | write(*,*)'2) If necessary atmosfera_LMD_nitr.dat (and others)' |
---|
513 | write(*,*)' can be obtained online on:' |
---|
514 | write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir' |
---|
515 | call abort_physic("inichim_newstart","missing input file",1) |
---|
516 | endif |
---|
517 | endif ! Of if(flagnitro) |
---|
518 | |
---|
519 | ! 2.2 read initialization files |
---|
520 | |
---|
521 | ! major species |
---|
522 | |
---|
523 | read(210,*) |
---|
524 | do l = 1,nalt |
---|
525 | read(210,*) dummy, tinit(l), pinit(l), densinit(l), & |
---|
526 | (vmrinit(l,n), n = 1,7) |
---|
527 | pinit(l) = pinit(l)*100. ! conversion in Pa |
---|
528 | pinit(l) = log(pinit(l)) ! for the vertical interpolation |
---|
529 | end do |
---|
530 | close(210) |
---|
531 | |
---|
532 | ! minor species |
---|
533 | |
---|
534 | read(220,*) |
---|
535 | do l = 1,nalt |
---|
536 | read(220,*) dummy, (vmrinit(l,n), n = 8,14) |
---|
537 | end do |
---|
538 | close(220) |
---|
539 | |
---|
540 | ! nitrogen species |
---|
541 | |
---|
542 | if (flagnitro) then |
---|
543 | read(230,*) |
---|
544 | do l = 1,nalt |
---|
545 | read(230,*) dummy, (vmrinit(l,n), n = 15,18) |
---|
546 | end do |
---|
547 | close(230) |
---|
548 | end if |
---|
549 | |
---|
550 | ! 3. initialization of tracers |
---|
551 | |
---|
552 | do i = 1,nbp_lon+1 |
---|
553 | do j = 1,nbp_lat |
---|
554 | do l = 1,nbp_lev |
---|
555 | |
---|
556 | pgcm = aps(l) + bps(l)*ps(i,j) ! gcm pressure |
---|
557 | pgcm = log(pgcm) ! for the vertical interpolation |
---|
558 | mmean(i,j,l) = 0. |
---|
559 | |
---|
560 | ! 3.1 vertical interpolation |
---|
561 | |
---|
562 | do n = 1,nspe |
---|
563 | call intrplf(pgcm,vmr,pinit,vmrinit(:,n),nalt) |
---|
564 | vmrint(n) = vmr |
---|
565 | iq = niq(n) |
---|
566 | mmean(i,j,l) = mmean(i,j,l) + vmrint(n)*mmol(iq) |
---|
567 | end do |
---|
568 | |
---|
569 | ! 3.2 attribute mixing ratio: - all layers or only thermosphere |
---|
570 | ! - with our without h2o |
---|
571 | |
---|
572 | if (flagthermo == 0 .or. (flagthermo == 1 .and. exp(pgcm) < 0.1)) then |
---|
573 | do n = 1,nspe |
---|
574 | iq = niq(n) |
---|
575 | if (iq /= igcm_h2o_vap .or. flagh2o == 1) then |
---|
576 | pq(i,j,l,iq) = vmrint(n)*mmol(iq)/mmean(i,j,l) |
---|
577 | end if |
---|
578 | end do |
---|
579 | end if |
---|
580 | |
---|
581 | end do |
---|
582 | end do |
---|
583 | end do |
---|
584 | |
---|
585 | ! set surface values of chemistry tracers to zero |
---|
586 | |
---|
587 | if (flagthermo == 0) then |
---|
588 | ! NB: no problem for "surface water vapour" tracer which is always 0 |
---|
589 | do n = 1,nspe |
---|
590 | iq = niq(n) |
---|
591 | qsurf(1:ngrid,iq) = 0. |
---|
592 | end do |
---|
593 | end if |
---|
594 | |
---|
595 | ! 3.3 initialization of tracers not contained in the initialization files |
---|
596 | |
---|
597 | ! methane : 10 ppbv |
---|
598 | |
---|
599 | if (igcm_ch4 /= 0) then |
---|
600 | vmr = 10.e-9 |
---|
601 | do i = 1,nbp_lon+1 |
---|
602 | do j = 1,nbp_lat |
---|
603 | do l = 1,nbp_lev |
---|
604 | pq(i,j,l,igcm_ch4) = vmr*mmol(igcm_ch4)/mmean(i,j,l) |
---|
605 | end do |
---|
606 | end do |
---|
607 | end do |
---|
608 | ! set surface value to zero |
---|
609 | qsurf(1:ngrid,igcm_ch4) = 0. |
---|
610 | end if |
---|
611 | |
---|
612 | ! ions: 0 |
---|
613 | |
---|
614 | if (igcm_co2plus /= 0) then |
---|
615 | !check that all required ions are in traceur.def |
---|
616 | if (igcm_o2plus == 0 .or. igcm_oplus == 0 .or. igcm_coplus == 0 & |
---|
617 | .or. igcm_cplus == 0 .or. igcm_nplus == 0 .or. igcm_noplus == 0 & |
---|
618 | .or. igcm_n2plus == 0 .or. igcm_hplus == 0 .or. igcm_hco2plus == 0 & |
---|
619 | .or. igcm_hcoplus == 0 .or. igcm_h2oplus == 0 & |
---|
620 | .or. igcm_h3oplus == 0.or. igcm_ohplus == 0 .or.igcm_elec == 0) then |
---|
621 | write(*,*)'inichim_newstart error:' |
---|
622 | write(*,*)'if co2plus is in traceur.def, all other ions must also be' |
---|
623 | write(*,*)'o2plus, oplus, coplus, cplus, nplus, noplus, n2plus' |
---|
624 | write(*,*)'hplus, hco2plus, hcoplus, h2oplus, h3oplus, ohplus, and elec' |
---|
625 | write(*,*)'stop' |
---|
626 | call abort_physic("inichim_newstart","missing ions in tracers",1) |
---|
627 | end if |
---|
628 | |
---|
629 | do i = 1,nbp_lon+1 |
---|
630 | do j = 1,nbp_lat |
---|
631 | do l = 1,nbp_lev |
---|
632 | ! all ions to 0 |
---|
633 | pq(i,j,l,igcm_co2plus) = 0. |
---|
634 | pq(i,j,l,igcm_o2plus) = 0. |
---|
635 | pq(i,j,l,igcm_oplus) = 0. |
---|
636 | pq(i,j,l,igcm_coplus) = 0. |
---|
637 | pq(i,j,l,igcm_cplus) = 0. |
---|
638 | pq(i,j,l,igcm_nplus) = 0. |
---|
639 | pq(i,j,l,igcm_noplus) = 0. |
---|
640 | pq(i,j,l,igcm_n2plus) = 0. |
---|
641 | pq(i,j,l,igcm_hplus) = 0. |
---|
642 | pq(i,j,l,igcm_hco2plus) = 0. |
---|
643 | pq(i,j,l,igcm_hcoplus) = 0. |
---|
644 | pq(i,j,l,igcm_h2oplus) = 0. |
---|
645 | pq(i,j,l,igcm_h3oplus) = 0. |
---|
646 | pq(i,j,l,igcm_ohplus) = 0. |
---|
647 | pq(i,j,l,igcm_elec) = 0. |
---|
648 | end do |
---|
649 | end do |
---|
650 | end do |
---|
651 | |
---|
652 | ! surface value to 0 |
---|
653 | |
---|
654 | qsurf(1:ngrid,igcm_co2plus) = 0. |
---|
655 | qsurf(1:ngrid,igcm_o2plus) = 0. |
---|
656 | qsurf(1:ngrid,igcm_oplus) = 0. |
---|
657 | qsurf(1:ngrid,igcm_coplus) = 0. |
---|
658 | qsurf(1:ngrid,igcm_cplus) = 0. |
---|
659 | qsurf(1:ngrid,igcm_nplus) = 0. |
---|
660 | qsurf(1:ngrid,igcm_noplus) = 0. |
---|
661 | qsurf(1:ngrid,igcm_n2plus) = 0. |
---|
662 | qsurf(1:ngrid,igcm_hplus) = 0. |
---|
663 | qsurf(1:ngrid,igcm_hco2plus) = 0. |
---|
664 | qsurf(1:ngrid,igcm_hcoplus) = 0. |
---|
665 | qsurf(1:ngrid,igcm_h2oplus) = 0. |
---|
666 | qsurf(1:ngrid,igcm_h3oplus) = 0. |
---|
667 | qsurf(1:ngrid,igcm_ohplus) = 0. |
---|
668 | qsurf(1:ngrid,igcm_elec) = 0. |
---|
669 | |
---|
670 | else |
---|
671 | |
---|
672 | if (igcm_o2plus /= 0 .or. igcm_oplus /= 0 .or. igcm_coplus /= 0 & |
---|
673 | .or. igcm_cplus /= 0 .or. igcm_nplus /= 0 .or. igcm_noplus /= 0 & |
---|
674 | .or. igcm_n2plus /= 0 .or. igcm_hplus /= 0 .or. igcm_hco2plus /= 0 & |
---|
675 | .or. igcm_hcoplus /= 0 .or. igcm_h2oplus /= 0 & |
---|
676 | .or. igcm_h3oplus /= 0 .or. igcm_ohplus /= 0 .or. igcm_elec /= 0) then |
---|
677 | write(*,*)'inichim_newstart error:' |
---|
678 | write(*,*)'some ions are in traceur.def, but not co2plus' |
---|
679 | write(*,*)'stop' |
---|
680 | call abort_physic("inichim_newstart","missing ions in tracers",1) |
---|
681 | end if |
---|
682 | end if ! of if(igcm_co2 /= 0) |
---|
683 | |
---|
684 | ! deallocations |
---|
685 | |
---|
686 | deallocate(niq) |
---|
687 | deallocate(vmrinit) |
---|
688 | deallocate(vmrint) |
---|
689 | |
---|
690 | end |
---|