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

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since dced13c was 462f74b, checked in by Sabine Eckhardt <sabine@…>, 8 years ago

first version of the backward scavenging

  • Property mode set to 100644
File size: 23.6 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
[462f74b]266    if (iout.eq.6) then !scavdep output
267      open(unitoutgrid,file=path(2)(1:length(2))//'grid_scav_'//adate// &
268           atime//'_'//anspec,form='unformatted')
269      write(unitoutgrid) itime
270     endif
271
[b069789]272    do kp=1,maxpointspec_act
273      do nage=1,nageclass
[e200b7a]274
[b069789]275        do jy=0,numygrid-1
276          do ix=0,numxgrid-1
[e200b7a]277
[b069789]278! WET DEPOSITION
279            if ((WETDEP).and.(ldirect.gt.0)) then
280              do l=1,nclassunc
281                auxgrid(l)=wetgridunc(ix,jy,ks,kp,l,nage)
282              end do
283              call mean(auxgrid,wetgrid(ix,jy), &
284                   wetgridsigma(ix,jy),nclassunc)
285! Multiply by number of classes to get total concentration
286              wetgrid(ix,jy)=wetgrid(ix,jy) &
287                   *nclassunc
288              wetgridtotal=wetgridtotal+wetgrid(ix,jy)
289! Calculate standard deviation of the mean
290              wetgridsigma(ix,jy)= &
291                   wetgridsigma(ix,jy)* &
292                   sqrt(real(nclassunc))
293              wetgridsigmatotal=wetgridsigmatotal+ &
294                   wetgridsigma(ix,jy)
295            endif
296
297! DRY DEPOSITION
298            if ((DRYDEP).and.(ldirect.gt.0)) then
299              do l=1,nclassunc
300                auxgrid(l)=drygridunc(ix,jy,ks,kp,l,nage)
301              end do
302              call mean(auxgrid,drygrid(ix,jy), &
303                   drygridsigma(ix,jy),nclassunc)
304! Multiply by number of classes to get total concentration
305              drygrid(ix,jy)=drygrid(ix,jy)* &
306                   nclassunc
307              drygridtotal=drygridtotal+drygrid(ix,jy)
308! Calculate standard deviation of the mean
309              drygridsigma(ix,jy)= &
310                   drygridsigma(ix,jy)* &
311                   sqrt(real(nclassunc))
312125           drygridsigmatotal=drygridsigmatotal+ &
313                   drygridsigma(ix,jy)
314            endif
315
316! CONCENTRATION OR MIXING RATIO
317            do kz=1,numzgrid
318              do l=1,nclassunc
319                auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
320              end do
321              call mean(auxgrid,grid(ix,jy,kz), &
322                   gridsigma(ix,jy,kz),nclassunc)
323! Multiply by number of classes to get total concentration
324              grid(ix,jy,kz)= &
325                   grid(ix,jy,kz)*nclassunc
326              gridtotal=gridtotal+grid(ix,jy,kz)
327! Calculate standard deviation of the mean
328              gridsigma(ix,jy,kz)= &
329                   gridsigma(ix,jy,kz)* &
330                   sqrt(real(nclassunc))
331              gridsigmatotal=gridsigmatotal+ &
332                   gridsigma(ix,jy,kz)
[e200b7a]333            end do
[b069789]334          end do
[e200b7a]335        end do
336
337
338
339
[b069789]340!*******************************************************************
341! Generate output: may be in concentration (ng/m3) or in mixing
342! ratio (ppt) or both
343! Output the position and the values alternated multiplied by
344! 1 or -1, first line is number of values, number of positions
345! For backward simulations, the unit is seconds, stored in grid_time
346!*******************************************************************
[e200b7a]347
[b069789]348! Concentration output
349!*********************
[462f74b]350        if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5).or.(iout.eq.6)) then
[e200b7a]351
[b069789]352! Wet deposition
353          sp_count_i=0
354          sp_count_r=0
355          sp_fact=-1.
356          sp_zer=.true.
357          if ((ldirect.eq.1).and.(WETDEP)) then
358            do jy=0,numygrid-1
359              do ix=0,numxgrid-1
360!oncentraion greater zero
361                if (wetgrid(ix,jy).gt.smallnum) then
362                  if (sp_zer.eqv..true.) then ! first non zero value
[e200b7a]363                    sp_count_i=sp_count_i+1
364                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
365                    sp_zer=.false.
366                    sp_fact=sp_fact*(-1.)
[b069789]367                  endif
368                  sp_count_r=sp_count_r+1
369                  sparse_dump_r(sp_count_r)= &
370                       sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
371!                sparse_dump_u(sp_count_r)=
372!+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
373                else ! concentration is zero
[e200b7a]374                  sp_zer=.true.
[b069789]375                endif
376              end do
[e200b7a]377            end do
[b069789]378          else
[e200b7a]379            sp_count_i=0
380            sp_count_r=0
[b069789]381          endif
382          write(unitoutgrid) sp_count_i
383          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
384          write(unitoutgrid) sp_count_r
385          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
386!       write(unitoutgrid) sp_count_u
387!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
388
389! Dry deposition
390          sp_count_i=0
391          sp_count_r=0
392          sp_fact=-1.
393          sp_zer=.true.
394          if ((ldirect.eq.1).and.(DRYDEP)) then
395            do jy=0,numygrid-1
396              do ix=0,numxgrid-1
397                if (drygrid(ix,jy).gt.smallnum) then
398                  if (sp_zer.eqv..true.) then ! first non zero value
[e200b7a]399                    sp_count_i=sp_count_i+1
400                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
401                    sp_zer=.false.
402                    sp_fact=sp_fact*(-1.)
[b069789]403                  endif
404                  sp_count_r=sp_count_r+1
405                  sparse_dump_r(sp_count_r)= &
406                       sp_fact* &
407                       1.e12*drygrid(ix,jy)/area(ix,jy)
408!                sparse_dump_u(sp_count_r)=
409!+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
410                else ! concentration is zero
[e200b7a]411                  sp_zer=.true.
[b069789]412                endif
413              end do
[e200b7a]414            end do
[b069789]415          else
[e200b7a]416            sp_count_i=0
417            sp_count_r=0
[b069789]418          endif
419          write(unitoutgrid) sp_count_i
420          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
421          write(unitoutgrid) sp_count_r
422          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
423!       write(*,*) sp_count_u
424!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
[e200b7a]425
426
427
[b069789]428! Concentrations
429          sp_count_i=0
430          sp_count_r=0
431          sp_fact=-1.
432          sp_zer=.true.
[e200b7a]433          do kz=1,numzgrid
434            do jy=0,numygrid-1
435              do ix=0,numxgrid-1
436                if (grid(ix,jy,kz).gt.smallnum) then
437                  if (sp_zer.eqv..true.) then ! first non zero value
438                    sp_count_i=sp_count_i+1
439                    sparse_dump_i(sp_count_i)= &
440                         ix+jy*numxgrid+kz*numxgrid*numygrid
441                    sp_zer=.false.
442                    sp_fact=sp_fact*(-1.)
[b069789]443                  endif
444                  sp_count_r=sp_count_r+1
445                  sparse_dump_r(sp_count_r)= &
446                       sp_fact* &
447                       grid(ix,jy,kz)* &
448                       factor3d(ix,jy,kz)/tot_mu(ks,kp)
449!                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
450!    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
451!                sparse_dump_u(sp_count_r)=
452!+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
453!+               factor(ix,jy,kz)/tot_mu(ks,kp)
454                else ! concentration is zero
[e200b7a]455                  sp_zer=.true.
[b069789]456                endif
[e200b7a]457              end do
458            end do
459          end do
[b069789]460          write(unitoutgrid) sp_count_i
461          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
462          write(unitoutgrid) sp_count_r
463          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
464!       write(unitoutgrid) sp_count_u
465!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
[e200b7a]466
467
468
[b069789]469        endif !  concentration output
[e200b7a]470
[b069789]471! Mixing ratio output
472!********************
[e200b7a]473
[b069789]474        if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
[e200b7a]475
[b069789]476! Wet deposition
477          sp_count_i=0
478          sp_count_r=0
479          sp_fact=-1.
480          sp_zer=.true.
481          if ((ldirect.eq.1).and.(WETDEP)) then
482            do jy=0,numygrid-1
483              do ix=0,numxgrid-1
[e200b7a]484                if (wetgrid(ix,jy).gt.smallnum) then
485                  if (sp_zer.eqv..true.) then ! first non zero value
486                    sp_count_i=sp_count_i+1
487                    sparse_dump_i(sp_count_i)= &
488                         ix+jy*numxgrid
489                    sp_zer=.false.
490                    sp_fact=sp_fact*(-1.)
[b069789]491                  endif
492                  sp_count_r=sp_count_r+1
493                  sparse_dump_r(sp_count_r)= &
494                       sp_fact* &
495                       1.e12*wetgrid(ix,jy)/area(ix,jy)
496!                sparse_dump_u(sp_count_r)=
497!    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
498                else ! concentration is zero
[e200b7a]499                  sp_zer=.true.
[b069789]500                endif
501              end do
[e200b7a]502            end do
[b069789]503          else
504            sp_count_i=0
505            sp_count_r=0
506          endif
507          write(unitoutgridppt) sp_count_i
508          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
509          write(unitoutgridppt) sp_count_r
510          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
511!       write(unitoutgridppt) sp_count_u
512!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
513
514
515! Dry deposition
516          sp_count_i=0
517          sp_count_r=0
518          sp_fact=-1.
519          sp_zer=.true.
520          if ((ldirect.eq.1).and.(DRYDEP)) then
521            do jy=0,numygrid-1
522              do ix=0,numxgrid-1
[e200b7a]523                if (drygrid(ix,jy).gt.smallnum) then
524                  if (sp_zer.eqv..true.) then ! first non zero value
525                    sp_count_i=sp_count_i+1
526                    sparse_dump_i(sp_count_i)= &
527                         ix+jy*numxgrid
528                    sp_zer=.false.
529                    sp_fact=sp_fact*(-1)
[b069789]530                  endif
531                  sp_count_r=sp_count_r+1
532                  sparse_dump_r(sp_count_r)= &
533                       sp_fact* &
534                       1.e12*drygrid(ix,jy)/area(ix,jy)
535!                sparse_dump_u(sp_count_r)=
536!    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
537                else ! concentration is zero
[e200b7a]538                  sp_zer=.true.
[b069789]539                endif
540              end do
[e200b7a]541            end do
[b069789]542          else
543            sp_count_i=0
544            sp_count_r=0
545          endif
546          write(unitoutgridppt) sp_count_i
547          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
548          write(unitoutgridppt) sp_count_r
549          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
550!       write(unitoutgridppt) sp_count_u
551!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
552
553
554! Mixing ratios
555          sp_count_i=0
556          sp_count_r=0
557          sp_fact=-1.
558          sp_zer=.true.
[e200b7a]559          do kz=1,numzgrid
560            do jy=0,numygrid-1
561              do ix=0,numxgrid-1
562                if (grid(ix,jy,kz).gt.smallnum) then
563                  if (sp_zer.eqv..true.) then ! first non zero value
564                    sp_count_i=sp_count_i+1
565                    sparse_dump_i(sp_count_i)= &
566                         ix+jy*numxgrid+kz*numxgrid*numygrid
567                    sp_zer=.false.
568                    sp_fact=sp_fact*(-1.)
[b069789]569                  endif
570                  sp_count_r=sp_count_r+1
571                  sparse_dump_r(sp_count_r)= &
572                       sp_fact* &
573                       1.e12*grid(ix,jy,kz) &
574                       /volume(ix,jy,kz)/outnum* &
575                       weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
576!                sparse_dump_u(sp_count_r)=
577!+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
578!+              outnum*weightair/weightmolar(ks)/
579!+              densityoutgrid(ix,jy,kz)
580                else ! concentration is zero
[e200b7a]581                  sp_zer=.true.
[b069789]582                endif
[e200b7a]583              end do
584            end do
585          end do
[b069789]586          write(unitoutgridppt) sp_count_i
587          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
588          write(unitoutgridppt) sp_count_r
589          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
590!       write(unitoutgridppt) sp_count_u
591!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
[e200b7a]592
[b069789]593        endif ! output for ppt
[e200b7a]594
[b069789]595      end do
596    end do
[e200b7a]597
598    close(unitoutgridppt)
599    close(unitoutgrid)
600
601  end do
602
603  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
604  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
605       wetgridtotal
606  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
607       drygridtotal
608
[b069789]609! Dump of receptor concentrations
[e200b7a]610
[b069789]611  if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
612    write(unitoutreceptppt) itime
613    do ks=1,nspec
614      write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
615           weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
616    end do
617  endif
[e200b7a]618
[b069789]619! Dump of receptor concentrations
[e200b7a]620
[b069789]621  if (numreceptor.gt.0) then
622    write(unitoutrecept) itime
623    do ks=1,nspec
624      write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
625           i=1,numreceptor)
626    end do
627  endif
[e200b7a]628
629
630
[b069789]631! Reinitialization of grid
632!*************************
[e200b7a]633
[db712a8]634  ! do ks=1,nspec
635  !   do kp=1,maxpointspec_act
636  !     do i=1,numreceptor
637  !       creceptor(i,ks)=0.
638  !     end do
639  !     do jy=0,numygrid-1
640  !       do ix=0,numxgrid-1
641  !         do l=1,nclassunc
642  !           do nage=1,nageclass
643  !             do kz=1,numzgrid
644  !               gridunc(ix,jy,kz,ks,kp,l,nage)=0.
645  !             end do
646  !           end do
647  !         end do
648  !       end do
649  !     end do
650  !   end do
651  ! end do
652  creceptor(:,:)=0.
653  gridunc(:,:,:,:,:,:,:)=0.
[e200b7a]654
655
656end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG