source: flexpart.git/src/concoutput.f90 @ 478e9e6

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 478e9e6 was 478e9e6, checked in by Espen Sollum ATMOS <eso@…>, 9 years ago

dates file is now created at program start, instead of appending to it if it already exists

  • Property mode set to 100644
File size: 22.5 KB
Line 
1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine concoutput(itime,outnum,gridtotalunc,wetgridtotalunc, &
23       drygridtotalunc)
24  !                        i     i          o             o
25  !       o
26  !*****************************************************************************
27  !                                                                            *
28  !     Output of the concentration grid and the receptor concentrations.      *
29  !                                                                            *
30  !     Author: A. Stohl                                                       *
31  !                                                                            *
32  !     24 May 1995                                                            *
33  !                                                                            *
34  !     13 April 1999, Major update: if output size is smaller, dump output    *
35  !                    in sparse matrix format; additional output of           *
36  !                    uncertainty                                             *
37  !                                                                            *
38  !     05 April 2000, Major update: output of age classes; output for backward*
39  !                    runs is time spent in grid cell times total mass of     *
40  !                    species.                                                *
41  !                                                                            *
42  !     17 February 2002, Appropriate dimensions for backward and forward runs *
43  !                       are now specified in file par_mod                    *
44  !                                                                            *
45  !     June 2006, write grid in sparse matrix with a single write command     *
46  !                in order to save disk space                                 *
47  !                                                                            *
48  !     2008 new sparse matrix format                                          *
49  !                                                                            *
50  !*****************************************************************************
51  !                                                                            *
52  ! Variables:                                                                 *
53  ! outnum          number of samples                                          *
54  ! ncells          number of cells with non-zero concentrations               *
55  ! sparse          .true. if in sparse matrix format, else .false.            *
56  ! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
57  !                                                                            *
58  !*****************************************************************************
59
60  use unc_mod
61  use point_mod
62  use outg_mod
63  use par_mod
64  use com_mod
65
66  implicit none
67
68  real(kind=dp) :: jul
69  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
70  integer :: sp_count_i,sp_count_r
71  real :: sp_fact
72  real :: outnum,densityoutrecept(maxreceptor),xl,yl
73
74  !real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
75  !    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
76  !    +    maxageclass)
77  !real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
78  !    +       maxageclass)
79  !real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
80  !    +       maxpointspec_act,maxageclass)
81  !real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
82  !    +       maxpointspec_act,maxageclass),
83  !    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
84  !    +     maxpointspec_act,maxageclass),
85  !    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
86  !    +     maxpointspec_act,maxageclass)
87  !real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
88  !real sparse_dump_r(numxgrid*numygrid*numzgrid)
89  !integer sparse_dump_i(numxgrid*numygrid*numzgrid)
90
91  !real sparse_dump_u(numxgrid*numygrid*numzgrid)
92  real :: auxgrid(nclassunc),gridtotal,gridsigmatotal,gridtotalunc
93  real :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
94  real :: drygridtotal,drygridsigmatotal,drygridtotalunc
95  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
96  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
97  real,parameter :: weightair=28.97
98  logical :: sp_zer
99  LOGICAL,save :: init=.true.
100  character :: adate*8,atime*6
101  character(len=3) :: anspec
102  integer :: mind
103! mind        eso:added to ensure identical results between 2&3-fields versions
104  CHARACTER(LEN=8),save :: file_stat='REPLACE'
105
106  ! Determine current calendar date, needed for the file name
107  !**********************************************************
108
109  jul=bdate+real(itime,kind=dp)/86400._dp
110  call caldate(jul,jjjjmmdd,ihmmss)
111  write(adate,'(i8.8)') jjjjmmdd
112  write(atime,'(i6.6)') ihmmss
113  OPEN(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND', STATUS=file_stat)
114  write(unitdates,'(a)') adate//atime
115  close(unitdates) 
116
117  ! Overwrite existing dates file on first call, later append to it
118  ! This fixes a bug where the dates file kept growing across multiple runs
119  ! TODO check if the 'always APPEND'-behaviour is useful in other scenarioes
120  ! e.g. (restart?)
121  IF (init) THEN
122    file_stat='OLD'
123    init=.false.
124  END IF
125
126
127
128  ! For forward simulations, output fields have dimension MAXSPEC,
129  ! for backward simulations, output fields have dimension MAXPOINT.
130  ! Thus, make loops either about nspec, or about numpoint
131  !*****************************************************************
132
133
134  if (ldirect.eq.1) then
135    do ks=1,nspec
136      do kp=1,maxpointspec_act
137        tot_mu(ks,kp)=1
138      end do
139    end do
140  else
141    do ks=1,nspec
142      do kp=1,maxpointspec_act
143        tot_mu(ks,kp)=xmass(kp,ks)
144      end do
145    end do
146  endif
147
148
149  !*******************************************************************
150  ! Compute air density: sufficiently accurate to take it
151  ! from coarse grid at some time
152  ! Determine center altitude of output layer, and interpolate density
153  ! data to that altitude
154  !*******************************************************************
155
156  mind=memind(2)
157  do kz=1,numzgrid
158    if (kz.eq.1) then
159      halfheight=outheight(1)/2.
160    else
161      halfheight=(outheight(kz)+outheight(kz-1))/2.
162    endif
163    do kzz=2,nz
164      if ((height(kzz-1).lt.halfheight).and. &
165           (height(kzz).gt.halfheight)) goto 46
166    end do
16746   kzz=max(min(kzz,nz),2)
168    dz1=halfheight-height(kzz-1)
169    dz2=height(kzz)-halfheight
170    dz=dz1+dz2
171    do jy=0,numygrid-1
172      do ix=0,numxgrid-1
173        xl=outlon0+real(ix)*dxout
174        yl=outlat0+real(jy)*dyout
175        xl=(xl-xlon0)/dx
176        yl=(yl-ylat0)/dy !v9.1.1
177        iix=max(min(nint(xl),nxmin1),0)
178        jjy=max(min(nint(yl),nymin1),0)
179        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
180        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
181        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
182             rho(iix,jjy,kzz-1,mind)*dz2)/dz
183      end do
184    end do
185  end do
186
187  do i=1,numreceptor
188    xl=xreceptor(i)
189    yl=yreceptor(i)
190    iix=max(min(nint(xl),nxmin1),0)
191    jjy=max(min(nint(yl),nymin1),0)
192    !densityoutrecept(i)=rho(iix,jjy,1,2)
193    densityoutrecept(i)=rho(iix,jjy,1,mind)
194  end do
195
196
197  ! Output is different for forward and backward simulations
198    do kz=1,numzgrid
199      do jy=0,numygrid-1
200        do ix=0,numxgrid-1
201          if (ldirect.eq.1) then
202            factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
203          else
204            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
205          endif
206        end do
207      end do
208    end do
209
210  !*********************************************************************
211  ! Determine the standard deviation of the mean concentration or mixing
212  ! ratio (uncertainty of the output) and the dry and wet deposition
213  !*********************************************************************
214
215  gridtotal=0.
216  gridsigmatotal=0.
217  gridtotalunc=0.
218  wetgridtotal=0.
219  wetgridsigmatotal=0.
220  wetgridtotalunc=0.
221  drygridtotal=0.
222  drygridsigmatotal=0.
223  drygridtotalunc=0.
224
225  do ks=1,nspec
226
227  write(anspec,'(i3.3)') ks
228  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
229    if (ldirect.eq.1) then
230      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
231           atime//'_'//anspec,form='unformatted')
232    else
233      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
234           atime//'_'//anspec,form='unformatted')
235    endif
236     write(unitoutgrid) itime
237   endif
238
239  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
240   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
241        atime//'_'//anspec,form='unformatted')
242
243    write(unitoutgridppt) itime
244  endif
245
246  do kp=1,maxpointspec_act
247  do nage=1,nageclass
248
249    do jy=0,numygrid-1
250      do ix=0,numxgrid-1
251
252  ! WET DEPOSITION
253        if ((WETDEP).and.(ldirect.gt.0)) then
254            do l=1,nclassunc
255              auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
256            end do
257            call mean(auxgrid,wetgrid(ix,jy), &
258                 wetgridsigma(ix,jy),nclassunc)
259  ! Multiply by number of classes to get total concentration
260            wetgrid(ix,jy)=wetgrid(ix,jy) &
261                 *nclassunc
262            wetgridtotal=wetgridtotal+wetgrid(ix,jy)
263  ! Calculate standard deviation of the mean
264            wetgridsigma(ix,jy)= &
265                 wetgridsigma(ix,jy)* &
266                 sqrt(real(nclassunc))
267            wetgridsigmatotal=wetgridsigmatotal+ &
268                 wetgridsigma(ix,jy)
269        endif
270
271  ! DRY DEPOSITION
272        if ((DRYDEP).and.(ldirect.gt.0)) then
273            do l=1,nclassunc
274              auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
275            end do
276            call mean(auxgrid,drygrid(ix,jy), &
277                 drygridsigma(ix,jy),nclassunc)
278  ! Multiply by number of classes to get total concentration
279            drygrid(ix,jy)=drygrid(ix,jy)* &
280                 nclassunc
281            drygridtotal=drygridtotal+drygrid(ix,jy)
282  ! Calculate standard deviation of the mean
283            drygridsigma(ix,jy)= &
284                 drygridsigma(ix,jy)* &
285                 sqrt(real(nclassunc))
286125         drygridsigmatotal=drygridsigmatotal+ &
287                 drygridsigma(ix,jy)
288        endif
289
290  ! CONCENTRATION OR MIXING RATIO
291        do kz=1,numzgrid
292            do l=1,nclassunc
293              auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
294            end do
295            call mean(auxgrid,grid(ix,jy,kz), &
296                 gridsigma(ix,jy,kz),nclassunc)
297  ! Multiply by number of classes to get total concentration
298            grid(ix,jy,kz)= &
299                 grid(ix,jy,kz)*nclassunc
300            gridtotal=gridtotal+grid(ix,jy,kz)
301  ! Calculate standard deviation of the mean
302            gridsigma(ix,jy,kz)= &
303                 gridsigma(ix,jy,kz)* &
304                 sqrt(real(nclassunc))
305            gridsigmatotal=gridsigmatotal+ &
306                 gridsigma(ix,jy,kz)
307        end do
308      end do
309    end do
310
311
312
313
314  !*******************************************************************
315  ! Generate output: may be in concentration (ng/m3) or in mixing
316  ! ratio (ppt) or both
317  ! Output the position and the values alternated multiplied by
318  ! 1 or -1, first line is number of values, number of positions
319  ! For backward simulations, the unit is seconds, stored in grid_time
320  !*******************************************************************
321
322  ! Concentration output
323  !*********************
324  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
325
326  ! Wet deposition
327         sp_count_i=0
328         sp_count_r=0
329         sp_fact=-1.
330         sp_zer=.true.
331         if ((ldirect.eq.1).and.(WETDEP)) then
332         do jy=0,numygrid-1
333            do ix=0,numxgrid-1
334  !oncentraion greater zero
335              if (wetgrid(ix,jy).gt.smallnum) then
336                 if (sp_zer.eqv..true.) then ! first non zero value
337                    sp_count_i=sp_count_i+1
338                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
339                    sp_zer=.false.
340                    sp_fact=sp_fact*(-1.)
341                 endif
342                 sp_count_r=sp_count_r+1
343                 sparse_dump_r(sp_count_r)= &
344                      sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
345  !                sparse_dump_u(sp_count_r)=
346  !+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
347              else ! concentration is zero
348                  sp_zer=.true.
349              endif
350            end do
351         end do
352         else
353            sp_count_i=0
354            sp_count_r=0
355         endif
356         write(unitoutgrid) sp_count_i
357         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
358         write(unitoutgrid) sp_count_r
359         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
360  !       write(unitoutgrid) sp_count_u
361  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
362
363  ! Dry deposition
364         sp_count_i=0
365         sp_count_r=0
366         sp_fact=-1.
367         sp_zer=.true.
368         if ((ldirect.eq.1).and.(DRYDEP)) then
369          do jy=0,numygrid-1
370            do ix=0,numxgrid-1
371              if (drygrid(ix,jy).gt.smallnum) then
372                 if (sp_zer.eqv..true.) then ! first non zero value
373                    sp_count_i=sp_count_i+1
374                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
375                    sp_zer=.false.
376                    sp_fact=sp_fact*(-1.)
377                 endif
378                 sp_count_r=sp_count_r+1
379                 sparse_dump_r(sp_count_r)= &
380                      sp_fact* &
381                      1.e12*drygrid(ix,jy)/area(ix,jy)
382  !                sparse_dump_u(sp_count_r)=
383  !+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
384              else ! concentration is zero
385                  sp_zer=.true.
386              endif
387            end do
388          end do
389         else
390            sp_count_i=0
391            sp_count_r=0
392         endif
393         write(unitoutgrid) sp_count_i
394         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
395         write(unitoutgrid) sp_count_r
396         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
397  !       write(*,*) sp_count_u
398  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
399
400
401
402  ! Concentrations
403         sp_count_i=0
404         sp_count_r=0
405         sp_fact=-1.
406         sp_zer=.true.
407          do kz=1,numzgrid
408            do jy=0,numygrid-1
409              do ix=0,numxgrid-1
410                if (grid(ix,jy,kz).gt.smallnum) then
411                  if (sp_zer.eqv..true.) then ! first non zero value
412                    sp_count_i=sp_count_i+1
413                    sparse_dump_i(sp_count_i)= &
414                         ix+jy*numxgrid+kz*numxgrid*numygrid
415                    sp_zer=.false.
416                    sp_fact=sp_fact*(-1.)
417                   endif
418                   sp_count_r=sp_count_r+1
419                   sparse_dump_r(sp_count_r)= &
420                        sp_fact* &
421                        grid(ix,jy,kz)* &
422                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
423  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
424  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
425  !                sparse_dump_u(sp_count_r)=
426  !+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
427  !+               factor(ix,jy,kz)/tot_mu(ks,kp)
428              else ! concentration is zero
429                  sp_zer=.true.
430              endif
431              end do
432            end do
433          end do
434         write(unitoutgrid) sp_count_i
435         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
436         write(unitoutgrid) sp_count_r
437         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
438  !       write(unitoutgrid) sp_count_u
439  !       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
440
441
442
443    endif !  concentration output
444
445  ! Mixing ratio output
446  !********************
447
448  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
449
450  ! Wet deposition
451         sp_count_i=0
452         sp_count_r=0
453         sp_fact=-1.
454         sp_zer=.true.
455         if ((ldirect.eq.1).and.(WETDEP)) then
456          do jy=0,numygrid-1
457            do ix=0,numxgrid-1
458                if (wetgrid(ix,jy).gt.smallnum) then
459                  if (sp_zer.eqv..true.) then ! first non zero value
460                    sp_count_i=sp_count_i+1
461                    sparse_dump_i(sp_count_i)= &
462                         ix+jy*numxgrid
463                    sp_zer=.false.
464                    sp_fact=sp_fact*(-1.)
465                 endif
466                 sp_count_r=sp_count_r+1
467                 sparse_dump_r(sp_count_r)= &
468                      sp_fact* &
469                      1.e12*wetgrid(ix,jy)/area(ix,jy)
470  !                sparse_dump_u(sp_count_r)=
471  !    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
472              else ! concentration is zero
473                  sp_zer=.true.
474              endif
475            end do
476          end do
477         else
478           sp_count_i=0
479           sp_count_r=0
480         endif
481         write(unitoutgridppt) sp_count_i
482         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
483         write(unitoutgridppt) sp_count_r
484         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
485  !       write(unitoutgridppt) sp_count_u
486  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
487
488
489  ! Dry deposition
490         sp_count_i=0
491         sp_count_r=0
492         sp_fact=-1.
493         sp_zer=.true.
494         if ((ldirect.eq.1).and.(DRYDEP)) then
495          do jy=0,numygrid-1
496            do ix=0,numxgrid-1
497                if (drygrid(ix,jy).gt.smallnum) then
498                  if (sp_zer.eqv..true.) then ! first non zero value
499                    sp_count_i=sp_count_i+1
500                    sparse_dump_i(sp_count_i)= &
501                         ix+jy*numxgrid
502                    sp_zer=.false.
503                    sp_fact=sp_fact*(-1)
504                 endif
505                 sp_count_r=sp_count_r+1
506                 sparse_dump_r(sp_count_r)= &
507                      sp_fact* &
508                      1.e12*drygrid(ix,jy)/area(ix,jy)
509  !                sparse_dump_u(sp_count_r)=
510  !    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
511              else ! concentration is zero
512                  sp_zer=.true.
513              endif
514            end do
515          end do
516         else
517           sp_count_i=0
518           sp_count_r=0
519         endif
520         write(unitoutgridppt) sp_count_i
521         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
522         write(unitoutgridppt) sp_count_r
523         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
524  !       write(unitoutgridppt) sp_count_u
525  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
526
527
528  ! Mixing ratios
529         sp_count_i=0
530         sp_count_r=0
531         sp_fact=-1.
532         sp_zer=.true.
533          do kz=1,numzgrid
534            do jy=0,numygrid-1
535              do ix=0,numxgrid-1
536                if (grid(ix,jy,kz).gt.smallnum) then
537                  if (sp_zer.eqv..true.) then ! first non zero value
538                    sp_count_i=sp_count_i+1
539                    sparse_dump_i(sp_count_i)= &
540                         ix+jy*numxgrid+kz*numxgrid*numygrid
541                    sp_zer=.false.
542                    sp_fact=sp_fact*(-1.)
543                 endif
544                 sp_count_r=sp_count_r+1
545                 sparse_dump_r(sp_count_r)= &
546                      sp_fact* &
547                      1.e12*grid(ix,jy,kz) &
548                      /volume(ix,jy,kz)/outnum* &
549                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
550  !                sparse_dump_u(sp_count_r)=
551  !+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
552  !+              outnum*weightair/weightmolar(ks)/
553  !+              densityoutgrid(ix,jy,kz)
554              else ! concentration is zero
555                  sp_zer=.true.
556              endif
557              end do
558            end do
559          end do
560         write(unitoutgridppt) sp_count_i
561         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
562         write(unitoutgridppt) sp_count_r
563         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
564  !       write(unitoutgridppt) sp_count_u
565  !       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
566
567      endif ! output for ppt
568
569  end do
570  end do
571
572    close(unitoutgridppt)
573    close(unitoutgrid)
574
575  end do
576
577  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
578  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
579       wetgridtotal
580  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
581       drygridtotal
582
583  ! Dump of receptor concentrations
584
585    if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
586      write(unitoutreceptppt) itime
587      do ks=1,nspec
588        write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
589             weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
590      end do
591    endif
592
593  ! Dump of receptor concentrations
594
595    if (numreceptor.gt.0) then
596      write(unitoutrecept) itime
597      do ks=1,nspec
598        write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
599             i=1,numreceptor)
600      end do
601    endif
602
603
604
605  ! Reinitialization of grid
606  !*************************
607
608  do ks=1,nspec
609  do kp=1,maxpointspec_act
610    do i=1,numreceptor
611      creceptor(i,ks)=0.
612    end do
613    do jy=0,numygrid-1
614      do ix=0,numxgrid-1
615        do l=1,nclassunc
616          do nage=1,nageclass
617            do kz=1,numzgrid
618              gridunc(ix,jy,kz,ks,kp,l,nage)=0.
619            end do
620          end do
621        end do
622      end do
623    end do
624  end do
625  end do
626
627
628end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG