source: flexpart.git/src/concoutput_mpi.f90 @ b069789

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

If the 'dates' file exists on simulation start, backup to 'dates.bak' and create new file 'dates' instead of appending to it

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