source: flexpart.git/src/concoutput_surf.f90 @ 8a65cb0

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

Added code, makefile for dev branch

  • Property mode set to 100644
File size: 22.9 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_surf(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  character :: adate*8,atime*6
100  character(len=3) :: anspec
101
102
103  if (verbosity.eq.1) then
104     print*,'inside concoutput_surf '
105     CALL SYSTEM_CLOCK(count_clock)
106     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
107  endif
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  !write(unitdates,'(a)') adate//atime
117
118    open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND')
119      write(unitdates,'(a)') adate//atime
120    close(unitdates)
121
122  ! For forward simulations, output fields have dimension MAXSPEC,
123  ! for backward simulations, output fields have dimension MAXPOINT.
124  ! Thus, make loops either about nspec, or about numpoint
125  !*****************************************************************
126
127
128  if (ldirect.eq.1) then
129    do ks=1,nspec
130      do kp=1,maxpointspec_act
131        tot_mu(ks,kp)=1
132      end do
133    end do
134  else
135    do ks=1,nspec
136      do kp=1,maxpointspec_act
137        tot_mu(ks,kp)=xmass(kp,ks)
138      end do
139    end do
140  endif
141
142
143  if (verbosity.eq.1) then
144     print*,'concoutput_surf 2'
145     CALL SYSTEM_CLOCK(count_clock)
146     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
147  endif
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  do kz=1,numzgrid
157    if (kz.eq.1) then
158      halfheight=outheight(1)/2.
159    else
160      halfheight=(outheight(kz)+outheight(kz-1))/2.
161    endif
162    do kzz=2,nz
163      if ((height(kzz-1).lt.halfheight).and. &
164           (height(kzz).gt.halfheight)) goto 46
165    end do
16646   kzz=max(min(kzz,nz),2)
167    dz1=halfheight-height(kzz-1)
168    dz2=height(kzz)-halfheight
169    dz=dz1+dz2
170    do jy=0,numygrid-1
171      do ix=0,numxgrid-1
172        xl=outlon0+real(ix)*dxout
173        yl=outlat0+real(jy)*dyout
174        xl=(xl-xlon0)/dx
175        yl=(yl-ylat0)/dy
176        iix=max(min(nint(xl),nxmin1),0)
177        jjy=max(min(nint(yl),nymin1),0)
178        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
179             rho(iix,jjy,kzz-1,2)*dz2)/dz
180      end do
181    end do
182  end do
183
184    do i=1,numreceptor
185      xl=xreceptor(i)
186      yl=yreceptor(i)
187      iix=max(min(nint(xl),nxmin1),0)
188      jjy=max(min(nint(yl),nymin1),0)
189      densityoutrecept(i)=rho(iix,jjy,1,2)
190    end do
191
192
193  ! Output is different for forward and backward simulations
194    do kz=1,numzgrid
195      do jy=0,numygrid-1
196        do ix=0,numxgrid-1
197          if (ldirect.eq.1) then
198            factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
199          else
200            factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
201          endif
202        end do
203      end do
204    end do
205
206  !*********************************************************************
207  ! Determine the standard deviation of the mean concentration or mixing
208  ! ratio (uncertainty of the output) and the dry and wet deposition
209  !*********************************************************************
210
211  if (verbosity.eq.1) then
212     print*,'concoutput_surf 3 (sd)'
213     CALL SYSTEM_CLOCK(count_clock)
214     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
215  endif
216  gridtotal=0.
217  gridsigmatotal=0.
218  gridtotalunc=0.
219  wetgridtotal=0.
220  wetgridsigmatotal=0.
221  wetgridtotalunc=0.
222  drygridtotal=0.
223  drygridsigmatotal=0.
224  drygridtotalunc=0.
225
226  do ks=1,nspec
227
228  write(anspec,'(i3.3)') ks
229  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
230    if (ldirect.eq.1) then
231      open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
232           atime//'_'//anspec,form='unformatted')
233    else
234      open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
235           atime//'_'//anspec,form='unformatted')
236    endif
237     write(unitoutgrid) itime
238   endif
239
240  if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
241   open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
242        atime//'_'//anspec,form='unformatted')
243
244    write(unitoutgridppt) itime
245  endif
246
247  do kp=1,maxpointspec_act
248  do nage=1,nageclass
249
250    do jy=0,numygrid-1
251      do ix=0,numxgrid-1
252
253  ! WET DEPOSITION
254        if ((WETDEP).and.(ldirect.gt.0)) then
255            do l=1,nclassunc
256              auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
257            end do
258            call mean(auxgrid,wetgrid(ix,jy), &
259                 wetgridsigma(ix,jy),nclassunc)
260  ! Multiply by number of classes to get total concentration
261            wetgrid(ix,jy)=wetgrid(ix,jy) &
262                 *nclassunc
263            wetgridtotal=wetgridtotal+wetgrid(ix,jy)
264  ! Calculate standard deviation of the mean
265            wetgridsigma(ix,jy)= &
266                 wetgridsigma(ix,jy)* &
267                 sqrt(real(nclassunc))
268            wetgridsigmatotal=wetgridsigmatotal+ &
269                 wetgridsigma(ix,jy)
270        endif
271
272  ! DRY DEPOSITION
273        if ((DRYDEP).and.(ldirect.gt.0)) then
274            do l=1,nclassunc
275              auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
276            end do
277            call mean(auxgrid,drygrid(ix,jy), &
278                 drygridsigma(ix,jy),nclassunc)
279  ! Multiply by number of classes to get total concentration
280            drygrid(ix,jy)=drygrid(ix,jy)* &
281                 nclassunc
282            drygridtotal=drygridtotal+drygrid(ix,jy)
283  ! Calculate standard deviation of the mean
284            drygridsigma(ix,jy)= &
285                 drygridsigma(ix,jy)* &
286                 sqrt(real(nclassunc))
287125         drygridsigmatotal=drygridsigmatotal+ &
288                 drygridsigma(ix,jy)
289        endif
290
291  ! CONCENTRATION OR MIXING RATIO
292        do kz=1,numzgrid
293            do l=1,nclassunc
294              auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
295            end do
296            call mean(auxgrid,grid(ix,jy,kz), &
297                 gridsigma(ix,jy,kz),nclassunc)
298  ! Multiply by number of classes to get total concentration
299            grid(ix,jy,kz)= &
300                 grid(ix,jy,kz)*nclassunc
301            gridtotal=gridtotal+grid(ix,jy,kz)
302  ! Calculate standard deviation of the mean
303            gridsigma(ix,jy,kz)= &
304                 gridsigma(ix,jy,kz)* &
305                 sqrt(real(nclassunc))
306            gridsigmatotal=gridsigmatotal+ &
307                 gridsigma(ix,jy,kz)
308        end do
309      end do
310    end do
311
312
313  !*******************************************************************
314  ! Generate output: may be in concentration (ng/m3) or in mixing
315  ! ratio (ppt) or both
316  ! Output the position and the values alternated multiplied by
317  ! 1 or -1, first line is number of values, number of positions
318  ! For backward simulations, the unit is seconds, stored in grid_time
319  !*******************************************************************
320
321  if (verbosity.eq.1) then
322     print*,'concoutput_surf 4 (output)'
323     CALL SYSTEM_CLOCK(count_clock)
324     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
325  endif
326
327  ! Concentration output
328  !*********************
329
330  if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
331
332  if (verbosity.eq.1) then
333     print*,'concoutput_surf (Wet deposition)'
334     CALL SYSTEM_CLOCK(count_clock)
335     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
336  endif
337
338  ! Wet deposition
339         sp_count_i=0
340         sp_count_r=0
341         sp_fact=-1.
342         sp_zer=.true.
343         if ((ldirect.eq.1).and.(WETDEP)) then
344         do jy=0,numygrid-1
345            do ix=0,numxgrid-1
346  ! concentraion greater zero
347              if (wetgrid(ix,jy).gt.smallnum) then
348                 if (sp_zer.eqv..true.) then ! first non zero value
349                    sp_count_i=sp_count_i+1
350                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
351                    sp_zer=.false.
352                    sp_fact=sp_fact*(-1.)
353                 endif
354                 sp_count_r=sp_count_r+1
355                 sparse_dump_r(sp_count_r)= &
356                      sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
357                 sparse_dump_u(sp_count_r)= &
358                      1.e12*wetgridsigma(ix,jy)/area(ix,jy)
359              else ! concentration is zero
360                  sp_zer=.true.
361              endif
362            end do
363         end do
364         else
365            sp_count_i=0
366            sp_count_r=0
367         endif
368         write(unitoutgrid) sp_count_i
369         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
370         write(unitoutgrid) sp_count_r
371         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
372!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
373
374  if (verbosity.eq.1) then
375     print*,'concoutput_surf (Dry deposition)'
376     CALL SYSTEM_CLOCK(count_clock)
377     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
378  endif
379  ! Dry deposition
380         sp_count_i=0
381         sp_count_r=0
382         sp_fact=-1.
383         sp_zer=.true.
384         if ((ldirect.eq.1).and.(DRYDEP)) then
385          do jy=0,numygrid-1
386            do ix=0,numxgrid-1
387              if (drygrid(ix,jy).gt.smallnum) then
388                 if (sp_zer.eqv..true.) then ! first non zero value
389                    sp_count_i=sp_count_i+1
390                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
391                    sp_zer=.false.
392                    sp_fact=sp_fact*(-1.)
393                 endif
394                 sp_count_r=sp_count_r+1
395                 sparse_dump_r(sp_count_r)= &
396                      sp_fact* &
397                      1.e12*drygrid(ix,jy)/area(ix,jy)
398                  sparse_dump_u(sp_count_r)= &
399                      1.e12*drygridsigma(ix,jy)/area(ix,jy)
400              else ! concentration is zero
401                  sp_zer=.true.
402              endif
403            end do
404          end do
405         else
406            sp_count_i=0
407            sp_count_r=0
408         endif
409         write(unitoutgrid) sp_count_i
410         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
411         write(unitoutgrid) sp_count_r
412         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
413!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
414
415  if (verbosity.eq.1) then
416     print*,'concoutput_surf (Concentrations)'
417     CALL SYSTEM_CLOCK(count_clock)
418     WRITE(*,*) 'SYSTEM_CLOCK',count_clock - count_clock0   
419  endif
420
421  ! Concentrations
422
423  ! surf_only write only 1st layer
424
425         sp_count_i=0
426         sp_count_r=0
427         sp_fact=-1.
428         sp_zer=.true.
429          do kz=1,1
430            do jy=0,numygrid-1
431              do ix=0,numxgrid-1
432                if (grid(ix,jy,kz).gt.smallnum) then
433                  if (sp_zer.eqv..true.) then ! first non zero value
434                    sp_count_i=sp_count_i+1
435                    sparse_dump_i(sp_count_i)= &
436                         ix+jy*numxgrid+kz*numxgrid*numygrid
437                    sp_zer=.false.
438                    sp_fact=sp_fact*(-1.)
439                  endif
440                   sp_count_r=sp_count_r+1
441                   sparse_dump_r(sp_count_r)= &
442                        sp_fact* &
443                        grid(ix,jy,kz)* &
444                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
445  !                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
446  !    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
447                   sparse_dump_u(sp_count_r)= &
448                        gridsigma(ix,jy,kz)* &
449                        factor3d(ix,jy,kz)/tot_mu(ks,kp)
450                else ! concentration is zero
451                  sp_zer=.true.
452                endif
453             end do
454           end do
455         end do
456         write(unitoutgrid) sp_count_i
457         write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
458         write(unitoutgrid) sp_count_r
459         write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
460!         write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
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)/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) (sparse_dump_u(i),i=1,sp_count_r)
505
506
507  ! Dry deposition
508         sp_count_i=0
509         sp_count_r=0
510         sp_fact=-1.
511         sp_zer=.true.
512         if ((ldirect.eq.1).and.(DRYDEP)) then
513          do jy=0,numygrid-1
514            do ix=0,numxgrid-1
515                if (drygrid(ix,jy).gt.smallnum) then
516                  if (sp_zer.eqv..true.) then ! first non zero value
517                    sp_count_i=sp_count_i+1
518                    sparse_dump_i(sp_count_i)= &
519                         ix+jy*numxgrid
520                    sp_zer=.false.
521                    sp_fact=sp_fact*(-1)
522                 endif
523                 sp_count_r=sp_count_r+1
524                 sparse_dump_r(sp_count_r)= &
525                      sp_fact* &
526                      1.e12*drygrid(ix,jy)/area(ix,jy)
527                 sparse_dump_u(sp_count_r)= &
528                      1.e12*drygridsigma(ix,jy)/area(ix,jy)
529              else ! concentration is zero
530                  sp_zer=.true.
531              endif
532            end do
533          end do
534         else
535           sp_count_i=0
536           sp_count_r=0
537         endif
538         write(unitoutgridppt) sp_count_i
539         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
540         write(unitoutgridppt) sp_count_r
541         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
542!         write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
543
544
545  ! Mixing ratios
546
547  ! surf_only write only 1st layer
548
549         sp_count_i=0
550         sp_count_r=0
551         sp_fact=-1.
552         sp_zer=.true.
553          do kz=1,1
554            do jy=0,numygrid-1
555              do ix=0,numxgrid-1
556                if (grid(ix,jy,kz).gt.smallnum) then
557                  if (sp_zer.eqv..true.) then ! first non zero value
558                    sp_count_i=sp_count_i+1
559                    sparse_dump_i(sp_count_i)= &
560                         ix+jy*numxgrid+kz*numxgrid*numygrid
561                    sp_zer=.false.
562                    sp_fact=sp_fact*(-1.)
563                 endif
564                 sp_count_r=sp_count_r+1
565                 sparse_dump_r(sp_count_r)= &
566                      sp_fact* &
567                      1.e12*grid(ix,jy,kz) &
568                      /volume(ix,jy,kz)/outnum* &
569                      weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
570                 sparse_dump_u(sp_count_r)= &
571                      1.e12*gridsigma(ix,jy,kz)/volume(ix,jy,kz)/ &
572                      outnum*weightair/weightmolar(ks)/ &
573                      densityoutgrid(ix,jy,kz)
574              else ! concentration is zero
575                  sp_zer=.true.
576              endif
577              end do
578            end do
579          end do
580         write(unitoutgridppt) sp_count_i
581         write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
582         write(unitoutgridppt) sp_count_r
583         write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
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_surf
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG