1 | ! SPDX-FileCopyrightText: FLEXPART 1998-2019, see flexpart_license.txt |
---|
2 | ! SPDX-License-Identifier: GPL-3.0-or-later |
---|
3 | |
---|
4 | subroutine concoutput_inversion_nest(itime,outnum) |
---|
5 | ! i i |
---|
6 | !***************************************************************************** |
---|
7 | ! * |
---|
8 | ! Output of the concentration grid and the receptor concentrations. * |
---|
9 | ! * |
---|
10 | ! Author: A. Stohl * |
---|
11 | ! * |
---|
12 | ! 24 May 1995 * |
---|
13 | ! * |
---|
14 | ! 13 April 1999, Major update: if output size is smaller, dump output * |
---|
15 | ! in sparse matrix format; additional output of * |
---|
16 | ! uncertainty * |
---|
17 | ! * |
---|
18 | ! 05 April 2000, Major update: output of age classes; output for backward* |
---|
19 | ! runs is time spent in grid cell times total mass of * |
---|
20 | ! species. * |
---|
21 | ! * |
---|
22 | ! 17 February 2002, Appropriate dimensions for backward and forward runs * |
---|
23 | ! are now specified in file par_mod * |
---|
24 | ! * |
---|
25 | ! June 2006, write grid in sparse matrix with a single write command * |
---|
26 | ! in order to save disk space * |
---|
27 | ! * |
---|
28 | ! 2008 new sparse matrix format * |
---|
29 | ! |
---|
30 | ! January 2017, Separate files by release but include all timesteps * |
---|
31 | ! * |
---|
32 | !***************************************************************************** |
---|
33 | ! * |
---|
34 | ! Variables: * |
---|
35 | ! outnum number of samples * |
---|
36 | ! ncells number of cells with non-zero concentrations * |
---|
37 | ! sparse .true. if in sparse matrix format, else .false. * |
---|
38 | ! tot_mu 1 for forward, initial mass mixing ration for backw. runs * |
---|
39 | ! * |
---|
40 | !***************************************************************************** |
---|
41 | |
---|
42 | use unc_mod |
---|
43 | use point_mod |
---|
44 | use outg_mod |
---|
45 | use par_mod |
---|
46 | use com_mod |
---|
47 | use mean_mod |
---|
48 | |
---|
49 | implicit none |
---|
50 | |
---|
51 | real(kind=dp) :: jul |
---|
52 | integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss |
---|
53 | integer :: sp_count_i,sp_count_r |
---|
54 | real :: sp_fact |
---|
55 | real :: outnum,densityoutrecept(maxreceptor),xl,yl |
---|
56 | ! RLT |
---|
57 | real :: densitydryrecept(maxreceptor) |
---|
58 | real :: factor_dryrecept(maxreceptor) |
---|
59 | |
---|
60 | !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid), |
---|
61 | ! +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act, |
---|
62 | ! + maxageclass) |
---|
63 | !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act, |
---|
64 | ! + maxageclass) |
---|
65 | !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec, |
---|
66 | ! + maxpointspec_act,maxageclass) |
---|
67 | !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec, |
---|
68 | ! + maxpointspec_act,maxageclass), |
---|
69 | ! + drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec, |
---|
70 | ! + maxpointspec_act,maxageclass), |
---|
71 | ! + wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec, |
---|
72 | ! + maxpointspec_act,maxageclass) |
---|
73 | !real factor(0:numxgrid-1,0:numygrid-1,numzgrid) |
---|
74 | !real sparse_dump_r(numxgrid*numygrid*numzgrid) |
---|
75 | !integer sparse_dump_i(numxgrid*numygrid*numzgrid) |
---|
76 | |
---|
77 | !real sparse_dump_u(numxgrid*numygrid*numzgrid) |
---|
78 | real(dep_prec) :: auxgrid(nclassunc) |
---|
79 | real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act) |
---|
80 | real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled |
---|
81 | real,parameter :: weightair=28.97 |
---|
82 | logical :: sp_zer |
---|
83 | logical,save :: lnstart=.true. |
---|
84 | logical,save,allocatable,dimension(:) :: lnstartrel |
---|
85 | character :: adate*8,atime*6 |
---|
86 | character(len=3) :: anspec |
---|
87 | logical :: lexist |
---|
88 | character :: areldate*8,areltime*6 |
---|
89 | |
---|
90 | if(lnstart) then |
---|
91 | allocate(lnstartrel(maxpointspec_act)) |
---|
92 | lnstartrel(:)=.true. |
---|
93 | endif |
---|
94 | print*, 'lnstartrel = ',lnstartrel |
---|
95 | |
---|
96 | ! Determine current calendar date, needed for the file name |
---|
97 | !********************************************************** |
---|
98 | |
---|
99 | jul=bdate+real(itime,kind=dp)/86400._dp |
---|
100 | call caldate(jul,jjjjmmdd,ihmmss) |
---|
101 | write(adate,'(i8.8)') jjjjmmdd |
---|
102 | write(atime,'(i6.6)') ihmmss |
---|
103 | |
---|
104 | print*, 'outnum:',outnum |
---|
105 | print*, 'datetime:',adate//atime |
---|
106 | |
---|
107 | ! For forward simulations, output fields have dimension MAXSPEC, |
---|
108 | ! for backward simulations, output fields have dimension MAXPOINT. |
---|
109 | ! Thus, make loops either about nspec, or about numpoint |
---|
110 | !***************************************************************** |
---|
111 | |
---|
112 | |
---|
113 | if (ldirect.eq.1) then |
---|
114 | do ks=1,nspec |
---|
115 | do kp=1,maxpointspec_act |
---|
116 | tot_mu(ks,kp)=1 |
---|
117 | end do |
---|
118 | end do |
---|
119 | else |
---|
120 | do ks=1,nspec |
---|
121 | do kp=1,maxpointspec_act |
---|
122 | tot_mu(ks,kp)=xmass(kp,ks) |
---|
123 | end do |
---|
124 | end do |
---|
125 | endif |
---|
126 | |
---|
127 | |
---|
128 | !******************************************************************* |
---|
129 | ! Compute air density: sufficiently accurate to take it |
---|
130 | ! from coarse grid at some time |
---|
131 | ! Determine center altitude of output layer, and interpolate density |
---|
132 | ! data to that altitude |
---|
133 | !******************************************************************* |
---|
134 | |
---|
135 | do kz=1,numzgrid |
---|
136 | if (kz.eq.1) then |
---|
137 | halfheight=outheight(1)/2. |
---|
138 | else |
---|
139 | halfheight=(outheight(kz)+outheight(kz-1))/2. |
---|
140 | endif |
---|
141 | do kzz=2,nz |
---|
142 | if ((height(kzz-1).lt.halfheight).and. & |
---|
143 | (height(kzz).gt.halfheight)) goto 46 |
---|
144 | end do |
---|
145 | 46 kzz=max(min(kzz,nz),2) |
---|
146 | dz1=halfheight-height(kzz-1) |
---|
147 | dz2=height(kzz)-halfheight |
---|
148 | dz=dz1+dz2 |
---|
149 | do jy=0,numygridn-1 |
---|
150 | do ix=0,numxgridn-1 |
---|
151 | xl=outlon0n+real(ix)*dxoutn |
---|
152 | yl=outlat0n+real(jy)*dyoutn |
---|
153 | xl=(xl-xlon0)/dx |
---|
154 | yl=(yl-ylat0)/dy |
---|
155 | iix=max(min(nint(xl),nxmin1),0) |
---|
156 | jjy=max(min(nint(yl),nymin1),0) |
---|
157 | densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ & |
---|
158 | rho(iix,jjy,kzz-1,2)*dz2)/dz |
---|
159 | ! RLT |
---|
160 | densitydrygrid(ix,jy,kz)=(rho_dry(iix,jjy,kzz,2)*dz1+ & |
---|
161 | rho_dry(iix,jjy,kzz-1,2)*dz2)/dz |
---|
162 | end do |
---|
163 | end do |
---|
164 | end do |
---|
165 | |
---|
166 | do i=1,numreceptor |
---|
167 | xl=xreceptor(i) |
---|
168 | yl=yreceptor(i) |
---|
169 | iix=max(min(nint(xl),nxmin1),0) |
---|
170 | jjy=max(min(nint(yl),nymin1),0) |
---|
171 | densityoutrecept(i)=rho(iix,jjy,1,2) |
---|
172 | ! RLT |
---|
173 | densitydryrecept(i)=rho_dry(iix,jjy,1,2) |
---|
174 | end do |
---|
175 | |
---|
176 | ! RLT |
---|
177 | ! conversion factor for output relative to dry air |
---|
178 | factor_drygrid=densityoutgrid/densitydrygrid |
---|
179 | factor_dryrecept=densityoutrecept/densitydryrecept |
---|
180 | |
---|
181 | ! Output is different for forward and backward simulations |
---|
182 | do kz=1,numzgrid |
---|
183 | do jy=0,numygridn-1 |
---|
184 | do ix=0,numxgridn-1 |
---|
185 | if (ldirect.eq.1) then |
---|
186 | factor3d(ix,jy,kz)=1.e12/volumen(ix,jy,kz)/outnum |
---|
187 | else |
---|
188 | factor3d(ix,jy,kz)=real(abs(loutaver))/outnum |
---|
189 | endif |
---|
190 | end do |
---|
191 | end do |
---|
192 | end do |
---|
193 | |
---|
194 | !********************************************************************* |
---|
195 | ! Determine the standard deviation of the mean concentration or mixing |
---|
196 | ! ratio (uncertainty of the output) and the dry and wet deposition |
---|
197 | !********************************************************************* |
---|
198 | |
---|
199 | do ks=1,nspec |
---|
200 | |
---|
201 | write(anspec,'(i3.3)') ks |
---|
202 | |
---|
203 | do kp=1,maxpointspec_act |
---|
204 | |
---|
205 | print*, 'itime = ',itime |
---|
206 | print*, 'lage(1) = ',lage(1) |
---|
207 | print*, 'ireleasestart(kp) = ',ireleasestart(kp) |
---|
208 | print*, 'ireleaseend(kp) = ',ireleaseend(kp) |
---|
209 | |
---|
210 | ! check itime is within release and backward trajectory length |
---|
211 | if (nageclass.eq.1) then |
---|
212 | if ((itime.gt.ireleaseend(kp)).or.(itime.lt.(ireleasestart(kp)-lage(1)))) then |
---|
213 | go to 10 |
---|
214 | endif |
---|
215 | endif |
---|
216 | |
---|
217 | ! calculate date of release |
---|
218 | jul=bdate+real(ireleasestart(kp),kind=dp)/86400._dp ! this is the current day |
---|
219 | call caldate(jul,jjjjmmdd,ihmmss) |
---|
220 | write(areldate,'(i8.8)') jjjjmmdd |
---|
221 | write(areltime,'(i6.6)') ihmmss |
---|
222 | print*, areldate//areltime |
---|
223 | |
---|
224 | ! calculate date of field |
---|
225 | jul=bdate+real(itime,kind=dp)/86400._dp |
---|
226 | call caldate(jul,jjjjmmdd,ihmmss) |
---|
227 | write(adate,'(i8.8)') jjjjmmdd |
---|
228 | write(atime,'(i6.6)') ihmmss |
---|
229 | print*, adate//atime |
---|
230 | |
---|
231 | if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then |
---|
232 | if (ldirect.eq.1) then |
---|
233 | ! concentrations |
---|
234 | inquire(file=path(2)(1:length(2))//'grid_conc_nest_'//areldate// & |
---|
235 | areltime//'_'//anspec,exist=lexist) |
---|
236 | if(lexist.and..not.lnstartrel(kp)) then |
---|
237 | ! open and append to existing file |
---|
238 | open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_nest_'//areldate// & |
---|
239 | areltime//'_'//anspec,form='unformatted',status='old',action='write',access='append') |
---|
240 | else |
---|
241 | ! open new file |
---|
242 | open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_nest_'//areldate// & |
---|
243 | areltime//'_'//anspec,form='unformatted',status='replace',action='write') |
---|
244 | endif |
---|
245 | else |
---|
246 | ! residence times |
---|
247 | inquire(file=path(2)(1:length(2))//'grid_time_nest_'//areldate// & |
---|
248 | areltime//'_'//anspec,exist=lexist) |
---|
249 | if(lexist.and..not.lnstartrel(kp)) then |
---|
250 | ! open and append to existing file |
---|
251 | open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_nest_'//areldate// & |
---|
252 | areltime//'_'//anspec,form='unformatted',status='old',action='write',access='append') |
---|
253 | else |
---|
254 | ! open new file |
---|
255 | open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_nest_'//areldate// & |
---|
256 | areltime//'_'//anspec,form='unformatted',status='replace',action='write') |
---|
257 | endif |
---|
258 | endif |
---|
259 | write(unitoutgrid) jjjjmmdd |
---|
260 | write(unitoutgrid) ihmmss |
---|
261 | endif |
---|
262 | |
---|
263 | if ((iout.eq.2).or.(iout.eq.3)) then |
---|
264 | ! mixing ratio |
---|
265 | inquire(file=path(2)(1:length(2))//'grid_pptv_nest_'//areldate// & |
---|
266 | areltime//'_'//anspec,exist=lexist) |
---|
267 | if(lexist.and..not.lnstartrel(kp)) then |
---|
268 | ! open and append to existing file |
---|
269 | open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_nest_'//areldate// & |
---|
270 | areltime//'_'//anspec,form='unformatted',status='old',action='write',access='append') |
---|
271 | else |
---|
272 | ! open new file |
---|
273 | open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_nest_'//areldate// & |
---|
274 | areltime//'_'//anspec,form='unformatted',status='replace',action='write') |
---|
275 | endif |
---|
276 | write(unitoutgridppt) jjjjmmdd |
---|
277 | write(unitoutgridppt) ihmmss |
---|
278 | endif |
---|
279 | |
---|
280 | lnstartrel(kp)=.false. |
---|
281 | |
---|
282 | do nage=1,nageclass |
---|
283 | |
---|
284 | do jy=0,numygridn-1 |
---|
285 | do ix=0,numxgridn-1 |
---|
286 | |
---|
287 | ! ! WET DEPOSITION |
---|
288 | ! if ((WETDEP).and.(ldirect.gt.0)) then |
---|
289 | ! do l=1,nclassunc |
---|
290 | ! auxgrid(l)=wetgriduncn(ix,jy,ks,kp,l,nage) |
---|
291 | ! end do |
---|
292 | ! call mean(auxgrid,wetgrid(ix,jy), & |
---|
293 | ! wetgridsigma(ix,jy),nclassunc) |
---|
294 | ! ! Multiply by number of classes to get total concentration |
---|
295 | ! wetgrid(ix,jy)=wetgrid(ix,jy) & |
---|
296 | ! *nclassunc |
---|
297 | ! ! Calculate standard deviation of the mean |
---|
298 | ! wetgridsigma(ix,jy)= & |
---|
299 | ! wetgridsigma(ix,jy)* & |
---|
300 | ! sqrt(real(nclassunc)) |
---|
301 | ! endif |
---|
302 | |
---|
303 | ! ! DRY DEPOSITION |
---|
304 | ! if ((DRYDEP).and.(ldirect.gt.0)) then |
---|
305 | ! do l=1,nclassunc |
---|
306 | ! auxgrid(l)=drygriduncn(ix,jy,ks,kp,l,nage) |
---|
307 | ! end do |
---|
308 | ! call mean(auxgrid,drygrid(ix,jy), & |
---|
309 | ! drygridsigma(ix,jy),nclassunc) |
---|
310 | ! ! Multiply by number of classes to get total concentration |
---|
311 | ! drygrid(ix,jy)=drygrid(ix,jy)* & |
---|
312 | ! nclassunc |
---|
313 | ! ! Calculate standard deviation of the mean |
---|
314 | ! drygridsigma(ix,jy)= & |
---|
315 | ! drygridsigma(ix,jy)* & |
---|
316 | ! sqrt(real(nclassunc)) |
---|
317 | ! endif |
---|
318 | |
---|
319 | ! CONCENTRATION OR MIXING RATIO |
---|
320 | do kz=1,numzgrid |
---|
321 | do l=1,nclassunc |
---|
322 | auxgrid(l)=griduncn(ix,jy,kz,ks,kp,l,nage) |
---|
323 | end do |
---|
324 | call mean(auxgrid,grid(ix,jy,kz), & |
---|
325 | gridsigma(ix,jy,kz),nclassunc) |
---|
326 | ! Multiply by number of classes to get total concentration |
---|
327 | grid(ix,jy,kz)= & |
---|
328 | grid(ix,jy,kz)*nclassunc |
---|
329 | ! Calculate standard deviation of the mean |
---|
330 | gridsigma(ix,jy,kz)= & |
---|
331 | gridsigma(ix,jy,kz)* & |
---|
332 | sqrt(real(nclassunc)) |
---|
333 | end do |
---|
334 | end do |
---|
335 | end do |
---|
336 | |
---|
337 | |
---|
338 | !******************************************************************* |
---|
339 | ! Generate output: may be in concentration (ng/m3) or in mixing |
---|
340 | ! ratio (ppt) or both |
---|
341 | ! Output the position and the values alternated multiplied by |
---|
342 | ! 1 or -1, first line is number of values, number of positions |
---|
343 | ! For backward simulations, the unit is seconds, stored in grid_time |
---|
344 | !******************************************************************* |
---|
345 | |
---|
346 | ! Concentration output |
---|
347 | !********************* |
---|
348 | |
---|
349 | if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then |
---|
350 | |
---|
351 | ! ! Wet deposition |
---|
352 | ! sp_count_i=0 |
---|
353 | ! sp_count_r=0 |
---|
354 | ! sp_fact=-1. |
---|
355 | ! sp_zer=.true. |
---|
356 | ! if ((ldirect.eq.1).and.(WETDEP)) then |
---|
357 | ! do jy=0,numygridn-1 |
---|
358 | ! do ix=0,numxgridn-1 |
---|
359 | ! ! concentration greater zero |
---|
360 | ! if (wetgrid(ix,jy).gt.smallnum) then |
---|
361 | ! if (sp_zer.eqv..true.) then ! first non zero value |
---|
362 | ! sp_count_i=sp_count_i+1 |
---|
363 | ! sparse_dump_i(sp_count_i)=ix+jy*numxgridn |
---|
364 | ! sp_zer=.false. |
---|
365 | ! sp_fact=sp_fact*(-1.) |
---|
366 | ! endif |
---|
367 | ! sp_count_r=sp_count_r+1 |
---|
368 | ! sparse_dump_r(sp_count_r)= & |
---|
369 | ! sp_fact*1.e12*wetgrid(ix,jy)/arean(ix,jy) |
---|
370 | ! sparse_dump_u(sp_count_r)= & |
---|
371 | ! 1.e12*wetgridsigma(ix,jy)/area(ix,jy) |
---|
372 | ! else ! concentration is zero |
---|
373 | ! sp_zer=.true. |
---|
374 | ! endif |
---|
375 | ! end do |
---|
376 | ! end do |
---|
377 | ! else |
---|
378 | ! sp_count_i=0 |
---|
379 | ! sp_count_r=0 |
---|
380 | ! endif |
---|
381 | ! write(unitoutgrid) sp_count_i |
---|
382 | ! write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i) |
---|
383 | ! write(unitoutgrid) sp_count_r |
---|
384 | ! write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r) |
---|
385 | ! write(unitoutgrid) sp_count_r |
---|
386 | ! write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r) |
---|
387 | |
---|
388 | ! ! Dry deposition |
---|
389 | ! sp_count_i=0 |
---|
390 | ! sp_count_r=0 |
---|
391 | ! sp_fact=-1. |
---|
392 | ! sp_zer=.true. |
---|
393 | ! if ((ldirect.eq.1).and.(DRYDEP)) then |
---|
394 | ! do jy=0,numygridn-1 |
---|
395 | ! do ix=0,numxgridn-1 |
---|
396 | ! if (drygrid(ix,jy).gt.smallnum) then |
---|
397 | ! if (sp_zer.eqv..true.) then ! first non zero value |
---|
398 | ! sp_count_i=sp_count_i+1 |
---|
399 | ! sparse_dump_i(sp_count_i)=ix+jy*numxgridn |
---|
400 | ! sp_zer=.false. |
---|
401 | ! sp_fact=sp_fact*(-1.) |
---|
402 | ! endif |
---|
403 | ! sp_count_r=sp_count_r+1 |
---|
404 | ! sparse_dump_r(sp_count_r)= & |
---|
405 | ! sp_fact* & |
---|
406 | ! 1.e12*drygrid(ix,jy)/arean(ix,jy) |
---|
407 | ! sparse_dump_u(sp_count_r)= & |
---|
408 | ! 1.e12*drygridsigma(ix,jy)/area(ix,jy) |
---|
409 | ! else ! concentration is zero |
---|
410 | ! sp_zer=.true. |
---|
411 | ! endif |
---|
412 | ! end do |
---|
413 | ! end do |
---|
414 | ! else |
---|
415 | ! sp_count_i=0 |
---|
416 | ! sp_count_r=0 |
---|
417 | ! endif |
---|
418 | ! write(unitoutgrid) sp_count_i |
---|
419 | ! write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i) |
---|
420 | ! write(unitoutgrid) sp_count_r |
---|
421 | ! write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r) |
---|
422 | ! write(unitoutgrid) sp_count_r |
---|
423 | ! write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r) |
---|
424 | ! |
---|
425 | |
---|
426 | ! Concentrations |
---|
427 | |
---|
428 | ! surf_only write only 1st layer |
---|
429 | |
---|
430 | sp_count_i=0 |
---|
431 | sp_count_r=0 |
---|
432 | sp_fact=-1. |
---|
433 | sp_zer=.true. |
---|
434 | do kz=1,1 |
---|
435 | do jy=0,numygridn-1 |
---|
436 | do ix=0,numxgridn-1 |
---|
437 | if (grid(ix,jy,kz).gt.smallnum) then |
---|
438 | if (sp_zer.eqv..true.) then ! first non zero value |
---|
439 | sp_count_i=sp_count_i+1 |
---|
440 | sparse_dump_i(sp_count_i)= & |
---|
441 | ix+jy*numxgridn+kz*numxgridn*numygridn |
---|
442 | sp_zer=.false. |
---|
443 | sp_fact=sp_fact*(-1.) |
---|
444 | endif |
---|
445 | sp_count_r=sp_count_r+1 |
---|
446 | sparse_dump_r(sp_count_r)= & |
---|
447 | sp_fact* & |
---|
448 | grid(ix,jy,kz)* & |
---|
449 | factor3d(ix,jy,kz)/tot_mu(ks,kp) |
---|
450 | ! if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0) |
---|
451 | ! + write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp |
---|
452 | sparse_dump_u(sp_count_r)= & |
---|
453 | gridsigma(ix,jy,kz)* & |
---|
454 | factor3d(ix,jy,kz)/tot_mu(ks,kp) |
---|
455 | else ! concentration is zero |
---|
456 | sp_zer=.true. |
---|
457 | endif |
---|
458 | end do |
---|
459 | end do |
---|
460 | end do |
---|
461 | write(unitoutgrid) sp_count_i |
---|
462 | write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i) |
---|
463 | write(unitoutgrid) sp_count_r |
---|
464 | write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r) |
---|
465 | ! write(unitoutgrid) sp_count_r |
---|
466 | ! write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r) |
---|
467 | |
---|
468 | endif ! concentration output |
---|
469 | |
---|
470 | ! Mixing ratio output |
---|
471 | !******************** |
---|
472 | |
---|
473 | if ((iout.eq.2).or.(iout.eq.3)) then ! mixing ratio |
---|
474 | |
---|
475 | ! ! Wet deposition |
---|
476 | ! sp_count_i=0 |
---|
477 | ! sp_count_r=0 |
---|
478 | ! sp_fact=-1. |
---|
479 | ! sp_zer=.true. |
---|
480 | ! if ((ldirect.eq.1).and.(WETDEP)) then |
---|
481 | ! do jy=0,numygridn-1 |
---|
482 | ! do ix=0,numxgridn-1 |
---|
483 | ! if (wetgrid(ix,jy).gt.smallnum) then |
---|
484 | ! if (sp_zer.eqv..true.) then ! first non zero value |
---|
485 | ! sp_count_i=sp_count_i+1 |
---|
486 | ! sparse_dump_i(sp_count_i)= & |
---|
487 | ! ix+jy*numxgridn |
---|
488 | ! sp_zer=.false. |
---|
489 | ! sp_fact=sp_fact*(-1.) |
---|
490 | ! endif |
---|
491 | ! sp_count_r=sp_count_r+1 |
---|
492 | ! sparse_dump_r(sp_count_r)= & |
---|
493 | ! sp_fact* & |
---|
494 | ! 1.e12*wetgrid(ix,jy)/arean(ix,jy) |
---|
495 | ! sparse_dump_u(sp_count_r)= & |
---|
496 | ! 1.e12*wetgridsigma(ix,jy)/area(ix,jy) |
---|
497 | ! else ! concentration is zero |
---|
498 | ! sp_zer=.true. |
---|
499 | ! endif |
---|
500 | ! end do |
---|
501 | ! end do |
---|
502 | ! else |
---|
503 | ! sp_count_i=0 |
---|
504 | ! sp_count_r=0 |
---|
505 | ! endif |
---|
506 | ! write(unitoutgridppt) sp_count_i |
---|
507 | ! write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i) |
---|
508 | ! write(unitoutgridppt) sp_count_r |
---|
509 | ! write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r) |
---|
510 | ! write(unitoutgridppt) sp_count_r |
---|
511 | ! write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r) |
---|
512 | ! |
---|
513 | |
---|
514 | ! ! Dry deposition |
---|
515 | ! sp_count_i=0 |
---|
516 | ! sp_count_r=0 |
---|
517 | ! sp_fact=-1. |
---|
518 | ! sp_zer=.true. |
---|
519 | ! if ((ldirect.eq.1).and.(DRYDEP)) then |
---|
520 | ! do jy=0,numygridn-1 |
---|
521 | ! do ix=0,numxgridn-1 |
---|
522 | ! if (drygrid(ix,jy).gt.smallnum) then |
---|
523 | ! if (sp_zer.eqv..true.) then ! first non zero value |
---|
524 | ! sp_count_i=sp_count_i+1 |
---|
525 | ! sparse_dump_i(sp_count_i)= & |
---|
526 | ! ix+jy*numxgridn |
---|
527 | ! sp_zer=.false. |
---|
528 | ! sp_fact=sp_fact*(-1) |
---|
529 | ! endif |
---|
530 | ! sp_count_r=sp_count_r+1 |
---|
531 | ! sparse_dump_r(sp_count_r)= & |
---|
532 | ! sp_fact* & |
---|
533 | ! 1.e12*drygrid(ix,jy)/arean(ix,jy) |
---|
534 | ! sparse_dump_u(sp_count_r)= & |
---|
535 | ! 1.e12*drygridsigma(ix,jy)/area(ix,jy) |
---|
536 | ! else ! concentration is zero |
---|
537 | ! sp_zer=.true. |
---|
538 | ! endif |
---|
539 | ! end do |
---|
540 | ! end do |
---|
541 | ! else |
---|
542 | ! sp_count_i=0 |
---|
543 | ! sp_count_r=0 |
---|
544 | ! endif |
---|
545 | ! write(unitoutgridppt) sp_count_i |
---|
546 | ! write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i) |
---|
547 | ! write(unitoutgridppt) sp_count_r |
---|
548 | ! write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r) |
---|
549 | ! write(unitoutgridppt) sp_count_r |
---|
550 | ! write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r) |
---|
551 | ! |
---|
552 | |
---|
553 | ! Mixing ratios |
---|
554 | |
---|
555 | ! surf_only write only 1st layer |
---|
556 | |
---|
557 | sp_count_i=0 |
---|
558 | sp_count_r=0 |
---|
559 | sp_fact=-1. |
---|
560 | sp_zer=.true. |
---|
561 | do kz=1,1 |
---|
562 | do jy=0,numygridn-1 |
---|
563 | do ix=0,numxgridn-1 |
---|
564 | if (grid(ix,jy,kz).gt.smallnum) then |
---|
565 | if (sp_zer.eqv..true.) then ! first non zero value |
---|
566 | sp_count_i=sp_count_i+1 |
---|
567 | sparse_dump_i(sp_count_i)= & |
---|
568 | ix+jy*numxgridn+kz*numxgridn*numygridn |
---|
569 | sp_zer=.false. |
---|
570 | sp_fact=sp_fact*(-1.) |
---|
571 | endif |
---|
572 | sp_count_r=sp_count_r+1 |
---|
573 | sparse_dump_r(sp_count_r)= & |
---|
574 | sp_fact* & |
---|
575 | 1.e12*grid(ix,jy,kz) & |
---|
576 | /volumen(ix,jy,kz)/outnum* & |
---|
577 | weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz) |
---|
578 | sparse_dump_u(sp_count_r)= & |
---|
579 | 1.e12*gridsigma(ix,jy,kz)/volumen(ix,jy,kz)/ & |
---|
580 | outnum*weightair/weightmolar(ks)/ & |
---|
581 | densityoutgrid(ix,jy,kz) |
---|
582 | else ! concentration is zero |
---|
583 | sp_zer=.true. |
---|
584 | endif |
---|
585 | end do |
---|
586 | end do |
---|
587 | end do |
---|
588 | write(unitoutgridppt) sp_count_i |
---|
589 | write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i) |
---|
590 | write(unitoutgridppt) sp_count_r |
---|
591 | write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r) |
---|
592 | ! write(unitoutgridppt) sp_count_r |
---|
593 | ! write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r) |
---|
594 | |
---|
595 | endif ! output for ppt |
---|
596 | |
---|
597 | end do ! nageclass |
---|
598 | |
---|
599 | close(unitoutgridppt) |
---|
600 | close(unitoutgrid) |
---|
601 | |
---|
602 | ! itime is outside range |
---|
603 | 10 continue |
---|
604 | |
---|
605 | end do ! maxpointspec_act |
---|
606 | |
---|
607 | end do ! nspec |
---|
608 | |
---|
609 | |
---|
610 | ! RLT Aug 2017 |
---|
611 | ! Write out conversion factor for dry air |
---|
612 | inquire(file=path(2)(1:length(2))//'factor_drygrid_nest',exist=lexist) |
---|
613 | if (lexist.and..not.lnstart) then |
---|
614 | ! open and append |
---|
615 | open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',& |
---|
616 | status='old',action='write',access='append') |
---|
617 | else |
---|
618 | ! create new |
---|
619 | open(unitoutfactor,file=path(2)(1:length(2))//'factor_drygrid_nest',form='unformatted',& |
---|
620 | status='replace',action='write') |
---|
621 | endif |
---|
622 | sp_count_i=0 |
---|
623 | sp_count_r=0 |
---|
624 | sp_fact=-1. |
---|
625 | sp_zer=.true. |
---|
626 | do kz=1,1 |
---|
627 | do jy=0,numygridn-1 |
---|
628 | do ix=0,numxgridn-1 |
---|
629 | if (factor_drygrid(ix,jy,kz).gt.(1.+smallnum).or.factor_drygrid(ix,jy,kz).lt.(1.-smallnum)) then |
---|
630 | if (sp_zer.eqv..true.) then ! first value not equal to one |
---|
631 | sp_count_i=sp_count_i+1 |
---|
632 | sparse_dump_i(sp_count_i)= & |
---|
633 | ix+jy*numxgridn+kz*numxgridn*numygridn |
---|
634 | sp_zer=.false. |
---|
635 | sp_fact=sp_fact*(-1.) |
---|
636 | endif |
---|
637 | sp_count_r=sp_count_r+1 |
---|
638 | sparse_dump_r(sp_count_r)= & |
---|
639 | sp_fact*factor_drygrid(ix,jy,kz) |
---|
640 | else ! factor is one |
---|
641 | sp_zer=.true. |
---|
642 | endif |
---|
643 | end do |
---|
644 | end do |
---|
645 | end do |
---|
646 | write(unitoutfactor) sp_count_i |
---|
647 | write(unitoutfactor) (sparse_dump_i(i),i=1,sp_count_i) |
---|
648 | write(unitoutfactor) sp_count_r |
---|
649 | write(unitoutfactor) (sparse_dump_r(i),i=1,sp_count_r) |
---|
650 | close(unitoutfactor) |
---|
651 | |
---|
652 | ! reset lnstart |
---|
653 | if (lnstart) then |
---|
654 | lnstart=.false. |
---|
655 | endif |
---|
656 | |
---|
657 | ! Reinitialization of grid |
---|
658 | !************************* |
---|
659 | |
---|
660 | do ks=1,nspec |
---|
661 | do kp=1,maxpointspec_act |
---|
662 | do i=1,numreceptor |
---|
663 | creceptor(i,ks)=0. |
---|
664 | end do |
---|
665 | do jy=0,numygridn-1 |
---|
666 | do ix=0,numxgridn-1 |
---|
667 | do l=1,nclassunc |
---|
668 | do nage=1,nageclass |
---|
669 | do kz=1,numzgrid |
---|
670 | griduncn(ix,jy,kz,ks,kp,l,nage)=0. |
---|
671 | end do |
---|
672 | end do |
---|
673 | end do |
---|
674 | end do |
---|
675 | end do |
---|
676 | end do |
---|
677 | end do |
---|
678 | |
---|
679 | |
---|
680 | end subroutine concoutput_inversion_nest |
---|
681 | |
---|