source: flexpart.git/src/concoutput.f90 @ b069789

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

If the 'dates' file exists on simulation start, backup to 'dates.bak' and create new file 'dates' instead of appending to it

  • Property mode set to 100644
File size: 23.3 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  logical :: ldates_file
106  integer :: ierr
107  character(LEN=100) :: dates_char
108
109! Determine current calendar date, needed for the file name
110!**********************************************************
111
112  jul=bdate+real(itime,kind=dp)/86400._dp
113  call caldate(jul,jjjjmmdd,ihmmss)
114  write(adate,'(i8.8)') jjjjmmdd
115  write(atime,'(i6.6)') ihmmss
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
120! If 'dates' file exists, make a backup
121  inquire(file=path(2)(1:length(2))//'dates', exist=ldates_file)
122  if (ldates_file.and.init) then
123    open(unit=unitdates, file=path(2)(1:length(2))//'dates',form='formatted', &
124         &access='sequential', status='old', action='read', iostat=ierr)
125    open(unit=unittmp, file=path(2)(1:length(2))//'dates.bak', access='sequential', &
126         &status='replace', action='write', form='formatted', iostat=ierr)
127    do while (.true.)
128      read(unitdates, '(a)', iostat=ierr) dates_char
129      if (ierr.ne.0) exit
130      write(unit=unittmp, fmt='(a)', iostat=ierr, advance='yes') trim(dates_char)
131!      if (ierr.ne.0) write(*,*) "Write error, ", ierr
132    end do
133    close(unit=unitdates)
134    close(unit=unittmp)
135  end if
136
137  open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND', STATUS=file_stat)
138  write(unitdates,'(a)') adate//atime
139  close(unitdates) 
140
141  IF (init) THEN
142    file_stat='OLD'
143    init=.false.
144  END IF
145
146
147! For forward simulations, output fields have dimension MAXSPEC,
148! for backward simulations, output fields have dimension MAXPOINT.
149! Thus, make loops either about nspec, or about numpoint
150!*****************************************************************
151
152
153  if (ldirect.eq.1) then
154    do ks=1,nspec
155      do kp=1,maxpointspec_act
156        tot_mu(ks,kp)=1
157      end do
158    end do
159  else
160    do ks=1,nspec
161      do kp=1,maxpointspec_act
162        tot_mu(ks,kp)=xmass(kp,ks)
163      end do
164    end do
165  endif
166
167
168!*******************************************************************
169! Compute air density: sufficiently accurate to take it
170! from coarse grid at some time
171! Determine center altitude of output layer, and interpolate density
172! data to that altitude
173!*******************************************************************
174
175  mind=memind(2)
176  do kz=1,numzgrid
177    if (kz.eq.1) then
178      halfheight=outheight(1)/2.
179    else
180      halfheight=(outheight(kz)+outheight(kz-1))/2.
181    endif
182    do kzz=2,nz
183      if ((height(kzz-1).lt.halfheight).and. &
184           (height(kzz).gt.halfheight)) goto 46
185    end do
18646  kzz=max(min(kzz,nz),2)
187    dz1=halfheight-height(kzz-1)
188    dz2=height(kzz)-halfheight
189    dz=dz1+dz2
190    do jy=0,numygrid-1
191      do ix=0,numxgrid-1
192        xl=outlon0+real(ix)*dxout
193        yl=outlat0+real(jy)*dyout
194        xl=(xl-xlon0)/dx
195        yl=(yl-ylat0)/dy !v9.1.1
196        iix=max(min(nint(xl),nxmin1),0)
197        jjy=max(min(nint(yl),nymin1),0)
198! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
199!      rho(iix,jjy,kzz-1,2)*dz2)/dz
200        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
201             rho(iix,jjy,kzz-1,mind)*dz2)/dz
202      end do
203    end do
204  end do
205
206  do i=1,numreceptor
207    xl=xreceptor(i)
208    yl=yreceptor(i)
209    iix=max(min(nint(xl),nxmin1),0)
210    jjy=max(min(nint(yl),nymin1),0)
211!densityoutrecept(i)=rho(iix,jjy,1,2)
212    densityoutrecept(i)=rho(iix,jjy,1,mind)
213  end do
214
215
216! Output is different for forward and backward simulations
217  do kz=1,numzgrid
218    do jy=0,numygrid-1
219      do ix=0,numxgrid-1
220        if (ldirect.eq.1) then
221          factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
222        else
223          factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
224        endif
225      end do
226    end do
227  end do
228
229!*********************************************************************
230! Determine the standard deviation of the mean concentration or mixing
231! ratio (uncertainty of the output) and the dry and wet deposition
232!*********************************************************************
233
234  gridtotal=0.
235  gridsigmatotal=0.
236  gridtotalunc=0.
237  wetgridtotal=0.
238  wetgridsigmatotal=0.
239  wetgridtotalunc=0.
240  drygridtotal=0.
241  drygridsigmatotal=0.
242  drygridtotalunc=0.
243
244  do ks=1,nspec
245
246    write(anspec,'(i3.3)') ks
247    if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
248      if (ldirect.eq.1) then
249        open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
250             atime//'_'//anspec,form='unformatted')
251      else
252        open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
253             atime//'_'//anspec,form='unformatted')
254      endif
255      write(unitoutgrid) itime
256    endif
257
258    if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
259      open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
260           atime//'_'//anspec,form='unformatted')
261
262      write(unitoutgridppt) itime
263    endif
264
265    do kp=1,maxpointspec_act
266      do nage=1,nageclass
267
268        do jy=0,numygrid-1
269          do ix=0,numxgrid-1
270
271! WET DEPOSITION
272            if ((WETDEP).and.(ldirect.gt.0)) then
273              do l=1,nclassunc
274                auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
275              end do
276              call mean(auxgrid,wetgrid(ix,jy), &
277                   wetgridsigma(ix,jy),nclassunc)
278! Multiply by number of classes to get total concentration
279              wetgrid(ix,jy)=wetgrid(ix,jy) &
280                   *nclassunc
281              wetgridtotal=wetgridtotal+wetgrid(ix,jy)
282! Calculate standard deviation of the mean
283              wetgridsigma(ix,jy)= &
284                   wetgridsigma(ix,jy)* &
285                   sqrt(real(nclassunc))
286              wetgridsigmatotal=wetgridsigmatotal+ &
287                   wetgridsigma(ix,jy)
288            endif
289
290! DRY DEPOSITION
291            if ((DRYDEP).and.(ldirect.gt.0)) then
292              do l=1,nclassunc
293                auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
294              end do
295              call mean(auxgrid,drygrid(ix,jy), &
296                   drygridsigma(ix,jy),nclassunc)
297! Multiply by number of classes to get total concentration
298              drygrid(ix,jy)=drygrid(ix,jy)* &
299                   nclassunc
300              drygridtotal=drygridtotal+drygrid(ix,jy)
301! Calculate standard deviation of the mean
302              drygridsigma(ix,jy)= &
303                   drygridsigma(ix,jy)* &
304                   sqrt(real(nclassunc))
305125           drygridsigmatotal=drygridsigmatotal+ &
306                   drygridsigma(ix,jy)
307            endif
308
309! CONCENTRATION OR MIXING RATIO
310            do kz=1,numzgrid
311              do l=1,nclassunc
312                auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
313              end do
314              call mean(auxgrid,grid(ix,jy,kz), &
315                   gridsigma(ix,jy,kz),nclassunc)
316! Multiply by number of classes to get total concentration
317              grid(ix,jy,kz)= &
318                   grid(ix,jy,kz)*nclassunc
319              gridtotal=gridtotal+grid(ix,jy,kz)
320! Calculate standard deviation of the mean
321              gridsigma(ix,jy,kz)= &
322                   gridsigma(ix,jy,kz)* &
323                   sqrt(real(nclassunc))
324              gridsigmatotal=gridsigmatotal+ &
325                   gridsigma(ix,jy,kz)
326            end do
327          end do
328        end do
329
330
331
332
333!*******************************************************************
334! Generate output: may be in concentration (ng/m3) or in mixing
335! ratio (ppt) or both
336! Output the position and the values alternated multiplied by
337! 1 or -1, first line is number of values, number of positions
338! For backward simulations, the unit is seconds, stored in grid_time
339!*******************************************************************
340
341! Concentration output
342!*********************
343        if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
344
345! Wet deposition
346          sp_count_i=0
347          sp_count_r=0
348          sp_fact=-1.
349          sp_zer=.true.
350          if ((ldirect.eq.1).and.(WETDEP)) then
351            do jy=0,numygrid-1
352              do ix=0,numxgrid-1
353!oncentraion greater zero
354                if (wetgrid(ix,jy).gt.smallnum) then
355                  if (sp_zer.eqv..true.) then ! first non zero value
356                    sp_count_i=sp_count_i+1
357                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
358                    sp_zer=.false.
359                    sp_fact=sp_fact*(-1.)
360                  endif
361                  sp_count_r=sp_count_r+1
362                  sparse_dump_r(sp_count_r)= &
363                       sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
364!                sparse_dump_u(sp_count_r)=
365!+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
366                else ! concentration is zero
367                  sp_zer=.true.
368                endif
369              end do
370            end do
371          else
372            sp_count_i=0
373            sp_count_r=0
374          endif
375          write(unitoutgrid) sp_count_i
376          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
377          write(unitoutgrid) sp_count_r
378          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
379!       write(unitoutgrid) sp_count_u
380!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
381
382! Dry deposition
383          sp_count_i=0
384          sp_count_r=0
385          sp_fact=-1.
386          sp_zer=.true.
387          if ((ldirect.eq.1).and.(DRYDEP)) then
388            do jy=0,numygrid-1
389              do ix=0,numxgrid-1
390                if (drygrid(ix,jy).gt.smallnum) then
391                  if (sp_zer.eqv..true.) then ! first non zero value
392                    sp_count_i=sp_count_i+1
393                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
394                    sp_zer=.false.
395                    sp_fact=sp_fact*(-1.)
396                  endif
397                  sp_count_r=sp_count_r+1
398                  sparse_dump_r(sp_count_r)= &
399                       sp_fact* &
400                       1.e12*drygrid(ix,jy)/area(ix,jy)
401!                sparse_dump_u(sp_count_r)=
402!+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
403                else ! concentration is zero
404                  sp_zer=.true.
405                endif
406              end do
407            end do
408          else
409            sp_count_i=0
410            sp_count_r=0
411          endif
412          write(unitoutgrid) sp_count_i
413          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
414          write(unitoutgrid) sp_count_r
415          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
416!       write(*,*) sp_count_u
417!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
418
419
420
421! Concentrations
422          sp_count_i=0
423          sp_count_r=0
424          sp_fact=-1.
425          sp_zer=.true.
426          do kz=1,numzgrid
427            do jy=0,numygrid-1
428              do ix=0,numxgrid-1
429                if (grid(ix,jy,kz).gt.smallnum) then
430                  if (sp_zer.eqv..true.) then ! first non zero value
431                    sp_count_i=sp_count_i+1
432                    sparse_dump_i(sp_count_i)= &
433                         ix+jy*numxgrid+kz*numxgrid*numygrid
434                    sp_zer=.false.
435                    sp_fact=sp_fact*(-1.)
436                  endif
437                  sp_count_r=sp_count_r+1
438                  sparse_dump_r(sp_count_r)= &
439                       sp_fact* &
440                       grid(ix,jy,kz)* &
441                       factor3d(ix,jy,kz)/tot_mu(ks,kp)
442!                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
443!    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
444!                sparse_dump_u(sp_count_r)=
445!+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
446!+               factor(ix,jy,kz)/tot_mu(ks,kp)
447                else ! concentration is zero
448                  sp_zer=.true.
449                endif
450              end do
451            end do
452          end do
453          write(unitoutgrid) sp_count_i
454          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
455          write(unitoutgrid) sp_count_r
456          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
457!       write(unitoutgrid) sp_count_u
458!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
459
460
461
462        endif !  concentration output
463
464! Mixing ratio output
465!********************
466
467        if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
468
469! Wet deposition
470          sp_count_i=0
471          sp_count_r=0
472          sp_fact=-1.
473          sp_zer=.true.
474          if ((ldirect.eq.1).and.(WETDEP)) then
475            do jy=0,numygrid-1
476              do ix=0,numxgrid-1
477                if (wetgrid(ix,jy).gt.smallnum) then
478                  if (sp_zer.eqv..true.) then ! first non zero value
479                    sp_count_i=sp_count_i+1
480                    sparse_dump_i(sp_count_i)= &
481                         ix+jy*numxgrid
482                    sp_zer=.false.
483                    sp_fact=sp_fact*(-1.)
484                  endif
485                  sp_count_r=sp_count_r+1
486                  sparse_dump_r(sp_count_r)= &
487                       sp_fact* &
488                       1.e12*wetgrid(ix,jy)/area(ix,jy)
489!                sparse_dump_u(sp_count_r)=
490!    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
491                else ! concentration is zero
492                  sp_zer=.true.
493                endif
494              end do
495            end do
496          else
497            sp_count_i=0
498            sp_count_r=0
499          endif
500          write(unitoutgridppt) sp_count_i
501          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
502          write(unitoutgridppt) sp_count_r
503          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
504!       write(unitoutgridppt) sp_count_u
505!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
506
507
508! Dry deposition
509          sp_count_i=0
510          sp_count_r=0
511          sp_fact=-1.
512          sp_zer=.true.
513          if ((ldirect.eq.1).and.(DRYDEP)) then
514            do jy=0,numygrid-1
515              do ix=0,numxgrid-1
516                if (drygrid(ix,jy).gt.smallnum) then
517                  if (sp_zer.eqv..true.) then ! first non zero value
518                    sp_count_i=sp_count_i+1
519                    sparse_dump_i(sp_count_i)= &
520                         ix+jy*numxgrid
521                    sp_zer=.false.
522                    sp_fact=sp_fact*(-1)
523                  endif
524                  sp_count_r=sp_count_r+1
525                  sparse_dump_r(sp_count_r)= &
526                       sp_fact* &
527                       1.e12*drygrid(ix,jy)/area(ix,jy)
528!                sparse_dump_u(sp_count_r)=
529!    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
530                else ! concentration is zero
531                  sp_zer=.true.
532                endif
533              end do
534            end do
535          else
536            sp_count_i=0
537            sp_count_r=0
538          endif
539          write(unitoutgridppt) sp_count_i
540          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
541          write(unitoutgridppt) sp_count_r
542          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
543!       write(unitoutgridppt) sp_count_u
544!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
545
546
547! Mixing ratios
548          sp_count_i=0
549          sp_count_r=0
550          sp_fact=-1.
551          sp_zer=.true.
552          do kz=1,numzgrid
553            do jy=0,numygrid-1
554              do ix=0,numxgrid-1
555                if (grid(ix,jy,kz).gt.smallnum) then
556                  if (sp_zer.eqv..true.) then ! first non zero value
557                    sp_count_i=sp_count_i+1
558                    sparse_dump_i(sp_count_i)= &
559                         ix+jy*numxgrid+kz*numxgrid*numygrid
560                    sp_zer=.false.
561                    sp_fact=sp_fact*(-1.)
562                  endif
563                  sp_count_r=sp_count_r+1
564                  sparse_dump_r(sp_count_r)= &
565                       sp_fact* &
566                       1.e12*grid(ix,jy,kz) &
567                       /volume(ix,jy,kz)/outnum* &
568                       weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
569!                sparse_dump_u(sp_count_r)=
570!+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
571!+              outnum*weightair/weightmolar(ks)/
572!+              densityoutgrid(ix,jy,kz)
573                else ! concentration is zero
574                  sp_zer=.true.
575                endif
576              end do
577            end do
578          end do
579          write(unitoutgridppt) sp_count_i
580          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
581          write(unitoutgridppt) sp_count_r
582          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
583!       write(unitoutgridppt) sp_count_u
584!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
585
586        endif ! output for ppt
587
588      end do
589    end do
590
591    close(unitoutgridppt)
592    close(unitoutgrid)
593
594  end do
595
596  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
597  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
598       wetgridtotal
599  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
600       drygridtotal
601
602! Dump of receptor concentrations
603
604  if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
605    write(unitoutreceptppt) itime
606    do ks=1,nspec
607      write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
608           weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
609    end do
610  endif
611
612! Dump of receptor concentrations
613
614  if (numreceptor.gt.0) then
615    write(unitoutrecept) itime
616    do ks=1,nspec
617      write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
618           i=1,numreceptor)
619    end do
620  endif
621
622
623
624! Reinitialization of grid
625!*************************
626
627  do ks=1,nspec
628    do kp=1,maxpointspec_act
629      do i=1,numreceptor
630        creceptor(i,ks)=0.
631      end do
632      do jy=0,numygrid-1
633        do ix=0,numxgrid-1
634          do l=1,nclassunc
635            do nage=1,nageclass
636              do kz=1,numzgrid
637                gridunc(ix,jy,kz,ks,kp,l,nage)=0.
638              end do
639            end do
640          end do
641        end do
642      end do
643    end do
644  end do
645
646
647end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG