source: flexpart.git/src/concoutput.f90 @ 7999df47

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

Completed handling of nested wind fields with cloud water (for wet deposition).

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