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

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

OH files converted from netcdf to binary format

  • 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  use mean_mod
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)
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
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
106  LOGICAL,save :: init=.true.
107  character :: adate*8,atime*6
108  character(len=3) :: anspec
109  integer :: mind
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
116! Measure execution time
117  if (mp_measure_time) call mpif_mtime('rootonly',0)
118
119! Determine current calendar date, needed for the file name
120!**********************************************************
121
122  jul=bdate+real(itime,kind=dp)/86400._dp
123  call caldate(jul,jjjjmmdd,ihmmss)
124  write(adate,'(i8.8)') jjjjmmdd
125  write(atime,'(i6.6)') ihmmss
126
127! Overwrite existing dates file on first call, later append to it
128! This fixes a bug where the dates file kept growing across multiple runs
129
130! If 'dates' file exists, make a backup
131  inquire(file=path(2)(1:length(2))//'dates', exist=ldates_file)
132  if (ldates_file.and.init) then
133    open(unit=unitdates, file=path(2)(1:length(2))//'dates',form='formatted', &
134         &access='sequential', status='old', action='read', iostat=ierr)
135    open(unit=unittmp, file=path(2)(1:length(2))//'dates.bak', access='sequential', &
136         &status='replace', action='write', form='formatted', iostat=ierr)
137    do while (.true.)
138      read(unitdates, '(a)', iostat=ierr) dates_char
139      if (ierr.ne.0) exit
140      write(unit=unittmp, fmt='(a)', iostat=ierr, advance='yes') trim(dates_char)
141!      if (ierr.ne.0) write(*,*) "Write error, ", ierr
142    end do
143    close(unit=unitdates)
144    close(unit=unittmp)
145  end if
146
147  open(unitdates,file=path(2)(1:length(2))//'dates', ACCESS='APPEND', STATUS=file_stat)
148  write(unitdates,'(a)') adate//atime
149  close(unitdates) 
150
151  ! Overwrite existing dates file on first call, later append to it
152  ! This fixes a bug where the dates file kept growing across multiple runs
153  IF (init) THEN
154    file_stat='OLD'
155    init=.false.
156  END IF
157
158
159! For forward simulations, output fields have dimension MAXSPEC,
160! for backward simulations, output fields have dimension MAXPOINT.
161! Thus, make loops either about nspec, or about numpoint
162!*****************************************************************
163
164
165  if (ldirect.eq.1) then
166    do ks=1,nspec
167      do kp=1,maxpointspec_act
168        tot_mu(ks,kp)=1
169      end do
170    end do
171  else
172    do ks=1,nspec
173      do kp=1,maxpointspec_act
174        tot_mu(ks,kp)=xmass(kp,ks)
175      end do
176    end do
177  endif
178
179
180!*******************************************************************
181! Compute air density: sufficiently accurate to take it
182! from coarse grid at some time
183! Determine center altitude of output layer, and interpolate density
184! data to that altitude
185!*******************************************************************
186
187  mind=memind(2)
188  do kz=1,numzgrid
189    if (kz.eq.1) then
190      halfheight=outheight(1)/2.
191    else
192      halfheight=(outheight(kz)+outheight(kz-1))/2.
193    endif
194    do kzz=2,nz
195      if ((height(kzz-1).lt.halfheight).and. &
196           (height(kzz).gt.halfheight)) goto 46
197    end do
19846  kzz=max(min(kzz,nz),2)
199    dz1=halfheight-height(kzz-1)
200    dz2=height(kzz)-halfheight
201    dz=dz1+dz2
202    do jy=0,numygrid-1
203      do ix=0,numxgrid-1
204        xl=outlon0+real(ix)*dxout
205        yl=outlat0+real(jy)*dyout
206        xl=(xl-xlon0)/dx
207        yl=(yl-ylat0)/dy !v9.1.1
208        iix=max(min(nint(xl),nxmin1),0)
209        jjy=max(min(nint(yl),nymin1),0)
210        ! densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,2)*dz1+ &
211        !      rho(iix,jjy,kzz-1,2)*dz2)/dz
212        densityoutgrid(ix,jy,kz)=(rho(iix,jjy,kzz,mind)*dz1+ &
213             rho(iix,jjy,kzz-1,mind)*dz2)/dz
214      end do
215    end do
216  end do
217
218  do i=1,numreceptor
219    xl=xreceptor(i)
220    yl=yreceptor(i)
221    iix=max(min(nint(xl),nxmin1),0)
222    jjy=max(min(nint(yl),nymin1),0)
223    !densityoutrecept(i)=rho(iix,jjy,1,2)
224    densityoutrecept(i)=rho(iix,jjy,1,mind)
225  end do
226
227
228! Output is different for forward and backward simulations
229  do kz=1,numzgrid
230    do jy=0,numygrid-1
231      do ix=0,numxgrid-1
232        if (ldirect.eq.1) then
233          factor3d(ix,jy,kz)=1.e12/volume(ix,jy,kz)/outnum
234        else
235          factor3d(ix,jy,kz)=real(abs(loutaver))/outnum
236        endif
237      end do
238    end do
239  end do
240
241!*********************************************************************
242! Determine the standard deviation of the mean concentration or mixing
243! ratio (uncertainty of the output) and the dry and wet deposition
244!*********************************************************************
245
246  gridtotal=0.
247  gridsigmatotal=0.
248  gridtotalunc=0.
249  wetgridtotal=0.
250  wetgridsigmatotal=0.
251  wetgridtotalunc=0.
252  drygridtotal=0.
253  drygridsigmatotal=0.
254  drygridtotalunc=0.
255
256  do ks=1,nspec
257
258    write(anspec,'(i3.3)') ks
259    if ((iout.eq.1).or.(iout.eq.3).or.(iout.eq.5)) then
260      if (ldirect.eq.1) then
261        open(unitoutgrid,file=path(2)(1:length(2))//'grid_conc_'//adate// &
262             atime//'_'//anspec,form='unformatted')
263      else
264        open(unitoutgrid,file=path(2)(1:length(2))//'grid_time_'//adate// &
265             atime//'_'//anspec,form='unformatted')
266      endif
267      write(unitoutgrid) itime
268    endif
269
270    if ((iout.eq.2).or.(iout.eq.3)) then      ! mixing ratio
271      open(unitoutgridppt,file=path(2)(1:length(2))//'grid_pptv_'//adate// &
272           atime//'_'//anspec,form='unformatted')
273
274      write(unitoutgridppt) itime
275    endif
276
277    do kp=1,maxpointspec_act
278      do nage=1,nageclass
279
280        do jy=0,numygrid-1
281          do ix=0,numxgrid-1
282
283! WET DEPOSITION
284            if ((WETDEP).and.(ldirect.gt.0)) then
285              do l=1,nclassunc
286                auxgrid(l)=wetgridunc0(ix,jy,ks,kp,l,nage)
287              end do
288              call mean(auxgrid,wetgrid(ix,jy), &
289                   wetgridsigma(ix,jy),nclassunc)
290! Multiply by number of classes to get total concentration
291              wetgrid(ix,jy)=wetgrid(ix,jy) &
292                   *nclassunc
293              wetgridtotal=wetgridtotal+wetgrid(ix,jy)
294! Calculate standard deviation of the mean
295              wetgridsigma(ix,jy)= &
296                   wetgridsigma(ix,jy)* &
297                   sqrt(real(nclassunc))
298              wetgridsigmatotal=wetgridsigmatotal+ &
299                   wetgridsigma(ix,jy)
300            endif
301
302! DRY DEPOSITION
303            if ((DRYDEP).and.(ldirect.gt.0)) then
304              do l=1,nclassunc
305                auxgrid(l)=drygridunc0(ix,jy,ks,kp,l,nage)
306              end do
307              call mean(auxgrid,drygrid(ix,jy), &
308                   drygridsigma(ix,jy),nclassunc)
309! Multiply by number of classes to get total concentration
310              drygrid(ix,jy)=drygrid(ix,jy)* &
311                   nclassunc
312              drygridtotal=drygridtotal+drygrid(ix,jy)
313! Calculate standard deviation of the mean
314              drygridsigma(ix,jy)= &
315                   drygridsigma(ix,jy)* &
316                   sqrt(real(nclassunc))
317              drygridsigmatotal=drygridsigmatotal+ &
318                   drygridsigma(ix,jy)
319            endif
320
321! CONCENTRATION OR MIXING RATIO
322            do kz=1,numzgrid
323              do l=1,nclassunc
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