1 | subroutine inichim_newstart(pq,ps,sig,nq,latfi,lonfi,airefi, |
---|
2 | $ thermo,qsurf) |
---|
3 | |
---|
4 | implicit none |
---|
5 | |
---|
6 | c======================================================================= |
---|
7 | c |
---|
8 | c subject: |
---|
9 | c -------- |
---|
10 | c |
---|
11 | c Initialization of the composition for use in the building of a new start file |
---|
12 | c used by program newstart.F |
---|
13 | c also used by program testphys1d.F |
---|
14 | c |
---|
15 | c Author: Sebastien Lebonnois (08/11/2002) |
---|
16 | c ------- |
---|
17 | c Revised 07/2003 by Monica Angelats-i-Coll to use input files |
---|
18 | c Modified 10/2008 Identify tracers by their names Ehouarn Millour |
---|
19 | c |
---|
20 | c Arguments: |
---|
21 | c ---------- |
---|
22 | c |
---|
23 | c pq(iip1,jjp1,llm,nqmx) Advected fields, ie chemical species here |
---|
24 | c sig = sigma vertical coordinate (interface of layers) |
---|
25 | c nq: number of tracers to initialize (used to evaluate if only |
---|
26 | c chemistry species are to be initialized or if water vapour |
---|
27 | c and water ice should also be initialized. |
---|
28 | c NB: number of chemistry tracers is ncomp (set in chimiedata.h) |
---|
29 | c======================================================================= |
---|
30 | |
---|
31 | c Declarations : |
---|
32 | c -------------- |
---|
33 | |
---|
34 | #include "dimensions.h" |
---|
35 | #include "dimphys.h" |
---|
36 | #include "paramet.h" |
---|
37 | #include "chimiedata.h" |
---|
38 | #include "tracer.h" |
---|
39 | #include "comcstfi.h" |
---|
40 | #include "comdiurn.h" |
---|
41 | #include "callkeys.h" |
---|
42 | #include "temps.h" |
---|
43 | #include "datafile.h" |
---|
44 | |
---|
45 | c Arguments : |
---|
46 | c ----------- |
---|
47 | |
---|
48 | real ps(iip1,jjp1) |
---|
49 | real pq(iip1,jjp1,llm,nqmx) ! Advected fields, ie chemical species |
---|
50 | real sig(llm+1) ! vertical coordinate (interface of layers) |
---|
51 | integer nq ! =nqmx |
---|
52 | ! or =nqmx-1 if H2O kept |
---|
53 | ! or =nqmx-2 if H2O and ice kept |
---|
54 | REAL latfi(ngridmx),lonfi(ngridmx),airefi(ngridmx) |
---|
55 | integer thermo ! flag for thermosphere initialisation only |
---|
56 | real :: qsurf(ngridmx,nqmx) ! surface tracers |
---|
57 | |
---|
58 | c Local variables : |
---|
59 | c ----------------- |
---|
60 | |
---|
61 | integer iq,i,j,l,n, nspecini,inil,iqmax,iiq |
---|
62 | INTEGER aa(nqmx) |
---|
63 | REAL qtot(iip1,jjp1,llm) |
---|
64 | real densconc,zgcm,zfile(252),vmrint,nt, mmean |
---|
65 | real nxf(252),zfilemin(252),sigfile(252) |
---|
66 | real vmrf(252,14),nf(252,14) |
---|
67 | real tfile(252),zzfile(252) |
---|
68 | real ni(14),niold(14) |
---|
69 | real mmolf(14) ! molecular mass in amu (species in files) |
---|
70 | data mmolf/44.,40.,28.,32.,28.,16.,2., ! majors |
---|
71 | & 1.,17.,33.,18.,34.,16.,48./ ! minors |
---|
72 | character*3 tmp ! temporary variable |
---|
73 | integer ierr,lnblnk |
---|
74 | external lnblnk |
---|
75 | |
---|
76 | logical :: oldnames ! =.true. if old tracer naming convention (q01,...) |
---|
77 | character(len=20) :: txt ! to store some text |
---|
78 | integer :: nqchem_start,count |
---|
79 | integer :: nqchem(nqmx) ! local chemistry tracer index array |
---|
80 | integer :: nbqchem ! total number of chemical species (water included) |
---|
81 | |
---|
82 | c----------------------------------------------------------------------- |
---|
83 | c Input/Output |
---|
84 | c ------------ |
---|
85 | ! read in 'callphys.def' |
---|
86 | call inichim_readcallphys() |
---|
87 | |
---|
88 | ! check if tracers have 'old' names |
---|
89 | oldnames=.false. |
---|
90 | count=0 |
---|
91 | do iq=1,nqmx |
---|
92 | txt=" " |
---|
93 | write(txt,'(a1,i2.2)') 'q',iq |
---|
94 | if (txt.eq.noms(iq)) then |
---|
95 | count=count+1 |
---|
96 | endif |
---|
97 | enddo ! of do iq=1,nqmx |
---|
98 | |
---|
99 | if (count.eq.nqmx) then |
---|
100 | write(*,*) "inichim_newstart: tracers seem to follow old ", |
---|
101 | & "naming convention (q01,q02,...)" |
---|
102 | write(*,*) " => will work for now ... " |
---|
103 | write(*,*) " but you should rename them" |
---|
104 | oldnames=.true. |
---|
105 | endif |
---|
106 | |
---|
107 | ! copy/set tracer names |
---|
108 | if (oldnames) then ! old names (derived from iq & option values) |
---|
109 | ! 1. dust: |
---|
110 | if (dustbin.ne.0) then ! transported dust |
---|
111 | do iq=1,dustbin |
---|
112 | txt=" " |
---|
113 | write(txt,'(a4,i2.2)') 'dust',iq |
---|
114 | noms(iq)=txt |
---|
115 | mmol(iq)=100. |
---|
116 | enddo |
---|
117 | ! special case if doubleq |
---|
118 | if (doubleq) then |
---|
119 | if (dustbin.ne.2) then |
---|
120 | write(*,*) 'initracer: error expected dustbin=2' |
---|
121 | else |
---|
122 | noms(1)='dust_mass' ! dust mass mixing ratio |
---|
123 | noms(2)='dust_number' ! dust number mixing ratio |
---|
124 | endif |
---|
125 | endif |
---|
126 | endif |
---|
127 | ! 2. water & ice |
---|
128 | if (water) then |
---|
129 | noms(nqmx)='h2o_vap' |
---|
130 | mmol(nqmx)=18. |
---|
131 | noms(nqmx-1)='h2o_ice' |
---|
132 | mmol(nqmx-1)=18. |
---|
133 | endif |
---|
134 | ! 3. Chemistry |
---|
135 | if (photochem .or. callthermos) then |
---|
136 | if (doubleq) then |
---|
137 | nqchem_start=3 |
---|
138 | else |
---|
139 | nqchem_start=dustbin+1 |
---|
140 | end if |
---|
141 | endif ! of if (photochem .or. callthermos) |
---|
142 | do iq=nqchem_start,nqchem_start+ncomp-1 |
---|
143 | noms(iq)=nomchem(iq-nqchem_start+1) |
---|
144 | mmol(iq)=mmolchem(iq-nqchem_start+1) |
---|
145 | enddo |
---|
146 | ! 4. Other tracers ???? |
---|
147 | if ((dustbin.eq.0).and.(.not.water)) then |
---|
148 | noms(1)='co2' |
---|
149 | mmol(1)=44 |
---|
150 | if (nqmx.eq.2) then |
---|
151 | noms(nqmx)='Ar_N2' |
---|
152 | mmol(nqmx)=30 |
---|
153 | endif |
---|
154 | endif |
---|
155 | else |
---|
156 | ! noms(:) already contain tracer names |
---|
157 | ! mmol(:) array is initialized later (see below) |
---|
158 | endif ! of if (oldnames) |
---|
159 | |
---|
160 | ! special modification when using 'old' tracers: |
---|
161 | ! qsurf(nqmx) was h2o ice whereas q(nqmx) was water vapour |
---|
162 | ! and (if iceparty) q(nqmx-1) was null whereas q(nqmx-1) was water ice |
---|
163 | if (oldnames.and.water) then |
---|
164 | write(*,*)'inichim_newstart: moving surface water ice to index ' |
---|
165 | & ,nqmx-1 |
---|
166 | qsurf(1:ngridmx,nqmx-1)=qsurf(1:ngridmx,nqmx) |
---|
167 | qsurf(1:ngridmx,nqmx)=0 |
---|
168 | endif |
---|
169 | |
---|
170 | c Dimension test: |
---|
171 | c --------------- |
---|
172 | |
---|
173 | if (water) then |
---|
174 | ! if ((nqchem_min+ncomp+1).ne.nqmx) then |
---|
175 | if ((dustbin+ncomp+2).ne.nqmx) then |
---|
176 | ! print*,'********* Dimension problem! ********' |
---|
177 | ! print*,"nqchem_min+ncomp+1).ne.nqmx" |
---|
178 | print*,'inichim_newstart: tracer dimension problem:' |
---|
179 | print*,"(dustbin+ncomp+2).ne.nqmx" |
---|
180 | print*,"ncomp: ",ncomp |
---|
181 | ! print*,"nqchem_min: ",nqchem_min |
---|
182 | print*,"nqmx: ",nqmx |
---|
183 | print*,'Change ncomp in chimiedata.h' |
---|
184 | STOP |
---|
185 | endif |
---|
186 | else |
---|
187 | ! if ((nqchem_min+ncomp).ne.nqmx) then |
---|
188 | if ((dustbin+ncomp+1).ne.nqmx) then |
---|
189 | ! print*,'********* Dimension problem! ********' |
---|
190 | ! print*,"nqchem_min+ncomp).ne.nqmx" |
---|
191 | print*,'initracer: tracer dimension problem:' |
---|
192 | print*,"(dustbin+ncomp+1).ne.nqmx" |
---|
193 | print*,"ncomp: ",ncomp |
---|
194 | ! print*,"nqchem_min: ",nqchem_min |
---|
195 | print*,"nqmx: ",nqmx |
---|
196 | print*,'Change ncomp in chimiedata.h' |
---|
197 | STOP |
---|
198 | endif |
---|
199 | endif |
---|
200 | |
---|
201 | c noms and mmol vectors: |
---|
202 | c ---------------------- |
---|
203 | ! |
---|
204 | ! if (iceparty) then |
---|
205 | ! do iq=nqchem_min,nqmx-2 |
---|
206 | ! noms(iq) = nomchem(iq-nqchem_min+1) |
---|
207 | ! mmol(iq) = mmolchem(iq-nqchem_min+1) |
---|
208 | ! enddo |
---|
209 | ! noms(nqmx-1) = 'ice' |
---|
210 | ! mmol(nqmx-1) = 18. |
---|
211 | ! noms(nqmx) = 'h2o' |
---|
212 | ! mmol(nqmx) = 18. |
---|
213 | ! else |
---|
214 | ! do iq=nqchem_min,nqmx-1 |
---|
215 | ! noms(iq) = nomchem(iq-nqchem_min+1) |
---|
216 | ! mmol(iq) = mmolchem(iq-nqchem_min+1) |
---|
217 | ! enddo |
---|
218 | ! noms(nqmx) = 'h2o' |
---|
219 | ! mmol(nqmx) = 18. |
---|
220 | ! endif |
---|
221 | ! if (nqchem_min.gt.1) then |
---|
222 | ! do iq=1,nqchem_min-1 |
---|
223 | ! noms(iq) = 'dust' |
---|
224 | ! mmol(iq) = 100. |
---|
225 | ! enddo |
---|
226 | ! endif |
---|
227 | |
---|
228 | |
---|
229 | ! Identify tracers by their names: (and set corresponding values of mmol) |
---|
230 | ! 0. initialize tracer indexes to zero: |
---|
231 | do iq=1,nqmx |
---|
232 | igcm_dustbin(iq)=0 |
---|
233 | enddo |
---|
234 | igcm_dust_mass=0 |
---|
235 | igcm_dust_number=0 |
---|
236 | igcm_h2o_vap=0 |
---|
237 | igcm_h2o_ice=0 |
---|
238 | igcm_co2=0 |
---|
239 | igcm_co=0 |
---|
240 | igcm_o=0 |
---|
241 | igcm_o1d=0 |
---|
242 | igcm_o2=0 |
---|
243 | igcm_o3=0 |
---|
244 | igcm_h=0 |
---|
245 | igcm_h2=0 |
---|
246 | igcm_oh=0 |
---|
247 | igcm_ho2=0 |
---|
248 | igcm_h2o2=0 |
---|
249 | igcm_n2=0 |
---|
250 | igcm_ar=0 |
---|
251 | igcm_ar_n2=0 |
---|
252 | |
---|
253 | ! 1. find dust tracers |
---|
254 | count=0 |
---|
255 | if (dustbin.gt.0) then |
---|
256 | do iq=1,nqmx |
---|
257 | txt=" " |
---|
258 | write(txt,'(a4,i2.2)')'dust',count+1 |
---|
259 | if (noms(iq).eq.txt) then |
---|
260 | count=count+1 |
---|
261 | igcm_dustbin(count)=iq |
---|
262 | mmol(iq)=100. |
---|
263 | endif |
---|
264 | enddo !do iq=1,nqmx |
---|
265 | endif ! of if (dustbin.gt.0) |
---|
266 | if (doubleq) then |
---|
267 | do iq=1,nqmx |
---|
268 | if (noms(iq).eq."dust_mass") then |
---|
269 | igcm_dust_mass=iq |
---|
270 | count=count+1 |
---|
271 | endif |
---|
272 | if (noms(iq).eq."dust_number") then |
---|
273 | igcm_dust_number=iq |
---|
274 | count=count+1 |
---|
275 | endif |
---|
276 | enddo |
---|
277 | endif ! of if (doubleq) |
---|
278 | ! 2. find chemistry and water tracers |
---|
279 | nbqchem=0 |
---|
280 | do iq=1,nqmx |
---|
281 | if (noms(iq).eq."co2") then |
---|
282 | igcm_co2=iq |
---|
283 | mmol(igcm_co2)=44. |
---|
284 | count=count+1 |
---|
285 | nbqchem=nbqchem+1 |
---|
286 | endif |
---|
287 | if (noms(iq).eq."co") then |
---|
288 | igcm_co=iq |
---|
289 | mmol(igcm_co)=28. |
---|
290 | count=count+1 |
---|
291 | nbqchem=nbqchem+1 |
---|
292 | endif |
---|
293 | if (noms(iq).eq."o") then |
---|
294 | igcm_o=iq |
---|
295 | mmol(igcm_o)=16. |
---|
296 | count=count+1 |
---|
297 | nbqchem=nbqchem+1 |
---|
298 | endif |
---|
299 | if (noms(iq).eq."o1d") then |
---|
300 | igcm_o1d=iq |
---|
301 | mmol(igcm_o1d)=16. |
---|
302 | count=count+1 |
---|
303 | nbqchem=nbqchem+1 |
---|
304 | endif |
---|
305 | if (noms(iq).eq."o2") then |
---|
306 | igcm_o2=iq |
---|
307 | mmol(igcm_o2)=32. |
---|
308 | count=count+1 |
---|
309 | nbqchem=nbqchem+1 |
---|
310 | endif |
---|
311 | if (noms(iq).eq."o3") then |
---|
312 | igcm_o3=iq |
---|
313 | mmol(igcm_o3)=48. |
---|
314 | count=count+1 |
---|
315 | nbqchem=nbqchem+1 |
---|
316 | endif |
---|
317 | if (noms(iq).eq."h") then |
---|
318 | igcm_h=iq |
---|
319 | mmol(igcm_h)=1. |
---|
320 | count=count+1 |
---|
321 | nbqchem=nbqchem+1 |
---|
322 | endif |
---|
323 | if (noms(iq).eq."h2") then |
---|
324 | igcm_h2=iq |
---|
325 | mmol(igcm_h2)=2. |
---|
326 | count=count+1 |
---|
327 | nbqchem=nbqchem+1 |
---|
328 | endif |
---|
329 | if (noms(iq).eq."oh") then |
---|
330 | igcm_oh=iq |
---|
331 | mmol(igcm_oh)=17. |
---|
332 | count=count+1 |
---|
333 | nbqchem=nbqchem+1 |
---|
334 | endif |
---|
335 | if (noms(iq).eq."ho2") then |
---|
336 | igcm_ho2=iq |
---|
337 | mmol(igcm_ho2)=33. |
---|
338 | count=count+1 |
---|
339 | nbqchem=nbqchem+1 |
---|
340 | endif |
---|
341 | if (noms(iq).eq."h2o2") then |
---|
342 | igcm_h2o2=iq |
---|
343 | mmol(igcm_h2o2)=34. |
---|
344 | count=count+1 |
---|
345 | nbqchem=nbqchem+1 |
---|
346 | endif |
---|
347 | if (noms(iq).eq."n2") then |
---|
348 | igcm_n2=iq |
---|
349 | mmol(igcm_n2)=28. |
---|
350 | count=count+1 |
---|
351 | nbqchem=nbqchem+1 |
---|
352 | endif |
---|
353 | if (noms(iq).eq."ar") then |
---|
354 | igcm_ar=iq |
---|
355 | mmol(igcm_ar)=40. |
---|
356 | count=count+1 |
---|
357 | nbqchem=nbqchem+1 |
---|
358 | endif |
---|
359 | if (noms(iq).eq."h2o_vap") then |
---|
360 | igcm_h2o_vap=iq |
---|
361 | mmol(igcm_h2o_vap)=18. |
---|
362 | count=count+1 |
---|
363 | nbqchem=nbqchem+1 |
---|
364 | endif |
---|
365 | if (noms(iq).eq."h2o_ice") then |
---|
366 | igcm_h2o_ice=iq |
---|
367 | mmol(igcm_h2o_ice)=18. |
---|
368 | count=count+1 |
---|
369 | nbqchem=nbqchem+1 |
---|
370 | endif |
---|
371 | ! Other stuff: e.g. for simulations using co2 + neutral gaz |
---|
372 | if (noms(iq).eq."Ar_N2") then |
---|
373 | igcm_ar_n2=iq |
---|
374 | mmol(igcm_ar_n2)=30. |
---|
375 | count=count+1 |
---|
376 | endif |
---|
377 | enddo ! of do iq=1,nqmx |
---|
378 | |
---|
379 | ! check that we identified all tracers: |
---|
380 | if (count.ne.nqmx) then |
---|
381 | write(*,*) "inichim_newstart: found only ",count," tracers" |
---|
382 | write(*,*) " expected ",nqmx |
---|
383 | do iq=1,count |
---|
384 | write(*,*)' ',iq,' ',trim(noms(iq)) |
---|
385 | enddo |
---|
386 | stop |
---|
387 | else |
---|
388 | write(*,*)"inichim_newstart: found all expected tracers, namely:" |
---|
389 | do iq=1,nqmx |
---|
390 | write(*,*)' ',iq,' ',trim(noms(iq)) |
---|
391 | enddo |
---|
392 | endif |
---|
393 | |
---|
394 | ! build local chemistry tracer index array nqchem(:) |
---|
395 | ! as in the old days, water vapour is last and water ice, |
---|
396 | ! if present, is just before water vapour |
---|
397 | nqchem(1)=igcm_co2 |
---|
398 | nqchem(2)=igcm_co |
---|
399 | nqchem(3)=igcm_o |
---|
400 | nqchem(4)=igcm_o1d |
---|
401 | nqchem(5)=igcm_o2 |
---|
402 | nqchem(6)=igcm_o3 |
---|
403 | nqchem(7)=igcm_h |
---|
404 | nqchem(8)=igcm_h2 |
---|
405 | nqchem(9)=igcm_oh |
---|
406 | nqchem(10)=igcm_ho2 |
---|
407 | nqchem(11)=igcm_h2o2 |
---|
408 | nqchem(12)=igcm_n2 |
---|
409 | nqchem(13)=igcm_ar |
---|
410 | nqchem(14)=igcm_h2o_ice |
---|
411 | nqchem(15)=igcm_h2o_vap |
---|
412 | |
---|
413 | ! Load in chemistry data for initialization: |
---|
414 | c---------------------------------------------------------------------- |
---|
415 | c Order of Major species in file: |
---|
416 | c |
---|
417 | c 1=CO2 |
---|
418 | c 2=AR |
---|
419 | c 3=N2 |
---|
420 | c 4=O2 |
---|
421 | c 5=CO |
---|
422 | c 6=O |
---|
423 | c 7=H2 |
---|
424 | c |
---|
425 | c Order of Minor species in files: |
---|
426 | c |
---|
427 | c 1=H |
---|
428 | c 2=OH |
---|
429 | c 3=HO2 |
---|
430 | c 4=H2O |
---|
431 | c 5=H2O2 |
---|
432 | c 6=O1D |
---|
433 | c 7=O3 |
---|
434 | c |
---|
435 | c---------------------------------------------------------------------- |
---|
436 | c Major species: |
---|
437 | aa(igcm_co2) = 1 |
---|
438 | aa(igcm_ar) = 2 |
---|
439 | aa(igcm_n2) = 3 |
---|
440 | aa(igcm_o2) = 4 |
---|
441 | aa(igcm_co) = 5 |
---|
442 | aa(igcm_o) = 6 |
---|
443 | aa(igcm_h2) = 7 |
---|
444 | c Minor species: |
---|
445 | aa(igcm_h) = 8 |
---|
446 | aa(igcm_oh) = 9 |
---|
447 | aa(igcm_ho2) = 10 |
---|
448 | aa(igcm_h2o_vap) = 11 |
---|
449 | aa(igcm_h2o2) = 12 |
---|
450 | aa(igcm_o1d) = 13 |
---|
451 | aa(igcm_o3) = 14 |
---|
452 | c---------------------------------------------------------------------- |
---|
453 | open(210, iostat=ierr,file= |
---|
454 | & datafile(1:lnblnk(datafile))//'/atmosfera_LMD_may.dat') |
---|
455 | if (ierr.ne.0) then |
---|
456 | write(*,*)'Error : cannot open file atmosfera_LMD_may.dat ' |
---|
457 | write(*,*)'(in aeronomars/inichim.F)' |
---|
458 | write(*,*)'It should be in :', datafile(1:lnblnk(datafile)),'/' |
---|
459 | write(*,*)'1) You can change this directory address in ' |
---|
460 | write(*,*)' file phymars/datafile.h' |
---|
461 | write(*,*)'2) If necessary atmosfera_LMD_may.dat (and others)' |
---|
462 | write(*,*)' can be obtained online on:' |
---|
463 | write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile' |
---|
464 | STOP |
---|
465 | endif |
---|
466 | open(220, iostat=ierr,file= |
---|
467 | & datafile(1:lnblnk(datafile))//'/atmosfera_LMD_min.dat') |
---|
468 | if (ierr.ne.0) then |
---|
469 | write(*,*)'Error : cannot open file atmosfera_LMD_min.dat ' |
---|
470 | write(*,*)'(in aeronomars/inichim.F)' |
---|
471 | write(*,*)'It should be in :', datafile(1:lnblnk(datafile)),'/' |
---|
472 | write(*,*)'1) You can change this directory address in ' |
---|
473 | write(*,*)' file phymars/datafile.h' |
---|
474 | write(*,*)'2) If necessary atmosfera_LMD_min.dat (and others)' |
---|
475 | write(*,*)' can be obtained online on:' |
---|
476 | write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile' |
---|
477 | STOP |
---|
478 | endif |
---|
479 | |
---|
480 | read(210,*) tmp |
---|
481 | read(220,*) tmp |
---|
482 | do l=1,252 |
---|
483 | read(210,112) tfile(l),zfile(l),nxf(l),(vmrf(l,n),n=1,7) |
---|
484 | zfile(l)=zfile(l)*100. !pressure (Pa) |
---|
485 | sigfile(l)=zfile(l)/zfile(1) |
---|
486 | enddo |
---|
487 | do l=1,252 |
---|
488 | read(220,113)zfilemin(l),(vmrf(l,n),n=8,14) |
---|
489 | zfilemin(l)=zfilemin(l)*1000. !height (m) |
---|
490 | do n=1,14 |
---|
491 | nf(l,n)=vmrf(l,n)*nxf(l) |
---|
492 | enddo |
---|
493 | enddo |
---|
494 | close(210) |
---|
495 | close(220) |
---|
496 | |
---|
497 | c flag thermo set to 1 or 0 in newstart |
---|
498 | c inil=33 for initialisation of thermosphere only |
---|
499 | c inil=1 for initialisation of all layers |
---|
500 | if (thermo.eq.1) then |
---|
501 | inil=33 |
---|
502 | else |
---|
503 | inil=1 |
---|
504 | endif |
---|
505 | |
---|
506 | ! Initialization of tracers |
---|
507 | |
---|
508 | do i=1,iip1 |
---|
509 | do j=1,jjp1 |
---|
510 | do l=inil,llm |
---|
511 | |
---|
512 | c zgcm=sig(l) |
---|
513 | zgcm=sig(l)*ps(i,j) |
---|
514 | densconc=0. |
---|
515 | nt=0. |
---|
516 | |
---|
517 | do n=1,14 |
---|
518 | call intrplf(zgcm,vmrint,zfile,nf(1,n),252) |
---|
519 | c call intrplf(zgcm,vmrint,sigfile,nf(1,n),252) |
---|
520 | ni(n)=vmrint |
---|
521 | |
---|
522 | densconc=densconc+ni(n)*mmolf(n) |
---|
523 | nt=nt+ni(n) |
---|
524 | enddo |
---|
525 | |
---|
526 | mmean=densconc/nt ! in amu |
---|
527 | |
---|
528 | if (nq.ne.nbqchem) then ! don't initialize water vapour |
---|
529 | do n=1,nq-1 |
---|
530 | pq(i,j,l,nqchem(n))= |
---|
531 | & ni(aa(nqchem(n)))/nt*mmol(nqchem(n))/mmean |
---|
532 | if(i.eq.iip1) pq(i,j,l,nqchem(n))=pq(1,j,l,nqchem(n)) |
---|
533 | enddo |
---|
534 | if (water .and. (nq.gt.nbqchem-2)) then |
---|
535 | pq(i,j,l,nqchem(nq)) = 0. |
---|
536 | if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq)) |
---|
537 | else |
---|
538 | pq(i,j,l,nqchem(nq))= |
---|
539 | & ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean |
---|
540 | if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq)) |
---|
541 | endif |
---|
542 | endif ! of if (nq.ne.nbqchem) |
---|
543 | |
---|
544 | if (nq.eq.nbqchem) then ! also initialize water vapour |
---|
545 | do n=1,nq-2 |
---|
546 | pq(i,j,l,nqchem(n))= |
---|
547 | & ni(aa(nqchem(n)))/nt*mmol(nqchem(n))/mmean |
---|
548 | if(i.eq.iip1) pq(i,j,l,nqchem(n))=pq(1,j,l,nqchem(n)) |
---|
549 | enddo |
---|
550 | if (water) then |
---|
551 | pq(i,j,l,nqchem(nq-1)) = 0. |
---|
552 | if(i.eq.iip1) then |
---|
553 | pq(i,j,l,nqchem(nq-1))=pq(1,j,l,nqchem(nq-1)) |
---|
554 | endif |
---|
555 | pq(i,j,l,nqchem(nq))= |
---|
556 | & ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean |
---|
557 | if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq)) |
---|
558 | else |
---|
559 | do n=nq-1,nq |
---|
560 | pq(i,j,l,nqchem(nq))= |
---|
561 | & ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean |
---|
562 | if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq)) |
---|
563 | enddo |
---|
564 | endif ! of if (water) |
---|
565 | endif ! of if (nq.eq.nqmx) |
---|
566 | enddo !nlayer |
---|
567 | enddo !ngrid_j |
---|
568 | enddo !ngrid_i |
---|
569 | |
---|
570 | cccccccccccccccccccccccccccccc |
---|
571 | c make sure that sum of q = 1 |
---|
572 | c dominent species is = 1 - sum(all other species) |
---|
573 | cccccccccccccccccccccccccccccc |
---|
574 | ! iqmax=nqchem_min |
---|
575 | iqmax=1 |
---|
576 | |
---|
577 | ! if ((nqmx-nqchem_min).gt.10) then |
---|
578 | if (nbqchem.gt.10) then |
---|
579 | do l=1,llm |
---|
580 | do j=1,jjp1 |
---|
581 | do i=1,iip1 |
---|
582 | ! do iq=nqchem_min,nqmx |
---|
583 | do iiq=1,nbqchem |
---|
584 | iq=nqchem(iiq) |
---|
585 | if ( (pq(i,j,l,iq).gt.pq(i,j,l,iqmax)).and. |
---|
586 | . (noms(iq).ne."ice") ) then |
---|
587 | iqmax=iq |
---|
588 | endif |
---|
589 | enddo |
---|
590 | pq(i,j,l,iqmax)=1. |
---|
591 | qtot(i,j,l)=0 |
---|
592 | ! do iq=nqchem_min,nqmx |
---|
593 | do iiq=1,nbqchem |
---|
594 | iq=nqchem(iiq) |
---|
595 | if ( (iq.ne.iqmax).and. |
---|
596 | . (noms(iq).ne."ice") ) then |
---|
597 | pq(i,j,l,iqmax)=pq(i,j,l,iqmax)-pq(i,j,l,iq) |
---|
598 | endif |
---|
599 | enddo !iq |
---|
600 | ! do iq=nqchem_min,nqmx |
---|
601 | do iiq=1,nbqchem |
---|
602 | iq=nqchem(iiq) |
---|
603 | if (noms(iq).ne."ice") then |
---|
604 | qtot(i,j,l)=qtot(i,j,l)+pq(i,j,l,iq) |
---|
605 | endif |
---|
606 | c if (i.eq.1.and.j.eq.1.and.l.eq.1) write(*,*) 'qtot(i,j,l)', |
---|
607 | c $ qtot(i,j,l) |
---|
608 | enddo !iq |
---|
609 | if (i.eq.1.and.j.eq.1) write(*,*) 'inichim_newstart: ', |
---|
610 | $ 'qtot(i,j,l)=',qtot(i,j,l) |
---|
611 | enddo !i |
---|
612 | enddo !j |
---|
613 | enddo !l |
---|
614 | endif ! of if (nbqchem.gt.10) |
---|
615 | ccccccccccccccccccccccccccccccc |
---|
616 | |
---|
617 | 66 format(i2,6(1x,e11.5)) |
---|
618 | 112 format(7x,f7.3,3x,e12.6,3x,e12.6,7(3x,e12.6)) |
---|
619 | 113 format(1x,f6.2,7(3x,e12.6)) |
---|
620 | |
---|
621 | end |
---|