source: flexpart.git/src/concoutput.f90 @ 54cbd6c

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bugunivie
Last change on this file since 54cbd6c was 54cbd6c, checked in by Sabine <sabine.eckhardt@…>, 8 years ago

all f90 files for dry/wet backward mode

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