source: flexpart.git/src/concoutput_mpi.f90 @ 02095e3

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

Added option to use double precision for calculating the deposition fields

  • Property mode set to 100644
File size: 24.0 KB
RevLine 
[8a65cb0]1!**********************************************************************
2! Copyright 1998,1999,2000,2001,2002,2005,2007,2008,2009,2010         *
3! Andreas Stohl, Petra Seibert, A. Frank, Gerhard Wotawa,             *
4! Caroline Forster, Sabine Eckhardt, John Burkhart, Harald Sodemann   *
5!                                                                     *
6! This file is part of FLEXPART.                                      *
7!                                                                     *
8! FLEXPART is free software: you can redistribute it and/or modify    *
9! it under the terms of the GNU General Public License as published by*
10! the Free Software Foundation, either version 3 of the License, or   *
11! (at your option) any later version.                                 *
12!                                                                     *
13! FLEXPART is distributed in the hope that it will be useful,         *
14! but WITHOUT ANY WARRANTY; without even the implied warranty of      *
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the       *
16! GNU General Public License for more details.                        *
17!                                                                     *
18! You should have received a copy of the GNU General Public License   *
19! along with FLEXPART.  If not, see <http://www.gnu.org/licenses/>.   *
20!**********************************************************************
21
22subroutine concoutput(itime,outnum,gridtotalunc,wetgridtotalunc, &
23     drygridtotalunc)
24!                        i     i          o             o
25!       o
26!*****************************************************************************
27!                                                                            *
28!     Output of the concentration grid and the receptor concentrations.      *
29!                                                                            *
30!     Author: A. Stohl                                                       *
31!                                                                            *
32!     24 May 1995                                                            *
33!                                                                            *
34!     13 April 1999, Major update: if output size is smaller, dump output    *
35!                    in sparse matrix format; additional output of           *
36!                    uncertainty                                             *
37!                                                                            *
38!     05 April 2000, Major update: output of age classes; output for backward*
39!                    runs is time spent in grid cell times total mass of     *
40!                    species.                                                *
41!                                                                            *
42!     17 February 2002, Appropriate dimensions for backward and forward runs *
43!                       are now specified in file par_mod                    *
44!                                                                            *
45!     June 2006, write grid in sparse matrix with a single write command     *
46!                in order to save disk space                                 *
47!                                                                            *
48!     2008 new sparse matrix format                                          *
49!                                                                            *
50!     Changes eso:                                                           *
51!     2014        MPI version    This routine is only called by root MPI     *
52!                                process (the other processes have sent      *
53!                                their fields to root)                       *
54!*****************************************************************************
55!                                                                            *
56! Variables:                                                                 *
57! outnum          number of samples                                          *
58! ncells          number of cells with non-zero concentrations               *
59! sparse          .true. if in sparse matrix format, else .false.            *
60! tot_mu          1 for forward, initial mass mixing ration for backw. runs  *
61!                                                                            *
62!*****************************************************************************
63
64  use unc_mod
65  use point_mod
66  use outg_mod
67  use par_mod
68  use com_mod
69  use mpi_mod
[6a678e3]70  use mean_mod
[8a65cb0]71
72  implicit none
73
74  real(kind=dp) :: jul
75  integer :: itime,i,ix,jy,kz,ks,kp,l,iix,jjy,kzz,nage,jjjjmmdd,ihmmss
76  integer :: sp_count_i,sp_count_r
77  real :: sp_fact
78  real :: outnum,densityoutrecept(maxreceptor),xl,yl
79
80!real densityoutgrid(0:numxgrid-1,0:numygrid-1,numzgrid),
81!    +grid(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,maxpointspec_act,
82!    +    maxageclass)
83!real wetgrid(0:numxgrid-1,0:numygrid-1,maxspec,maxpointspec_act,
84!    +       maxageclass)
85!real drygrid(0:numxgrid-1,0:numygrid-1,maxspec,
86!    +       maxpointspec_act,maxageclass)
87!real gridsigma(0:numxgrid-1,0:numygrid-1,numzgrid,maxspec,
88!    +       maxpointspec_act,maxageclass),
89!    +     drygridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
90!    +     maxpointspec_act,maxageclass),
91!    +     wetgridsigma(0:numxgrid-1,0:numygrid-1,maxspec,
92!    +     maxpointspec_act,maxageclass)
93!real factor(0:numxgrid-1,0:numygrid-1,numzgrid)
94!real sparse_dump_r(numxgrid*numygrid*numzgrid)
95!integer sparse_dump_i(numxgrid*numygrid*numzgrid)
96
97!real sparse_dump_u(numxgrid*numygrid*numzgrid)
[6a678e3]98  real(dep_prec) :: auxgrid(nclassunc)
99  real(sp) :: gridtotal,gridsigmatotal,gridtotalunc
100  real(dep_prec) :: wetgridtotal,wetgridsigmatotal,wetgridtotalunc
101  real(dep_prec) :: drygridtotal,drygridsigmatotal,drygridtotalunc
[8a65cb0]102  real :: halfheight,dz,dz1,dz2,tot_mu(maxspec,maxpointspec_act)
103  real,parameter :: smallnum = tiny(0.0) ! smallest number that can be handled
104  real,parameter :: weightair=28.97
105  logical :: sp_zer
[478e9e6]106  LOGICAL,save :: init=.true.
[8a65cb0]107  character :: adate*8,atime*6
108  character(len=3) :: anspec
[5f9d14a]109  integer :: mind
[b069789]110! mind        eso: added to ensure identical results between 2&3-fields versions
111  character(LEN=8),save :: file_stat='REPLACE'
112  logical :: ldates_file
113  integer :: ierr
114  character(LEN=100) :: dates_char
115!  character :: dates_char
[8a65cb0]116
117! Measure execution time
118  if (mp_measure_time) call mpif_mtime('rootonly',0)
119
120! Determine current calendar date, needed for the file name
121!**********************************************************
122
123  jul=bdate+real(itime,kind=dp)/86400._dp
124  call caldate(jul,jjjjmmdd,ihmmss)
125  write(adate,'(i8.8)') jjjjmmdd
126  write(atime,'(i6.6)') ihmmss
[b069789]127
128! Overwrite existing dates file on first call, later append to it
129! This fixes a bug where the dates file kept growing across multiple runs
130
131! If 'dates' file exists, make a backup
132  inquire(file=path(2)(1:length(2))//'dates', exist=ldates_file)
133  if (ldates_file.and.init) then
134    open(unit=unitdates, file=path(2)(1:length(2))//'dates',form='formatted', &
135         &access='sequential', status='old', action='read', iostat=ierr)
136    open(unit=unittmp, file=path(2)(1:length(2))//'dates.bak', access='sequential', &
137         &status='replace', action='write', form='formatted', iostat=ierr)
138    do while (.true.)
139      read(unitdates, '(a)', iostat=ierr) dates_char
140      if (ierr.ne.0) exit
141      write(unit=unittmp, fmt='(a)', iostat=ierr, advance='yes') trim(dates_char)
142!      if (ierr.ne.0) write(*,*) "Write error, ", ierr
143    end do
144    close(unit=unitdates)
145    close(unit=unittmp)
146  end if
147
148  open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND', STATUS=file_stat)
[8a65cb0]149  write(unitdates,'(a)') adate//atime
150  close(unitdates) 
151
[478e9e6]152  ! Overwrite existing dates file on first call, later append to it
153  ! This fixes a bug where the dates file kept growing across multiple runs
154  IF (init) THEN
155    file_stat='OLD'
156    init=.false.
157  END IF
158
159
[8a65cb0]160! For forward simulations, output fields have dimension MAXSPEC,
161! for backward simulations, output fields have dimension MAXPOINT.
162! Thus, make loops either about nspec, or about numpoint
163!*****************************************************************
164
165
166  if (ldirect.eq.1) then
167    do ks=1,nspec
168      do kp=1,maxpointspec_act
169        tot_mu(ks,kp)=1
170      end do
171    end do
172  else
173    do ks=1,nspec
174      do kp=1,maxpointspec_act
175        tot_mu(ks,kp)=xmass(kp,ks)
176      end do
177    end do
178  endif
179
180
181!*******************************************************************
182! Compute air density: sufficiently accurate to take it
183! from coarse grid at some time
184! Determine center altitude of output layer, and interpolate density
185! data to that altitude
186!*******************************************************************
187
[5f9d14a]188  mind=memind(2)
[8a65cb0]189  do kz=1,numzgrid
190    if (kz.eq.1) then
191      halfheight=outheight(1)/2.
192    else
193      halfheight=(outheight(kz)+outheight(kz-1))/2.
194    endif
195    do kzz=2,nz
196      if ((height(kzz-1).lt.halfheight).and. &
197           (height(kzz).gt.halfheight)) goto 46
198    end do
19946  kzz=max(min(kzz,nz),2)
200    dz1=halfheight-height(kzz-1)
201    dz2=height(kzz)-halfheight
202    dz=dz1+dz2
203    do jy=0,numygrid-1
204      do ix=0,numxgrid-1
205        xl=outlon0+real(ix)*dxout
206        yl=outlat0+real(jy)*dyout
207        xl=(xl-xlon0)/dx
208        yl=(yl-ylat0)/dy !v9.1.1
209        iix=max(min(nint(xl),nxmin1),0)
210        jjy=max(min(nint(yl),nymin1),0)
[5f9d14a]211        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
212        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
213        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
214             rho(iix,jjy,kzz-1,mind)*dz2)/dz
[8a65cb0]215      end do
216    end do
217  end do
218
219  do i=1,numreceptor
220    xl=xreceptor(i)
221    yl=yreceptor(i)
222    iix=max(min(nint(xl),nxmin1),0)
223    jjy=max(min(nint(yl),nymin1),0)
[5f9d14a]224    !densityoutrecept(i)=rho(iix,jjy,1,2)
225    densityoutrecept(i)=rho(iix,jjy,1,mind)
[8a65cb0]226  end do
227
228
229! Output is different for forward and backward simulations
230  do kz=1,numzgrid
231    do jy=0,numygrid-1
232      do ix=0,numxgrid-1
233        if (ldirect.eq.1) then
234          factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
235        else
236          factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
237        endif
238      end do
239    end do
240  end do
241
242!*********************************************************************
243! Determine the standard deviation of the mean concentration or mixing
244! ratio (uncertainty of the output) and the dry and wet deposition
245!*********************************************************************
246
247  gridtotal=0.
248  gridsigmatotal=0.
249  gridtotalunc=0.
250  wetgridtotal=0.
251  wetgridsigmatotal=0.
252  wetgridtotalunc=0.
253  drygridtotal=0.
254  drygridsigmatotal=0.
255  drygridtotalunc=0.
256
257  do ks=1,nspec
258
259    write(anspec,'(i3.3)') ks
260    if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
261      if (ldirect.eq.1) then
262        open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
263             atime//'_'//anspec,form='unformatted')
264      else
265        open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
266             atime//'_'//anspec,form='unformatted')
267      endif
268      write(unitoutgrid) itime
269    endif
270
271    if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
272      open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
273           atime//'_'//anspec,form='unformatted')
274
275      write(unitoutgridppt) itime
276    endif
277
278    do kp=1,maxpointspec_act
279      do nage=1,nageclass
280
281        do jy=0,numygrid-1
282          do ix=0,numxgrid-1
283
284! WET DEPOSITION
285            if ((WETDEP).and.(ldirect.gt.0)) then
286              do l=1,nclassunc
287                auxgrid(l)=wetgridunc0(ix,jy,ks,kp,l,nage)
288              end do
289              call mean(auxgrid,wetgrid(ix,jy), &
290                   wetgridsigma(ix,jy),nclassunc)
291! Multiply by number of classes to get total concentration
292              wetgrid(ix,jy)=wetgrid(ix,jy) &
293                   *nclassunc
294              wetgridtotal=wetgridtotal+wetgrid(ix,jy)
295! Calculate standard deviation of the mean
296              wetgridsigma(ix,jy)= &
297                   wetgridsigma(ix,jy)* &
298                   sqrt(real(nclassunc))
299              wetgridsigmatotal=wetgridsigmatotal+ &
300                   wetgridsigma(ix,jy)
301            endif
302
303! DRY DEPOSITION
304            if ((DRYDEP).and.(ldirect.gt.0)) then
305              do l=1,nclassunc
306                auxgrid(l)=drygridunc0(ix,jy,ks,kp,l,nage)
307              end do
308              call mean(auxgrid,drygrid(ix,jy), &
309                   drygridsigma(ix,jy),nclassunc)
310! Multiply by number of classes to get total concentration
311              drygrid(ix,jy)=drygrid(ix,jy)* &
312                   nclassunc
313              drygridtotal=drygridtotal+drygrid(ix,jy)
314! Calculate standard deviation of the mean
315              drygridsigma(ix,jy)= &
316                   drygridsigma(ix,jy)* &
317                   sqrt(real(nclassunc))
318              drygridsigmatotal=drygridsigmatotal+ &
319                   drygridsigma(ix,jy)
320            endif
321
322! CONCENTRATION OR MIXING RATIO
323            do kz=1,numzgrid
324              do l=1,nclassunc
325                auxgrid(l)=gridunc(ix,jy,kz,ks,kp,l,nage)
326              end do
327              call mean(auxgrid,grid(ix,jy,kz), &
328                   gridsigma(ix,jy,kz),nclassunc)
329! Multiply by number of classes to get total concentration
330              grid(ix,jy,kz)= &
331                   grid(ix,jy,kz)*nclassunc
332              gridtotal=gridtotal+grid(ix,jy,kz)
333! Calculate standard deviation of the mean
334              gridsigma(ix,jy,kz)= &
335                   gridsigma(ix,jy,kz)* &
336                   sqrt(real(nclassunc))
337              gridsigmatotal=gridsigmatotal+ &
338                   gridsigma(ix,jy,kz)
339            end do
340          end do
341        end do
342
343
344
345
346!*******************************************************************
347! Generate output: may be in concentration (ng/m3) or in mixing
348! ratio (ppt) or both
349! Output the position and the values alternated multiplied by
350! 1 or -1, first line is number of values, number of positions
351! For backward simulations, the unit is seconds, stored in grid_time
352!*******************************************************************
353
354! Concentration output
355!*********************
356        if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
357
358! Wet deposition
359          sp_count_i=0
360          sp_count_r=0
361          sp_fact=-1.
362          sp_zer=.true.
363          if ((ldirect.eq.1).and.(WETDEP)) then
364            do jy=0,numygrid-1
365              do ix=0,numxgrid-1
366!oncentraion greater zero
367                if (wetgrid(ix,jy).gt.smallnum) then
368                  if (sp_zer.eqv..true.) then ! first non zero value
369                    sp_count_i=sp_count_i+1
370                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
371                    sp_zer=.false.
372                    sp_fact=sp_fact*(-1.)
373                  endif
374                  sp_count_r=sp_count_r+1
375                  sparse_dump_r(sp_count_r)= &
376                       sp_fact*1.e12*wetgrid(ix,jy)/area(ix,jy)
377!                sparse_dump_u(sp_count_r)=
378!+                1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
379                else ! concentration is zero
380                  sp_zer=.true.
381                endif
382              end do
383            end do
384          else
385            sp_count_i=0
386            sp_count_r=0
387          endif
388          write(unitoutgrid) sp_count_i
389          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
390          write(unitoutgrid) sp_count_r
391          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
392!       write(unitoutgrid) sp_count_u
393!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
394
395! Dry deposition
396          sp_count_i=0
397          sp_count_r=0
398          sp_fact=-1.
399          sp_zer=.true.
400          if ((ldirect.eq.1).and.(DRYDEP)) then
401            do jy=0,numygrid-1
402              do ix=0,numxgrid-1
403                if (drygrid(ix,jy).gt.smallnum) then
404                  if (sp_zer.eqv..true.) then ! first non zero value
405                    sp_count_i=sp_count_i+1
406                    sparse_dump_i(sp_count_i)=ix+jy*numxgrid
407                    sp_zer=.false.
408                    sp_fact=sp_fact*(-1.)
409                  endif
410                  sp_count_r=sp_count_r+1
411                  sparse_dump_r(sp_count_r)= &
412                       sp_fact* &
413                       1.e12*drygrid(ix,jy)/area(ix,jy)
414!                sparse_dump_u(sp_count_r)=
415!+                1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
416                else ! concentration is zero
417                  sp_zer=.true.
418                endif
419              end do
420            end do
421          else
422            sp_count_i=0
423            sp_count_r=0
424          endif
425          write(unitoutgrid) sp_count_i
426          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
427          write(unitoutgrid) sp_count_r
428          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
429!       write(*,*) sp_count_u
430!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
431
432
433
434! Concentrations
435          sp_count_i=0
436          sp_count_r=0
437          sp_fact=-1.
438          sp_zer=.true.
439          do kz=1,numzgrid
440            do jy=0,numygrid-1
441              do ix=0,numxgrid-1
442                if (grid(ix,jy,kz).gt.smallnum) then
443                  if (sp_zer.eqv..true.) then ! first non zero value
444                    sp_count_i=sp_count_i+1
445                    sparse_dump_i(sp_count_i)= &
446                         ix+jy*numxgrid+kz*numxgrid*numygrid
447                    sp_zer=.false.
448                    sp_fact=sp_fact*(-1.)
449                  endif
450                  sp_count_r=sp_count_r+1
451                  sparse_dump_r(sp_count_r)= &
452                       sp_fact* &
453                       grid(ix,jy,kz)* &
454                       factor3d(ix,jy,kz)/tot_mu(ks,kp)
455!                 if ((factor(ix,jy,kz)/tot_mu(ks,kp)).eq.0)
456!    +              write (*,*) factor(ix,jy,kz),tot_mu(ks,kp),ks,kp
457!                sparse_dump_u(sp_count_r)=
458!+               ,gridsigma(ix,jy,kz,ks,kp,nage)*
459!+               factor(ix,jy,kz)/tot_mu(ks,kp)
460                else ! concentration is zero
461                  sp_zer=.true.
462                endif
463              end do
464            end do
465          end do
466          write(unitoutgrid) sp_count_i
467          write(unitoutgrid) (sparse_dump_i(i),i=1,sp_count_i)
468          write(unitoutgrid) sp_count_r
469          write(unitoutgrid) (sparse_dump_r(i),i=1,sp_count_r)
470!       write(unitoutgrid) sp_count_u
471!       write(unitoutgrid) (sparse_dump_u(i),i=1,sp_count_r)
472
473
474
475        endif !  concentration output
476
477! Mixing ratio output
478!********************
479
480        if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
481
482! Wet deposition
483          sp_count_i=0
484          sp_count_r=0
485          sp_fact=-1.
486          sp_zer=.true.
487          if ((ldirect.eq.1).and.(WETDEP)) then
488            do jy=0,numygrid-1
489              do ix=0,numxgrid-1
490                if (wetgrid(ix,jy).gt.smallnum) then
491                  if (sp_zer.eqv..true.) then ! first non zero value
492                    sp_count_i=sp_count_i+1
493                    sparse_dump_i(sp_count_i)= &
494                         ix+jy*numxgrid
495                    sp_zer=.false.
496                    sp_fact=sp_fact*(-1.)
497                  endif
498                  sp_count_r=sp_count_r+1
499                  sparse_dump_r(sp_count_r)= &
500                       sp_fact* &
501                       1.e12*wetgrid(ix,jy)/area(ix,jy)
502!                sparse_dump_u(sp_count_r)=
503!    +            ,1.e12*wetgridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
504                else ! concentration is zero
505                  sp_zer=.true.
506                endif
507              end do
508            end do
509          else
510            sp_count_i=0
511            sp_count_r=0
512          endif
513          write(unitoutgridppt) sp_count_i
514          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
515          write(unitoutgridppt) sp_count_r
516          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
517!       write(unitoutgridppt) sp_count_u
518!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
519
520
521! Dry deposition
522          sp_count_i=0
523          sp_count_r=0
524          sp_fact=-1.
525          sp_zer=.true.
526          if ((ldirect.eq.1).and.(DRYDEP)) then
527            do jy=0,numygrid-1
528              do ix=0,numxgrid-1
529                if (drygrid(ix,jy).gt.smallnum) then
530                  if (sp_zer.eqv..true.) then ! first non zero value
531                    sp_count_i=sp_count_i+1
532                    sparse_dump_i(sp_count_i)= &
533                         ix+jy*numxgrid
534                    sp_zer=.false.
535                    sp_fact=sp_fact*(-1)
536                  endif
537                  sp_count_r=sp_count_r+1
538                  sparse_dump_r(sp_count_r)= &
539                       sp_fact* &
540                       1.e12*drygrid(ix,jy)/area(ix,jy)
541!                sparse_dump_u(sp_count_r)=
542!    +            ,1.e12*drygridsigma(ix,jy,ks,kp,nage)/area(ix,jy)
543                else ! concentration is zero
544                  sp_zer=.true.
545                endif
546              end do
547            end do
548          else
549            sp_count_i=0
550            sp_count_r=0
551          endif
552          write(unitoutgridppt) sp_count_i
553          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
554          write(unitoutgridppt) sp_count_r
555          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
556!       write(unitoutgridppt) sp_count_u
557!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
558
559
560! Mixing ratios
561          sp_count_i=0
562          sp_count_r=0
563          sp_fact=-1.
564          sp_zer=.true.
565          do kz=1,numzgrid
566            do jy=0,numygrid-1
567              do ix=0,numxgrid-1
568                if (grid(ix,jy,kz).gt.smallnum) then
569                  if (sp_zer.eqv..true.) then ! first non zero value
570                    sp_count_i=sp_count_i+1
571                    sparse_dump_i(sp_count_i)= &
572                         ix+jy*numxgrid+kz*numxgrid*numygrid
573                    sp_zer=.false.
574                    sp_fact=sp_fact*(-1.)
575                  endif
576                  sp_count_r=sp_count_r+1
577                  sparse_dump_r(sp_count_r)= &
578                       sp_fact* &
579                       1.e12*grid(ix,jy,kz) &
580                       /volume(ix,jy,kz)/outnum* &
581                       weightair/weightmolar(ks)/densityoutgrid(ix,jy,kz)
582!                sparse_dump_u(sp_count_r)=
583!+              ,1.e12*gridsigma(ix,jy,kz,ks,kp,nage)/volume(ix,jy,kz)/
584!+              outnum*weightair/weightmolar(ks)/
585!+              densityoutgrid(ix,jy,kz)
586                else ! concentration is zero
587                  sp_zer=.true.
588                endif
589              end do
590            end do
591          end do
592          write(unitoutgridppt) sp_count_i
593          write(unitoutgridppt) (sparse_dump_i(i),i=1,sp_count_i)
594          write(unitoutgridppt) sp_count_r
595          write(unitoutgridppt) (sparse_dump_r(i),i=1,sp_count_r)
596!       write(unitoutgridppt) sp_count_u
597!       write(unitoutgridppt) (sparse_dump_u(i),i=1,sp_count_r)
598
599        endif ! output for ppt
600
601      end do
602    end do
603
604    close(unitoutgridppt)
605    close(unitoutgrid)
606
607  end do
608
609  if (gridtotal.gt.0.) gridtotalunc=gridsigmatotal/gridtotal
610  if (wetgridtotal.gt.0.) wetgridtotalunc=wetgridsigmatotal/ &
611       wetgridtotal
612  if (drygridtotal.gt.0.) drygridtotalunc=drygridsigmatotal/ &
613       drygridtotal
614
615! Dump of receptor concentrations
616
617  if (numreceptor.gt.0 .and. (iout.eq.2 .or. iout.eq.3)  ) then
618    write(unitoutreceptppt) itime
619    do ks=1,nspec
620      write(unitoutreceptppt) (1.e12*creceptor(i,ks)/outnum* &
621           weightair/weightmolar(ks)/densityoutrecept(i),i=1,numreceptor)
622    end do
623  endif
624
625! Dump of receptor concentrations
626
627  if (numreceptor.gt.0) then
628    write(unitoutrecept) itime
629    do ks=1,nspec
630      write(unitoutrecept) (1.e12*creceptor(i,ks)/outnum, &
631           i=1,numreceptor)
632    end do
633  endif
634
635
636
637! Reinitialization of grid
638!*************************
639
640  do ks=1,nspec
641    do kp=1,maxpointspec_act
642      do i=1,numreceptor
643        creceptor(i,ks)=0.
644      end do
645      do jy=0,numygrid-1
646        do ix=0,numxgrid-1
647          do l=1,nclassunc
648            do nage=1,nageclass
649              do kz=1,numzgrid
650                gridunc(ix,jy,kz,ks,kp,l,nage)=0.
651              end do
652            end do
653          end do
654        end do
655      end do
656    end do
657  end do
658
659  if (mp_measure_time) call mpif_mtime('rootonly',1)
660 
661end subroutine concoutput
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG